Loading...
Searching...
No Matches
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-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
38namespace Foam
39{
41}
42
43const Foam::scalar Foam::twoDPointCorrector::edgeOrthogonalityTol = 1.0 - 1e-4;
44
45
46// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
48void Foam::twoDPointCorrector::calcAddressing() const
49{
50 // Find geometry normal
51 planeNormalPtr_ = std::make_unique<vector>(0, 0, 0);
52 auto& 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 {
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_ = std::make_unique<labelList>(mesh_.nEdges());
123 auto& 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
176void Foam::twoDPointCorrector::clearAddressing() const
177{
178 planeNormalPtr_.reset(nullptr);
179 normalEdgeIndicesPtr_.reset(nullptr);
180}
181
182
183void 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();
193 p = A + sign(n & p)*pDash;
194}
195
196
197// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
198
199Foam::twoDPointCorrector::twoDPointCorrector(const polyMesh& mesh)
200:
201 MeshObject_type(mesh),
202 required_(mesh_.nGeometricD() == 2),
203 isWedge_(false),
204 wedgeAxis_(Zero),
205 wedgeAngle_(0.0)
206{}
207
208
209// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
210
213 clearAddressing();
214}
215
216
217// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218
220{
221 const vector& pn = planeNormal();
222
223 if (mag(pn.x()) >= edgeOrthogonalityTol)
224 {
225 return vector::X;
226 }
227 else if (mag(pn.y()) >= edgeOrthogonalityTol)
228 {
229 return vector::Y;
230 }
231 else if (mag(pn.z()) >= edgeOrthogonalityTol)
232 {
233 return vector::Z;
234 }
235
237 << "Plane normal not aligned with the coordinate system" << nl
238 << " pn = " << pn
239 << abort(FatalError);
240
241 return vector::Z;
242}
243
244
246{
247 if (!planeNormalPtr_)
248 {
249 calcAddressing();
250 }
251
252 return *planeNormalPtr_;
253}
254
255
257{
258 if (!normalEdgeIndicesPtr_)
259 {
260 calcAddressing();
261 }
262
263 return *normalEdgeIndicesPtr_;
264}
265
266
268{
269 if (!required_) return;
270
271 // Algorithm:
272 // Loop through all edges. Calculate the average point position A for
273 // the front and the back. Correct the position of point P (in two planes)
274 // such that vectors AP and planeNormal are parallel
275
276 // Get reference to edges
277 const edgeList& meshEdges = mesh_.edges();
278
279 const labelList& neIndices = normalEdgeIndices();
280 const vector& pn = planeNormal();
281
282 for (const label edgei : neIndices)
283 {
284 point& pStart = p[meshEdges[edgei].start()];
285
286 point& pEnd = p[meshEdges[edgei].end()];
287
288 // calculate average point position
289 point A = 0.5*(pStart + pEnd);
291
292 if (isWedge_)
293 {
294 snapToWedge(pn, A, pStart);
295 snapToWedge(pn, A, pEnd);
296 }
297 else
298 {
299 // correct point locations
300 pStart = A + pn*(pn & (pStart - A));
301 pEnd = A + pn*(pn & (pEnd - A));
302 }
303 }
304}
305
306
308(
309 const pointField& p,
310 vectorField& disp
311) const
312{
313 if (!required_) return;
314
315 // Algorithm:
316 // Loop through all edges. Calculate the average point position A for
317 // the front and the back. Correct the position of point P (in two planes)
318 // such that vectors AP and planeNormal are parallel
319
320 // Get reference to edges
321 const edgeList& meshEdges = mesh_.edges();
322
323 const labelList& neIndices = normalEdgeIndices();
324 const vector& pn = planeNormal();
325
326 for (const label edgei : neIndices)
327 {
328 const edge& e = meshEdges[edgei];
329
330 label startPointi = e.start();
331 point pStart = p[startPointi] + disp[startPointi];
332
333 label endPointi = e.end();
334 point pEnd = p[endPointi] + disp[endPointi];
335
336 // calculate average point position
337 point A = 0.5*(pStart + pEnd);
339
340 if (isWedge_)
341 {
342 snapToWedge(pn, A, pStart);
343 snapToWedge(pn, A, pEnd);
344 disp[startPointi] = pStart - p[startPointi];
345 disp[endPointi] = pEnd - p[endPointi];
346 }
347 else
348 {
349 // correct point locations
350 disp[startPointi] = (A + pn*(pn & (pStart - A))) - p[startPointi];
351 disp[endPointi] = (A + pn*(pn & (pEnd - A))) - p[endPointi];
352 }
353 }
354}
355
358{
359 clearAddressing();
360}
361
362
364{
365 return true;
366}
367
368
369// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
label n
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition UListI.H:454
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
Class applies a two-dimensional correction to mesh motion point field.
bool movePoints()
Correct weighting factors for moving mesh.
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.
Definition className.H:142
volScalarField & p
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
#define WarningInFunction
Report a warning using Foam::Warning.
void constrainToMeshCentre(const polyMesh &mesh, point &pt)
Set the constrained components of position to mesh centre.
Definition meshTools.C:622
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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.
Vector< scalar > vector
Definition vector.H:57
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299