twoDPointCorrector.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "twoDPointCorrector.H"
30 #include "polyMesh.H"
31 #include "wedgePolyPatch.H"
32 #include "emptyPolyPatch.H"
33 #include "SubField.H"
34 #include "meshTools.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(twoDPointCorrector, 0);
41 }
42 
43 const Foam::scalar Foam::twoDPointCorrector::edgeOrthogonalityTol = 1.0 - 1e-4;
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::twoDPointCorrector::calcAddressing() const
49 {
50  // Find geometry normal
51  planeNormalPtr_ = new vector(0, 0, 0);
52  vector& pn = *planeNormalPtr_;
53 
54  // Algorithm:
55  // Attempt to find wedge patch and work out the normal from it.
56  // If not found, find an empty patch with faces in it and use the
57  // normal of the first face. If neither is found, declare an
58  // error.
59 
60  // Try and find a wedge patch
61  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
62 
63  for (const polyPatch& p : patches)
64  {
65  if (isA<wedgePolyPatch>(p))
66  {
67  isWedge_ = true;
68 
69  const wedgePolyPatch& wp = refCast<const wedgePolyPatch>(p);
70 
71  pn = wp.centreNormal();
72 
73  wedgeAxis_ = wp.axis();
74  wedgeAngle_ = mag(acos(wp.cosAngle()));
75 
76  if (polyMesh::debug)
77  {
78  Pout<< "Found normal from wedge patch " << p.index() << nl;
79  }
80 
81  break;
82  }
83  }
84 
85  // Try to find an empty patch with faces
86  if (!isWedge_)
87  {
88  for (const polyPatch& p : patches)
89  {
90  if (isA<emptyPolyPatch>(p) && p.size())
91  {
92  pn = p.faceAreas()[0];
93 
94  if (polyMesh::debug)
95  {
96  Pout<< "Found normal from empty patch " << p.index() << nl;
97  }
98 
99  break;
100  }
101  }
102  }
103 
104 
105  if (mag(pn) < VSMALL)
106  {
108  << "Cannot determine normal vector from patches."
109  << abort(FatalError);
110  }
111  else
112  {
113  pn /= mag(pn);
114  }
115 
116  if (polyMesh::debug)
117  {
118  Pout<< " twoDPointCorrector normal: " << pn << nl;
119  }
120 
121  // Select edges to be included in check.
122  normalEdgeIndicesPtr_ = new labelList(mesh_.nEdges());
123  labelList& neIndices = *normalEdgeIndicesPtr_;
124 
125  const edgeList& meshEdges = mesh_.edges();
126 
127  const pointField& meshPoints = mesh_.points();
128 
129  label nNormalEdges = 0;
130 
131  forAll(meshEdges, edgeI)
132  {
133  const edge& e = meshEdges[edgeI];
134 
135  vector edgeVector = e.unitVec(meshPoints);
136 
137  if (mag(edgeVector & pn) > edgeOrthogonalityTol)
138  {
139  // this edge is normal to the plane. Add it to the list
140  neIndices[nNormalEdges++] = edgeI;
141  }
142  }
143 
144  neIndices.setSize(nNormalEdges);
145 
146  // Construction check: number of points in a read 2-D or wedge geometry
147  // should be odd and the number of edges normal to the plane should be
148  // exactly half the number of points
149  if (!isWedge_)
150  {
151  if (meshPoints.size() % 2 != 0)
152  {
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."
157  << endl;
158  }
159 
160  if (2*nNormalEdges != meshPoints.size())
161  {
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()
170  << endl;
171  }
172  }
173 }
174 
175 
176 void Foam::twoDPointCorrector::clearAddressing() const
177 {
178  deleteDemandDrivenData(planeNormalPtr_);
179  deleteDemandDrivenData(normalEdgeIndicesPtr_);
180 }
181 
182 
183 void Foam::twoDPointCorrector::snapToWedge
184 (
185  const vector& n,
186  const point& A,
187  point& p
188 ) const
189 {
190  scalar ADash = mag(A - wedgeAxis_*(wedgeAxis_ & A));
191  vector pDash = ADash*tan(wedgeAngle_)*planeNormal();
192 
193  p = A + sign(n & p)*pDash;
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
198 
199 Foam::twoDPointCorrector::twoDPointCorrector(const polyMesh& mesh)
200 :
202  required_(mesh_.nGeometricD() == 2),
203  planeNormalPtr_(nullptr),
204  normalEdgeIndicesPtr_(nullptr),
205  isWedge_(false),
206  wedgeAxis_(Zero),
207  wedgeAngle_(0.0)
208 {}
209 
210 
211 
212 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
213 
215 {
216  clearAddressing();
217 }
218 
219 
220 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 
223 {
224  const vector& pn = planeNormal();
225 
226  if (mag(pn.x()) >= edgeOrthogonalityTol)
227  {
228  return vector::X;
229  }
230  else if (mag(pn.y()) >= edgeOrthogonalityTol)
231  {
232  return vector::Y;
233  }
234  else if (mag(pn.z()) >= edgeOrthogonalityTol)
235  {
236  return vector::Z;
237  }
238 
240  << "Plane normal not aligned with the coordinate system" << nl
241  << " pn = " << pn
242  << abort(FatalError);
243 
244  return vector::Z;
245 }
246 
247 
249 {
250  if (!planeNormalPtr_)
251  {
252  calcAddressing();
253  }
254 
255  return *planeNormalPtr_;
256 }
257 
258 
260 {
261  if (!normalEdgeIndicesPtr_)
262  {
263  calcAddressing();
264  }
265 
266  return *normalEdgeIndicesPtr_;
267 }
268 
269 
271 {
272  if (!required_) return;
273 
274  // Algorithm:
275  // Loop through all edges. Calculate the average point position A for
276  // the front and the back. Correct the position of point P (in two planes)
277  // such that vectors AP and planeNormal are parallel
278 
279  // Get reference to edges
280  const edgeList& meshEdges = mesh_.edges();
281 
282  const labelList& neIndices = normalEdgeIndices();
283  const vector& pn = planeNormal();
284 
285  for (const label edgei : neIndices)
286  {
287  point& pStart = p[meshEdges[edgei].start()];
288 
289  point& pEnd = p[meshEdges[edgei].end()];
290 
291  // calculate average point position
292  point A = 0.5*(pStart + pEnd);
294 
295  if (isWedge_)
296  {
297  snapToWedge(pn, A, pStart);
298  snapToWedge(pn, A, pEnd);
299  }
300  else
301  {
302  // correct point locations
303  pStart = A + pn*(pn & (pStart - A));
304  pEnd = A + pn*(pn & (pEnd - A));
305  }
306  }
307 }
308 
309 
311 (
312  const pointField& p,
313  vectorField& disp
314 ) const
315 {
316  if (!required_) return;
317 
318  // Algorithm:
319  // Loop through all edges. Calculate the average point position A for
320  // the front and the back. Correct the position of point P (in two planes)
321  // such that vectors AP and planeNormal are parallel
322 
323  // Get reference to edges
324  const edgeList& meshEdges = mesh_.edges();
325 
326  const labelList& neIndices = normalEdgeIndices();
327  const vector& pn = planeNormal();
328 
329  for (const label edgei : neIndices)
330  {
331  const edge& e = meshEdges[edgei];
332 
333  label startPointi = e.start();
334  point pStart = p[startPointi] + disp[startPointi];
335 
336  label endPointi = e.end();
337  point pEnd = p[endPointi] + disp[endPointi];
338 
339  // calculate average point position
340  point A = 0.5*(pStart + pEnd);
342 
343  if (isWedge_)
344  {
345  snapToWedge(pn, A, pStart);
346  snapToWedge(pn, A, pEnd);
347  disp[startPointi] = pStart - p[startPointi];
348  disp[endPointi] = pEnd - p[endPointi];
349  }
350  else
351  {
352  // correct point locations
353  disp[startPointi] = (A + pn*(pn & (pStart - A))) - p[startPointi];
354  disp[endPointi] = (A + pn*(pn & (pEnd - A))) - p[endPointi];
355  }
356  }
357 }
358 
359 
361 {
362  clearAddressing();
363 }
364 
365 
367 {
368  return true;
369 }
370 
371 
372 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
meshTools.H
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:78
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Vector< scalar >::Z
Definition: Vector.H:80
Foam::tan
dimensionedScalar tan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:266
Foam::Vector< scalar >::Y
Definition: Vector.H:80
Foam::edgeList
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
SubField.H
Foam::twoDPointCorrector::planeNormal
const vector & planeNormal() const
Return plane normal.
Definition: twoDPointCorrector.C:248
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::twoDPointCorrector::correctDisplacement
void correctDisplacement(const pointField &p, vectorField &disp) const
Correct motion displacements.
Definition: twoDPointCorrector.C:311
wedgePolyPatch.H
Foam::twoDPointCorrector::~twoDPointCorrector
~twoDPointCorrector()
Destructor.
Definition: twoDPointCorrector.C:214
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::primitiveMesh::nEdges
label nEdges() const
Number of mesh edges.
Definition: primitiveMeshI.H:67
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
polyMesh.H
Foam::UpdateableMeshObject
Definition: MeshObject.H:241
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:90
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::twoDPointCorrector
Class applies a two-dimensional correction to mesh motion point field.
Definition: twoDPointCorrector.H:63
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< vector >
Foam::Vector< scalar >::X
Definition: Vector.H:80
Foam::MeshObject< polyMesh, UpdateableMeshObject, twoDPointCorrector >::mesh_
const polyMesh & mesh_
Definition: MeshObject.H:96
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
emptyPolyPatch.H
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:84
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::Vector< scalar >
Foam::List< label >
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::direction
uint8_t direction
Definition: direction.H:47
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::meshTools::constrainToMeshCentre
void constrainToMeshCentre(const polyMesh &mesh, point &pt)
Set the constrained components of position to mesh centre.
Definition: meshTools.C:629
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
Foam::MeshObject
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:88
twoDPointCorrector.H
Foam::twoDPointCorrector::movePoints
bool movePoints()
Correct weighting factors for moving mesh.
Definition: twoDPointCorrector.C:366
Foam::twoDPointCorrector::normalEdgeIndices
const labelList & normalEdgeIndices() const
Return indices of normal edges.
Definition: twoDPointCorrector.C:259
Foam::twoDPointCorrector::correctPoints
void correctPoints(pointField &p) const
Correct motion points.
Definition: twoDPointCorrector.C:270
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::twoDPointCorrector::normalDir
direction normalDir() const
Return direction normal to plane.
Definition: twoDPointCorrector.C:222
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::twoDPointCorrector::updateMesh
void updateMesh(const mapPolyMesh &)
Update topology.
Definition: twoDPointCorrector.C:360