Loading...
Searching...
No Matches
basicFvGeometryScheme.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) 2020-2021 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 "surfaceFields.H"
31#include "volFields.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38 defineTypeNameAndDebug(basicFvGeometryScheme, 0);
39 addToRunTimeSelectionTable(fvGeometryScheme, basicFvGeometryScheme, dict);
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45Foam::basicFvGeometryScheme::basicFvGeometryScheme
46(
47 const fvMesh& mesh,
48 const dictionary& dict
49)
51 fvGeometryScheme(mesh, dict)
52{}
53
54
55// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56
58{
60
61 if (debug)
62 {
63 Pout<< "basicFvGeometryScheme::movePoints() : "
64 << "recalculating primitiveMesh centres" << endl;
65 }
66
67 // Use lower level to calculate the geometry
69}
70
71
73{
74 if (debug)
75 {
76 Pout<< "basicFvGeometryScheme::weights() : "
77 << "Constructing weighting factors for face interpolation"
78 << endl;
79 }
80
81 auto tweights =
83 (
85 (
86 "weights",
87 mesh_.pointsInstance(),
88 mesh_,
92 ),
93 mesh_,
95 );
96
97 auto& weights = tweights.ref();
98 weights.setOriented();
99
100 // Set local references to mesh data
101 // Note that we should not use fvMesh sliced fields at this point yet
102 // since this causes a loop when generating weighting factors in
103 // coupledFvPatchField evaluation phase
104 const labelUList& owner = mesh_.owner();
105 const labelUList& neighbour = mesh_.neighbour();
106
107 const vectorField& Cf = mesh_.faceCentres();
108 const vectorField& C = mesh_.cellCentres();
109 const vectorField& Sf = mesh_.faceAreas();
110
111 // ... and reference to the internal field of the weighting factors
112 scalarField& w = weights.primitiveFieldRef();
113
114 forAll(owner, facei)
115 {
116 // Note: mag in the dot-product.
117 // For all valid meshes, the non-orthogonality will be less than
118 // 90 deg and the dot-product will be positive. For invalid
119 // meshes (d & s <= 0), this will stabilise the calculation
120 // but the result will be poor.
121 scalar SfdOwn = mag(Sf[facei] & (Cf[facei] - C[owner[facei]]));
122 scalar SfdNei = mag(Sf[facei] & (C[neighbour[facei]] - Cf[facei]));
123
124 if (mag(SfdOwn + SfdNei) > ROOTVSMALL)
125 {
126 w[facei] = SfdNei/(SfdOwn + SfdNei);
127 }
128 else
129 {
130 w[facei] = 0.5;
131 }
132 }
133
134 auto& wBf = weights.boundaryFieldRef();
135
136 forAll(mesh_.boundary(), patchi)
137 {
138 mesh_.boundary()[patchi].makeWeights(wBf[patchi]);
139 }
140
141 if (debug)
142 {
143 Pout<< "basicFvGeometryScheme::weights : "
144 << "Finished constructing weighting factors for face interpolation"
146 }
147 return tweights;
148}
149
150
153{
154 if (debug)
155 {
156 Pout<< "basicFvGeometryScheme::deltaCoeffs() : "
157 << "Constructing differencing factors array for face gradient"
158 << endl;
159 }
160
161 // Force the construction of the weighting factors
162 // needed to make sure deltaCoeffs are calculated for parallel runs.
163 (void)mesh_.weights();
164
165 auto tdeltaCoeffs =
167 (
169 (
170 "deltaCoeffs",
171 mesh_.pointsInstance(),
172 mesh_,
176 ),
177 mesh_,
179 );
180
181 auto& deltaCoeffs = tdeltaCoeffs.ref();
182 deltaCoeffs.setOriented();
183
184
185 // Set local references to mesh data
186 const vectorField& C = mesh_.cellCentres();
187 const labelUList& owner = mesh_.owner();
188 const labelUList& neighbour = mesh_.neighbour();
189
190 forAll(owner, facei)
191 {
192 deltaCoeffs[facei] = 1.0/mag(C[neighbour[facei]] - C[owner[facei]]);
193 }
194
195 auto& deltaCoeffsBf = deltaCoeffs.boundaryFieldRef();
196
197 forAll(deltaCoeffsBf, patchi)
198 {
199 const fvPatch& p = mesh_.boundary()[patchi];
200 deltaCoeffsBf[patchi] = 1.0/mag(p.delta());
201
202 // Optionally correct
203 p.makeDeltaCoeffs(deltaCoeffsBf[patchi]);
205
206 return tdeltaCoeffs;
207}
208
209
212{
213 if (debug)
214 {
215 Pout<< "basicFvGeometryScheme::nonOrthDeltaCoeffs() : "
216 << "Constructing differencing factors array for face gradient"
217 << endl;
218 }
219
220 // Force the construction of the weighting factors
221 // needed to make sure deltaCoeffs are calculated for parallel runs.
222 weights();
223
224 auto tnonOrthDeltaCoeffs =
226 (
228 (
229 "nonOrthDeltaCoeffs",
230 mesh_.pointsInstance(),
231 mesh_,
235 ),
236 mesh_,
238 );
239
240 auto& nonOrthDeltaCoeffs = tnonOrthDeltaCoeffs.ref();
241 nonOrthDeltaCoeffs.setOriented();
242
243
244 // Set local references to mesh data
245 const volVectorField& C = mesh_.C();
246 const labelUList& owner = mesh_.owner();
247 const labelUList& neighbour = mesh_.neighbour();
248 const surfaceVectorField& Sf = mesh_.Sf();
249 const surfaceScalarField& magSf = mesh_.magSf();
250
251 forAll(owner, facei)
252 {
253 vector delta = C[neighbour[facei]] - C[owner[facei]];
254 vector unitArea = Sf[facei]/magSf[facei];
255
256 // Standard cell-centre distance form
257 //NonOrthDeltaCoeffs[facei] = (unitArea & delta)/magSqr(delta);
258
259 // Slightly under-relaxed form
260 //NonOrthDeltaCoeffs[facei] = 1.0/mag(delta);
261
262 // More under-relaxed form
263 //NonOrthDeltaCoeffs[facei] = 1.0/(mag(unitArea & delta) + VSMALL);
264
265 // Stabilised form for bad meshes
266 nonOrthDeltaCoeffs[facei] = 1.0/max(unitArea & delta, 0.05*mag(delta));
267 }
268
269 auto& nonOrthDeltaCoeffsBf = nonOrthDeltaCoeffs.boundaryFieldRef();
270
271 forAll(nonOrthDeltaCoeffsBf, patchi)
272 {
273 fvsPatchScalarField& patchDeltaCoeffs = nonOrthDeltaCoeffsBf[patchi];
274
275 const fvPatch& p = patchDeltaCoeffs.patch();
276
277 const vectorField patchDeltas(mesh_.boundary()[patchi].delta());
278
279 forAll(p, patchFacei)
280 {
281 vector unitArea =
282 Sf.boundaryField()[patchi][patchFacei]
283 /magSf.boundaryField()[patchi][patchFacei];
284
285 const vector& delta = patchDeltas[patchFacei];
286
287 patchDeltaCoeffs[patchFacei] =
288 1.0/max(unitArea & delta, 0.05*mag(delta));
289 }
290
291 // Optionally correct
292 p.makeNonOrthoDeltaCoeffs(patchDeltaCoeffs);
293 }
294 return tnonOrthDeltaCoeffs;
295}
296
297
300{
301 if (debug)
302 {
303 Pout<< "surfaceInterpolation::makeNonOrthCorrectionVectors() : "
304 << "Constructing non-orthogonal correction vectors"
305 << endl;
306 }
307
308 auto tnonOrthCorrectionVectors =
310 (
312 (
313 "nonOrthCorrectionVectors",
314 mesh_.pointsInstance(),
315 mesh_,
319 ),
320 mesh_,
321 dimless
322 );
323
324 auto& corrVecs = tnonOrthCorrectionVectors.ref();
325 corrVecs.setOriented();
326
327 // Set local references to mesh data
328 const volVectorField& C = mesh_.C();
329 const labelUList& owner = mesh_.owner();
330 const labelUList& neighbour = mesh_.neighbour();
331 const surfaceVectorField& Sf = mesh_.Sf();
332 const surfaceScalarField& magSf = mesh_.magSf();
333 tmp<surfaceScalarField> tNonOrthDeltaCoeffs(nonOrthDeltaCoeffs());
334 const surfaceScalarField& NonOrthDeltaCoeffs = tNonOrthDeltaCoeffs();
335
336 forAll(owner, facei)
337 {
338 vector unitArea(Sf[facei]/magSf[facei]);
339 vector delta(C[neighbour[facei]] - C[owner[facei]]);
340
341 corrVecs[facei] = unitArea - delta*NonOrthDeltaCoeffs[facei];
342 }
343
344 // Boundary correction vectors set to zero for boundary patches
345 // and calculated consistently with internal corrections for
346 // coupled patches
347
348 auto& corrVecsBf = corrVecs.boundaryFieldRef();
349
350 forAll(corrVecsBf, patchi)
351 {
352 fvsPatchVectorField& patchCorrVecs = corrVecsBf[patchi];
353
354 const fvPatch& p = patchCorrVecs.patch();
355
356 if (!patchCorrVecs.coupled())
357 {
358 patchCorrVecs = Zero;
359 }
360 else
361 {
362 const auto& patchNonOrthDeltaCoeffs =
363 NonOrthDeltaCoeffs.boundaryField()[patchi];
364
365 const vectorField patchDeltas(mesh_.boundary()[patchi].delta());
366
367 forAll(p, patchFacei)
368 {
369 vector unitArea =
370 Sf.boundaryField()[patchi][patchFacei]
371 /magSf.boundaryField()[patchi][patchFacei];
372
373 const vector& delta = patchDeltas[patchFacei];
374
375 patchCorrVecs[patchFacei] =
376 unitArea - delta*patchNonOrthDeltaCoeffs[patchFacei];
377 }
378 }
379
380 // Optionally correct
381 p.makeNonOrthoCorrVectors(patchCorrVecs);
382 }
383
384 if (debug)
385 {
386 Pout<< "surfaceInterpolation::makeNonOrthCorrectionVectors() : "
387 << "Finished constructing non-orthogonal correction vectors"
388 << endl;
389 }
390 return tnonOrthCorrectionVectors;
391}
392
393
395(
396 const pointField& points,
397 const refPtr<pointField>& oldPoints,
398 pointField& faceCentres,
399 vectorField& faceAreas,
400 pointField& cellCentres,
401 scalarField& cellVolumes
402) const
403{
405 (
406 mesh_,
407 points,
408 faceCentres,
409 faceAreas
410 );
411
413 (
414 mesh_,
415 faceCentres,
416 faceAreas,
417 cellCentres,
418 cellVolumes
419 );
420
421 // Assume something has changed
422 return true;
423}
424
425
426// ************************************************************************* //
scalar delta
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Graphite solid properties.
Definition C.H:49
C() noexcept
Default construct.
Definition C.C:36
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Default geometry calculation scheme. Slight stabilisation for bad meshes.
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.
virtual tmp< surfaceVectorField > nonOrthCorrectionVectors() const
Return non-orthogonality correction vectors.
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.
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.
virtual void movePoints()
Update basic geometric properties from provided points.
const fvMesh & mesh() const
Return mesh reference.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual bool coupled() const
True if the patch field is coupled.
const fvPatch & patch() const noexcept
Return the patch.
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.
virtual void updateGeom()
Update all geometric data.
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
const pointField & points
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvsPatchField< vector > fvsPatchVectorField
fvsPatchField< scalar > fvsPatchScalarField
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.