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-2020 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 #include "demandDrivenData.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(twoDPointCorrector, 0);
42 }
43 
44 const Foam::scalar Foam::twoDPointCorrector::edgeOrthogonalityTol = 1.0 - 1e-4;
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::twoDPointCorrector::calcAddressing() const
50 {
51  // Find geometry normal
52  planeNormalPtr_ = new vector(0, 0, 0);
53  vector& pn = *planeNormalPtr_;
54 
55  // Algorithm:
56  // Attempt to find wedge patch and work out the normal from it.
57  // If not found, find an empty patch with faces in it and use the
58  // normal of the first face. If neither is found, declare an
59  // error.
60 
61  // Try and find a wedge patch
62  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
63 
64  for (const polyPatch& p : patches)
65  {
66  if (isA<wedgePolyPatch>(p))
67  {
68  isWedge_ = true;
69 
70  const wedgePolyPatch& wp = refCast<const wedgePolyPatch>(p);
71 
72  pn = wp.centreNormal();
73 
74  wedgeAxis_ = wp.axis();
75  wedgeAngle_ = mag(acos(wp.cosAngle()));
76 
77  if (polyMesh::debug)
78  {
79  Pout<< "Found normal from wedge patch " << p.index() << nl;
80  }
81 
82  break;
83  }
84  }
85 
86  // Try to find an empty patch with faces
87  if (!isWedge_)
88  {
89  for (const polyPatch& p : patches)
90  {
91  if (isA<emptyPolyPatch>(p) && p.size())
92  {
93  pn = p.faceAreas()[0];
94 
95  if (polyMesh::debug)
96  {
97  Pout<< "Found normal from empty patch " << p.index() << nl;
98  }
99 
100  break;
101  }
102  }
103  }
104 
105 
106  if (mag(pn) < VSMALL)
107  {
109  << "Cannot determine normal vector from patches."
110  << abort(FatalError);
111  }
112  else
113  {
114  pn /= mag(pn);
115  }
116 
117  if (polyMesh::debug)
118  {
119  Pout<< " twoDPointCorrector normal: " << pn << nl;
120  }
121 
122  // Select edges to be included in check.
123  normalEdgeIndicesPtr_ = new labelList(mesh_.nEdges());
124  labelList& neIndices = *normalEdgeIndicesPtr_;
125 
126  const edgeList& meshEdges = mesh_.edges();
127 
128  const pointField& meshPoints = mesh_.points();
129 
130  label nNormalEdges = 0;
131 
132  forAll(meshEdges, edgeI)
133  {
134  const edge& e = meshEdges[edgeI];
135 
136  vector edgeVector = e.unitVec(meshPoints);
137 
138  if (mag(edgeVector & pn) > edgeOrthogonalityTol)
139  {
140  // this edge is normal to the plane. Add it to the list
141  neIndices[nNormalEdges++] = edgeI;
142  }
143  }
144 
145  neIndices.setSize(nNormalEdges);
146 
147  // Construction check: number of points in a read 2-D or wedge geometry
148  // should be odd and the number of edges normal to the plane should be
149  // exactly half the number of points
150  if (!isWedge_)
151  {
152  if (meshPoints.size() % 2 != 0)
153  {
155  << "the number of vertices in the geometry "
156  << "is odd - this should not be the case for a 2-D case. "
157  << "Please check the geometry."
158  << endl;
159  }
160 
161  if (2*nNormalEdges != meshPoints.size())
162  {
164  << "The number of points in the mesh is "
165  << "not equal to twice the number of edges normal to the plane "
166  << "- this may be OK only for wedge geometries.\n"
167  << " Please check the geometry or adjust "
168  << "the orthogonality tolerance.\n" << endl
169  << "Number of normal edges: " << nNormalEdges
170  << " number of points: " << meshPoints.size()
171  << endl;
172  }
173  }
174 }
175 
176 
177 void Foam::twoDPointCorrector::clearAddressing() const
178 {
179  deleteDemandDrivenData(planeNormalPtr_);
180  deleteDemandDrivenData(normalEdgeIndicesPtr_);
181 }
182 
183 
184 void Foam::twoDPointCorrector::snapToWedge
185 (
186  const vector& n,
187  const point& A,
188  point& p
189 ) const
190 {
191  scalar ADash = mag(A - wedgeAxis_*(wedgeAxis_ & A));
192  vector pDash = ADash*tan(wedgeAngle_)*planeNormal();
193 
194  p = A + sign(n & p)*pDash;
195 }
196 
197 
198 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
199 
200 Foam::twoDPointCorrector::twoDPointCorrector(const polyMesh& mesh)
201 :
203  required_(mesh_.nGeometricD() == 2),
204  planeNormalPtr_(nullptr),
205  normalEdgeIndicesPtr_(nullptr),
206  isWedge_(false),
207  wedgeAxis_(Zero),
208  wedgeAngle_(0.0)
209 {}
210 
211 
212 
213 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
214 
216 {
217  clearAddressing();
218 }
219 
220 
221 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222 
224 {
225  const vector& pn = planeNormal();
226 
227  if (mag(pn.x()) >= edgeOrthogonalityTol)
228  {
229  return vector::X;
230  }
231  else if (mag(pn.y()) >= edgeOrthogonalityTol)
232  {
233  return vector::Y;
234  }
235  else if (mag(pn.z()) >= edgeOrthogonalityTol)
236  {
237  return vector::Z;
238  }
239 
241  << "Plane normal not aligned with the coordinate system" << nl
242  << " pn = " << pn
243  << abort(FatalError);
244 
245  return vector::Z;
246 }
247 
248 
250 {
251  if (!planeNormalPtr_)
252  {
253  calcAddressing();
254  }
255 
256  return *planeNormalPtr_;
257 }
258 
259 
261 {
262  if (!normalEdgeIndicesPtr_)
263  {
264  calcAddressing();
265  }
266 
267  return *normalEdgeIndicesPtr_;
268 }
269 
270 
272 {
273  if (!required_) return;
274 
275  // Algorithm:
276  // Loop through all edges. Calculate the average point position A for
277  // the front and the back. Correct the position of point P (in two planes)
278  // such that vectors AP and planeNormal are parallel
279 
280  // Get reference to edges
281  const edgeList& meshEdges = mesh_.edges();
282 
283  const labelList& neIndices = normalEdgeIndices();
284  const vector& pn = planeNormal();
285 
286  for (const label edgei : neIndices)
287  {
288  point& pStart = p[meshEdges[edgei].start()];
289 
290  point& pEnd = p[meshEdges[edgei].end()];
291 
292  // calculate average point position
293  point A = 0.5*(pStart + pEnd);
295 
296  if (isWedge_)
297  {
298  snapToWedge(pn, A, pStart);
299  snapToWedge(pn, A, pEnd);
300  }
301  else
302  {
303  // correct point locations
304  pStart = A + pn*(pn & (pStart - A));
305  pEnd = A + pn*(pn & (pEnd - A));
306  }
307  }
308 }
309 
310 
312 (
313  const pointField& p,
314  vectorField& disp
315 ) const
316 {
317  if (!required_) return;
318 
319  // Algorithm:
320  // Loop through all edges. Calculate the average point position A for
321  // the front and the back. Correct the position of point P (in two planes)
322  // such that vectors AP and planeNormal are parallel
323 
324  // Get reference to edges
325  const edgeList& meshEdges = mesh_.edges();
326 
327  const labelList& neIndices = normalEdgeIndices();
328  const vector& pn = planeNormal();
329 
330  for (const label edgei : neIndices)
331  {
332  const edge& e = meshEdges[edgei];
333 
334  label startPointi = e.start();
335  point pStart = p[startPointi] + disp[startPointi];
336 
337  label endPointi = e.end();
338  point pEnd = p[endPointi] + disp[endPointi];
339 
340  // calculate average point position
341  point A = 0.5*(pStart + pEnd);
343 
344  if (isWedge_)
345  {
346  snapToWedge(pn, A, pStart);
347  snapToWedge(pn, A, pEnd);
348  disp[startPointi] = pStart - p[startPointi];
349  disp[endPointi] = pEnd - p[endPointi];
350  }
351  else
352  {
353  // correct point locations
354  disp[startPointi] = (A + pn*(pn & (pStart - A))) - p[startPointi];
355  disp[endPointi] = (A + pn*(pn & (pEnd - A))) - p[endPointi];
356  }
357  }
358 }
359 
360 
362 {
363  clearAddressing();
364 }
365 
366 
368 {
369  return true;
370 }
371 
372 
373 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
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:1069
meshTools.H
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Vector< scalar >::Z
Definition: Vector.H:81
Foam::tan
dimensionedScalar tan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:266
Foam::Vector< scalar >::Y
Definition: Vector.H:81
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:249
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::twoDPointCorrector::correctDisplacement
void correctDisplacement(const pointField &p, vectorField &disp) const
Correct motion displacements.
Definition: twoDPointCorrector.C:312
wedgePolyPatch.H
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::twoDPointCorrector::~twoDPointCorrector
~twoDPointCorrector()
Destructor.
Definition: twoDPointCorrector.C:215
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:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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:85
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:296
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::Field< vector >
Foam::Vector< scalar >::X
Definition: Vector.H:81
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:144
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:79
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
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:52
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:161
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:367
Foam::twoDPointCorrector::normalEdgeIndices
const labelList & normalEdgeIndices() const
Return indices of normal edges.
Definition: twoDPointCorrector.C:260
Foam::twoDPointCorrector::correctPoints
void correctPoints(pointField &p) const
Correct motion points.
Definition: twoDPointCorrector.C:271
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:223
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::twoDPointCorrector::updateMesh
void updateMesh(const mapPolyMesh &)
Update topology.
Definition: twoDPointCorrector.C:361