43const Foam::scalar Foam::twoDPointCorrector::edgeOrthogonalityTol = 1.0 - 1
e-4;
48void Foam::twoDPointCorrector::calcAddressing()
const
51 planeNormalPtr_ = std::make_unique<vector>(0, 0, 0);
52 auto& pn = *planeNormalPtr_;
71 pn = wp.centreNormal();
73 wedgeAxis_ = wp.axis();
74 wedgeAngle_ =
mag(
acos(wp.cosAngle()));
78 Pout<<
"Found normal from wedge patch " <<
p.index() <<
nl;
92 pn =
p.faceAreas()[0];
96 Pout<<
"Found normal from empty patch " <<
p.index() <<
nl;
105 if (
mag(pn) < VSMALL)
108 <<
"Cannot determine normal vector from patches."
118 Pout<<
" twoDPointCorrector normal: " << pn <<
nl;
122 normalEdgeIndicesPtr_ = std::make_unique<labelList>(
mesh_.nEdges());
123 auto& neIndices = *normalEdgeIndicesPtr_;
129 label nNormalEdges = 0;
133 const edge&
e = meshEdges[edgeI];
135 vector edgeVector =
e.unitVec(meshPoints);
137 if (
mag(edgeVector & pn) > edgeOrthogonalityTol)
140 neIndices[nNormalEdges++] = edgeI;
144 neIndices.setSize(nNormalEdges);
151 if (meshPoints.size() % 2 != 0)
154 <<
"the number of vertices in the geometry "
155 <<
"is odd - this should not be the case for a 2-D case. "
156 <<
"Please check the geometry."
160 if (2*nNormalEdges != meshPoints.size())
163 <<
"The number of points in the mesh is "
164 <<
"not equal to twice the number of edges normal to the plane "
165 <<
"- this may be OK only for wedge geometries.\n"
166 <<
" Please check the geometry or adjust "
167 <<
"the orthogonality tolerance.\n" <<
endl
168 <<
"Number of normal edges: " << nNormalEdges
169 <<
" number of points: " << meshPoints.size()
176void Foam::twoDPointCorrector::clearAddressing()
const
178 planeNormalPtr_.reset(
nullptr);
179 normalEdgeIndicesPtr_.reset(
nullptr);
183void Foam::twoDPointCorrector::snapToWedge
190 scalar ADash =
mag(
A - wedgeAxis_*(wedgeAxis_ &
A));
191 vector pDash = ADash*
tan(wedgeAngle_)*planeNormal();
199Foam::twoDPointCorrector::twoDPointCorrector(
const polyMesh&
mesh)
201 MeshObject_type(
mesh),
202 required_(mesh_.nGeometricD() == 2),
221 const vector& pn = planeNormal();
223 if (
mag(pn.
x()) >= edgeOrthogonalityTol)
227 else if (
mag(pn.
y()) >= edgeOrthogonalityTol)
231 else if (
mag(pn.
z()) >= edgeOrthogonalityTol)
237 <<
"Plane normal not aligned with the coordinate system" <<
nl
247 if (!planeNormalPtr_)
252 return *planeNormalPtr_;
258 if (!normalEdgeIndicesPtr_)
263 return *normalEdgeIndicesPtr_;
269 if (!required_)
return;
277 const edgeList& meshEdges = mesh_.edges();
279 const labelList& neIndices = normalEdgeIndices();
280 const vector& pn = planeNormal();
282 for (
const label edgei : neIndices)
284 point& pStart =
p[meshEdges[edgei].start()];
286 point& pEnd =
p[meshEdges[edgei].end()];
289 point A = 0.5*(pStart + pEnd);
294 snapToWedge(pn,
A, pStart);
295 snapToWedge(pn,
A, pEnd);
300 pStart =
A + pn*(pn & (pStart -
A));
301 pEnd =
A + pn*(pn & (pEnd -
A));
313 if (!required_)
return;
321 const edgeList& meshEdges = mesh_.edges();
323 const labelList& neIndices = normalEdgeIndices();
324 const vector& pn = planeNormal();
326 for (
const label edgei : neIndices)
328 const edge&
e = meshEdges[edgei];
330 label startPointi =
e.start();
331 point pStart =
p[startPointi] + disp[startPointi];
333 label endPointi =
e.
end();
334 point pEnd =
p[endPointi] + disp[endPointi];
337 point A = 0.5*(pStart + pEnd);
342 snapToWedge(pn,
A, pStart);
343 snapToWedge(pn,
A, pEnd);
344 disp[startPointi] = pStart -
p[startPointi];
345 disp[endPointi] = pEnd -
p[endPointi];
350 disp[startPointi] = (
A + pn*(pn & (pStart -
A))) -
p[startPointi];
351 disp[endPointi] = (
A + pn*(pn & (pEnd -
A))) -
p[endPointi];
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
const polyMesh & mesh() const noexcept
iterator end() noexcept
Return an iterator to end traversing the UList.
const Cmpt & x() const noexcept
Access to the vector x component.
const Cmpt & z() const noexcept
Access to the vector z component.
const Cmpt & y() const noexcept
Access to the vector y component.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Mesh consisting of general polyhedral cells.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Class applies a two-dimensional correction to mesh motion point field.
bool movePoints()
Correct weighting factors for moving mesh.
~twoDPointCorrector()
Destructor.
void correctDisplacement(const pointField &p, vectorField &disp) const
Correct motion displacements.
direction normalDir() const
Return direction normal to plane.
void correctPoints(pointField &p) const
Correct motion points.
void updateMesh(const mapPolyMesh &)
Update topology.
const labelList & normalEdgeIndices() const
Return indices of normal edges.
const vector & planeNormal() const
Return plane normal.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const polyBoundaryMesh & patches
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
List< edge > edgeList
List of edge.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.