Loading...
Searching...
No Matches
parallelFvGeometryScheme.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) 2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "fvMesh.H"
31#include "syncTools.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
41 (
44 dict
45 );
46}
47
48// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49
50void Foam::parallelFvGeometryScheme::adjustGeometry
51(
52 pointField& faceCentres,
53 vectorField& faceAreas
54) const
55{
56 pointField syncedBCentres;
57 syncedBCentres = SubField<vector>
58 (
59 faceCentres,
62 );
64 (
65 mesh_,
66 syncedBCentres
67 );
68
69 vectorField syncedBAreas;
70 syncedBAreas = SubField<vector>
71 (
72 faceAreas,
75 );
77 (
78 mesh_,
79 syncedBAreas,
80 eqOp<vector>(),
81 transformOriented()
82 );
83
84 const auto& pbm = mesh_.boundaryMesh();
85 for (const auto& pp : pbm)
86 {
87 const auto* ppp = isA<coupledPolyPatch>(pp);
88
89 //if (ppp)
90 //{
91 // Pout<< "For patch:" << ppp->name()
92 // << " size:" << ppp->size()
93 // << endl;
94 // forAll(*ppp, i)
95 // {
96 // const label facei = ppp->start()+i;
97 // Pout<< " Face:" << facei << nl
98 // << " meshFc:" << faceCentres[facei] << nl
99 // << " meshFa:" << faceAreas[facei] << nl
100 // << endl;
101 // }
102 //}
103
104 if (ppp && !ppp->owner())
105 {
106 SubField<point> patchFc
107 (
108 faceCentres,
109 ppp->size(),
110 ppp->start()
111 );
112 patchFc = SubField<vector>
113 (
114 syncedBCentres,
115 ppp->size(),
116 ppp->offset()
117 );
118
119 SubField<vector> patchArea
120 (
121 faceAreas,
122 ppp->size(),
123 ppp->start()
124 );
125 patchArea = SubField<vector>
126 (
127 syncedBAreas,
128 ppp->size(),
129 ppp->offset()
130 );
131 }
132 }
133}
134
135
136void Foam::parallelFvGeometryScheme::adjustGeometry()
137{
138 // Swap face centres and areas
139 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
140
141 pointField faceCentres(mesh_.faceCentres());
142 vectorField faceAreas(mesh_.faceAreas());
143 {
144 pointField syncedBCentres;
145 syncedBCentres = SubField<vector>
146 (
147 mesh_.faceCentres(),
148 mesh_.nBoundaryFaces(),
149 mesh_.nInternalFaces()
150 );
152 (
153 mesh_,
154 syncedBCentres
155 );
156
157 vectorField syncedBAreas;
158 syncedBAreas = SubField<vector>
159 (
160 mesh_.faceAreas(),
161 mesh_.nBoundaryFaces(),
162 mesh_.nInternalFaces()
163 );
165 (
166 mesh_,
167 syncedBAreas,
168 eqOp<vector>(),
169 transformOriented()
170 );
171
172 const auto& pbm = mesh_.boundaryMesh();
173 for (const auto& pp : pbm)
174 {
175 const auto* ppp = isA<coupledPolyPatch>(pp);
176
177 //if (ppp)
178 //{
179 // Pout<< "For patch:" << ppp->name()
180 // << " size:" << ppp->size()
181 // << endl;
182 // forAll(*ppp, i)
183 // {
184 // const label facei = ppp->start()+i;
185 // Pout<< " Face:" << facei << nl
186 // << " meshFc:" << faceCentres[facei] << nl
187 // << " meshFa:" << faceAreas[facei] << nl
188 // << endl;
189 // }
190 //}
191
192 if (ppp && !ppp->owner())
193 {
194 // Delete old cached geometry
195 const_cast<coupledPolyPatch&>
196 (
197 *ppp
199
200 SubField<point> patchFc
201 (
202 faceCentres,
203 ppp->size(),
204 ppp->start()
205 );
206 patchFc = SubField<vector>
207 (
208 syncedBCentres,
209 ppp->size(),
210 ppp->offset()
211 );
212
213 SubField<vector> patchArea
214 (
215 faceAreas,
216 ppp->size(),
217 ppp->start()
218 );
219 patchArea = SubField<vector>
220 (
221 syncedBAreas,
222 ppp->size(),
223 ppp->offset()
224 );
225 }
226 }
227 }
228
229
230 // Force recalculating cell properties
231 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
232
233 pointField cellCentres(mesh_.nCells());
234 scalarField cellVolumes(mesh_.nCells());
236 (
237 mesh_,
238 faceCentres,
239 faceAreas,
240 cellCentres,
241 cellVolumes
242 );
243
244
245 label nFaces = 0;
246 forAll(faceCentres, facei)
247 {
248 const auto& oldFc = mesh_.faceCentres()[facei];
249 const auto& newFc = faceCentres[facei];
250 const auto& oldFa = mesh_.faceAreas()[facei];
251 const auto& newFa = faceAreas[facei];
252
253 if (oldFc != newFc || oldFa != newFa)
254 {
255 nFaces++;
256
257 if (mesh_.isInternalFace(facei))
258 {
260 << "Different geometry for internal face:" << facei
261 << " oldFc:" << oldFc << nl
262 << " newFc:" << newFc << nl
263 << " oldFa:" << oldFa << nl
264 << " newFa:" << newFa << nl
265 << endl;
266 }
267 }
268 }
269
270 label nCells = 0;
271 forAll(cellCentres, celli)
272 {
273 const auto& oldCc = mesh_.cellCentres()[celli];
274 const auto& newCc = cellCentres[celli];
275 const auto& oldCv = mesh_.cellVolumes()[celli];
276 const auto& newCv = cellVolumes[celli];
277
278 if (oldCc != newCc || oldCv != newCv)
279 {
280 nCells++;
281 //Pout<< "cell:" << celli << nl
282 // << " oldCc:" << oldCc << nl
283 // << " newCc:" << newCc << nl
284 // << " oldCv:" << oldCv << nl
285 // << " newCv:" << newCv << nl
286 // << endl;
287 }
288 }
289
290 Info<< "parallelFvGeometryScheme::movePoints() : "
291 << "adjusted geometry of faces:"
292 << returnReduce(nFaces, sumOp<label>())
293 << " of cells:"
294 << returnReduce(nCells, sumOp<label>())
295 << endl;
296
297 const_cast<fvMesh&>(mesh_).primitiveMesh::resetGeometry
298 (
299 std::move(faceCentres),
300 std::move(faceAreas),
301 std::move(cellCentres),
302 std::move(cellVolumes)
303 );
304}
305
306
307// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
308
309Foam::parallelFvGeometryScheme::parallelFvGeometryScheme
310(
311 const fvMesh& mesh,
312 const dictionary& dict
313)
314:
315 fvGeometryScheme(mesh, dict),
316 dict_(dict.subOrEmptyDict("geometry"))
317{}
318
319
320// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
321
323{
324 if (!geometryPtr_)
325 {
326 if (debug)
327 {
328 Pout<< "parallelFvGeometryScheme::geometry() : "
329 << "constructing underlying scheme from " << dict_
330 << endl;
331 }
332 geometryPtr_ = fvGeometryScheme::New
333 (
334 mesh_,
335 dict_,
336 basicFvGeometryScheme::typeName
337 );
338 }
339 return geometryPtr_();
340}
341
342
344{
345 if (debug)
346 {
347 Pout<< "parallelFvGeometryScheme::movePoints() : "
348 << "recalculating primitiveMesh centres" << endl;
349 }
350
351 // Do any primitive geometry calculation
352 const_cast<fvGeometryScheme&>(geometry()).movePoints();
353
354 // Do in-place overriding processor boundary geometry
355 adjustGeometry();
356}
357
358
360{
361 return const_cast<fvGeometryScheme&>(geometry()).updateMesh(mpm);
362}
363
364
367{
368 return geometry().weights();
369}
370
371
374{
375 return geometry().deltaCoeffs();
376}
377
378
381{
382 return geometry().nonOrthDeltaCoeffs();
383}
384
385
388{
389 return geometry().nonOrthCorrectionVectors();
390}
391
392
394(
395 const pointField& points,
396 const refPtr<pointField>& oldPoints,
397 pointField& faceCentres,
398 vectorField& faceAreas,
399 pointField& cellCentres,
400 scalarField& cellVolumes
401) const
402{
404 (
405 mesh_,
406 points,
407 faceCentres,
408 cellCentres
409 );
410
411 // Parallel consistency
412 adjustGeometry(faceCentres, faceAreas);
413
415 (
416 mesh_,
417 faceCentres,
418 faceAreas,
419 cellCentres,
420 cellVolumes
421 );
422
423 return true;
424}
425
426
427// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition SubField.H:58
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base class for geometry calculation schemes.
const fvMesh & mesh_
Hold reference to mesh.
const fvMesh & mesh() const
Return mesh reference.
static tmp< fvGeometryScheme > New(const fvMesh &mesh, const dictionary &dict, const word &defaultScheme)
Return new tmp interpolation scheme.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Geometry calculation scheme with explicit sync of face-geometry across processor patches.
tmp< fvGeometryScheme > geometryPtr_
Demand-driven construction of underlying scheme.
virtual bool updateGeom(const pointField &points, const refPtr< pointField > &oldPoints, pointField &faceCentres, vectorField &faceAreas, pointField &cellCentres, scalarField &cellVolumes) const
Calculate geometry quantities using mesh topology and provided points. If oldPoints provided only doe...
virtual tmp< surfaceScalarField > nonOrthDeltaCoeffs() const
Return non-orthogonal cell-centre difference coefficients.
const fvGeometryScheme & geometry() const
Construct underlying fvGeometryScheme.
virtual tmp< surfaceVectorField > nonOrthCorrectionVectors() const
Return non-orthogonality correction vectors.
dictionary dict_
Dictionary for underlying scheme.
virtual void movePoints()
Do what is necessary if the mesh has moved.
virtual tmp< surfaceScalarField > weights() const
Return linear difference weighting factors.
virtual tmp< surfaceScalarField > deltaCoeffs() const
Return cell-centre difference coefficients.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh for topology changes.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
static void makeFaceCentresAndAreas(const UList< face > &faces, const pointField &p, vectorField &fCtrs, vectorField &fAreas)
Calculate face centres and areas for specified faces.
static void makeCellCentresAndVols(const primitiveMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, vectorField &cellCtrs, scalarField &cellVols)
Calculate cell centres and volumes from face properties.
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces).
label nInternalFaces() const noexcept
Number of internal faces.
void resetGeometry(pointField &&faceCentres, pointField &&faceAreas, pointField &&cellCentres, scalarField &&cellVolumes)
Reset the local geometry.
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
static void syncBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const TransformOp &top, const bool parRun=UPstream::parRun())
Synchronize values on boundary faces only.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled boundary face values. Uses eqOp.
Definition syncTools.H:524
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
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
Field< vector > vectorField
Specialisation of Field<T> for vector.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.