40 const Foam::scalar Foam::twoDPointCorrector::edgeOrthogonalityTol = 1.0 - 1
e-4;
45 void Foam::twoDPointCorrector::calcAddressing()
const
48 planeNormalPtr_ =
new vector(0, 0, 0);
49 vector& pn = *planeNormalPtr_;
66 const wedgePolyPatch& wp =
69 pn = wp.centreNormal();
71 wedgeAxis_ = wp.axis();
72 wedgeAngle_ =
mag(
acos(wp.cosAngle()));
76 Pout<<
"Found normal from wedge patch " <<
patchi;
94 Pout<<
"Found normal from empty patch " <<
patchi;
103 if (
mag(pn) < vSmall)
106 <<
"Cannot determine normal vector from patches."
116 Pout<<
" twoDPointCorrector normal: " << pn <<
endl;
121 labelList& neIndices = *normalEdgeIndicesPtr_;
127 label nNormalEdges = 0;
131 const edge&
e = meshEdges[edgeI];
133 vector edgeVector =
e.vec(meshPoints)/(
e.mag(meshPoints) + vSmall);
135 if (
mag(edgeVector & pn) > edgeOrthogonalityTol)
138 neIndices[nNormalEdges++] = edgeI;
142 neIndices.setSize(nNormalEdges);
149 if (meshPoints.size() % 2 != 0)
152 <<
"the number of vertices in the geometry "
153 <<
"is odd - this should not be the case for a 2-D case. "
154 <<
"Please check the geometry."
158 if (2*nNormalEdges != meshPoints.size())
161 <<
"The number of points in the mesh is "
162 <<
"not equal to twice the number of edges normal to the plane "
163 <<
"- this may be OK only for wedge geometries.\n"
164 <<
" Please check the geometry or adjust "
165 <<
"the orthogonality tolerance.\n" <<
endl
166 <<
"Number of normal edges: " << nNormalEdges
167 <<
" number of points: " << meshPoints.size()
174 void Foam::twoDPointCorrector::clearAddressing()
const
181 void Foam::twoDPointCorrector::snapToWedge
188 scalar ADash =
mag(
A - wedgeAxis_*(wedgeAxis_ &
A));
189 vector pDash = ADash*
tan(wedgeAngle_)*planeNormal();
205 required_(mesh.nGeometricD() == 2),
206 planeNormalPtr_(nullptr),
207 normalEdgeIndicesPtr_(nullptr),
227 const vector& pn = planeNormal();
229 if (
mag(pn.
x()) >= edgeOrthogonalityTol)
233 else if (
mag(pn.
y()) >= edgeOrthogonalityTol)
237 else if (
mag(pn.
z()) >= edgeOrthogonalityTol)
244 <<
"Plane normal not aligned with the coordinate system" <<
nl
255 if (!planeNormalPtr_)
260 return *planeNormalPtr_;
266 if (!normalEdgeIndicesPtr_)
271 return *normalEdgeIndicesPtr_;
277 if (!required_)
return;
285 const edgeList& meshEdges = mesh().edges();
287 const labelList& neIndices = normalEdgeIndices();
288 const vector& pn = planeNormal();
292 point& pStart =
p[meshEdges[neIndices[edgeI]].start()];
294 point& pEnd =
p[meshEdges[neIndices[edgeI]].
end()];
297 point A = 0.5*(pStart + pEnd);
302 snapToWedge(pn,
A, pStart);
303 snapToWedge(pn,
A, pEnd);
308 pStart =
A + pn*(pn & (pStart -
A));
309 pEnd =
A + pn*(pn & (pEnd -
A));
321 if (!required_)
return;
329 const edgeList& meshEdges = mesh().edges();
331 const labelList& neIndices = normalEdgeIndices();
332 const vector& pn = planeNormal();
336 const edge&
e = meshEdges[neIndices[edgeI]];
338 label startPointi =
e.start();
339 point pStart =
p[startPointi] + disp[startPointi];
341 label endPointi =
e.end();
342 point pEnd =
p[endPointi] + disp[endPointi];
345 point A = 0.5*(pStart + pEnd);
350 snapToWedge(pn,
A, pStart);
351 snapToWedge(pn,
A, pEnd);
352 disp[startPointi] = pStart -
p[startPointi];
353 disp[endPointi] = pEnd -
p[endPointi];
358 disp[startPointi] = (
A + pn*(pn & (pStart -
A))) -
p[startPointi];
359 disp[endPointi] = (
A + pn*(pn & (pEnd -
A))) -
p[endPointi];
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
#define forAll(list, i)
Loop across all elements in list.
Templated abstract base-class for demand-driven mesh objects used to automate their allocation to the...
const polyMesh & mesh() const
iterator end()
Return an iterator to end traversing the UList.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Mesh consisting of general polyhedral cells.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
virtual const pointField & points() const
Return raw points.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Class applies a two-dimensional correction to mesh motion point field.
virtual bool movePoints()
Correct weighting factors for moving mesh.
virtual void topoChange(const polyTopoChangeMap &)
Update topology.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
~twoDPointCorrector()
Destructor.
void correctDisplacement(const pointField &p, vectorField &disp) const
Correct motion displacements.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
twoDPointCorrector(const polyMesh &mesh)
Construct from mesh.
direction normalDir() const
Return direction normal to plane.
void correctPoints(pointField &p) const
Correct motion points.
const labelList & normalEdgeIndices() const
Return indices of normal edges.
const vector & planeNormal() const
Return plane normal.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar e
Elementary charge.
dimensionedScalar tan(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
dimensionedScalar sign(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void deleteDemandDrivenData(DataType *&dataPtr)
errorManip< error > abort(error &err)
vectorField pointField
pointField is a vectorField.
vector point
Point is a vector.
Vector< scalar > vector
A scalar version of the templated Vector.
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout, "Pout")
dimensionedScalar acos(const dimensionedScalar &ds)