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"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39}
40
41
42// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43
45(
46 const fvMesh& mesh,
47 const dictionary& dict
48)
49:
51{}
52
53
54// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55
57{
59
60 if (debug)
61 {
62 Pout<< "basicFvGeometryScheme::movePoints() : "
63 << "recalculating primitiveMesh centres" << endl;
64 }
65
66 // Use lower level to calculate the geometry
67 const_cast<fvMesh&>(mesh_).primitiveMesh::updateGeom();
68}
69
70
72{
73 if (debug)
74 {
75 Pout<< "basicFvGeometryScheme::weights() : "
76 << "Constructing weighting factors for face interpolation"
77 << endl;
78 }
79
80 auto tweights =
82 (
84 (
85 "weights",
86 mesh_.pointsInstance(),
87 mesh_,
90 false // Do not register
91 ),
92 mesh_,
94 );
95
96 auto& weights = tweights.ref();
97 weights.setOriented();
98
99 // Set local references to mesh data
100 // Note that we should not use fvMesh sliced fields at this point yet
101 // since this causes a loop when generating weighting factors in
102 // coupledFvPatchField evaluation phase
103 const labelUList& owner = mesh_.owner();
104 const labelUList& neighbour = mesh_.neighbour();
105
106 const vectorField& Cf = mesh_.faceCentres();
107 const vectorField& C = mesh_.cellCentres();
108 const vectorField& Sf = mesh_.faceAreas();
109
110 // ... and reference to the internal field of the weighting factors
111 scalarField& w = weights.primitiveFieldRef();
112
113 forAll(owner, facei)
114 {
115 // Note: mag in the dot-product.
116 // For all valid meshes, the non-orthogonality will be less than
117 // 90 deg and the dot-product will be positive. For invalid
118 // meshes (d & s <= 0), this will stabilise the calculation
119 // but the result will be poor.
120 scalar SfdOwn = mag(Sf[facei] & (Cf[facei] - C[owner[facei]]));
121 scalar SfdNei = mag(Sf[facei] & (C[neighbour[facei]] - Cf[facei]));
122
123 if (mag(SfdOwn + SfdNei) > ROOTVSMALL)
124 {
125 w[facei] = SfdNei/(SfdOwn + SfdNei);
126 }
127 else
128 {
129 w[facei] = 0.5;
130 }
131 }
132
133 auto& wBf = weights.boundaryFieldRef();
134
135 forAll(mesh_.boundary(), patchi)
136 {
137 mesh_.boundary()[patchi].makeWeights(wBf[patchi]);
138 }
139
140 if (debug)
141 {
142 Pout<< "basicFvGeometryScheme::weights : "
143 << "Finished constructing weighting factors for face interpolation"
144 << endl;
145 }
146 return tweights;
147}
148
149
152{
153 if (debug)
154 {
155 Pout<< "basicFvGeometryScheme::deltaCoeffs() : "
156 << "Constructing differencing factors array for face gradient"
157 << endl;
158 }
159
160 // Force the construction of the weighting factors
161 // needed to make sure deltaCoeffs are calculated for parallel runs.
162 (void)mesh_.weights();
163
164 auto tdeltaCoeffs =
166 (
168 (
169 "deltaCoeffs",
170 mesh_.pointsInstance(),
171 mesh_,
174 false // Do not register
175 ),
176 mesh_,
178 );
179
180 auto& deltaCoeffs = tdeltaCoeffs.ref();
181 deltaCoeffs.setOriented();
182
183
184 // Set local references to mesh data
185 const vectorField& C = mesh_.cellCentres();
186 const labelUList& owner = mesh_.owner();
187 const labelUList& neighbour = mesh_.neighbour();
188
189 forAll(owner, facei)
190 {
191 deltaCoeffs[facei] = 1.0/mag(C[neighbour[facei]] - C[owner[facei]]);
192 }
193
194 auto& deltaCoeffsBf = deltaCoeffs.boundaryFieldRef();
195
196 forAll(deltaCoeffsBf, patchi)
197 {
198 const fvPatch& p = mesh_.boundary()[patchi];
199 deltaCoeffsBf[patchi] = 1.0/mag(p.delta());
200
201 // Optionally correct
202 p.makeDeltaCoeffs(deltaCoeffsBf[patchi]);
203 }
204
205 return tdeltaCoeffs;
206}
207
208
211{
212 if (debug)
213 {
214 Pout<< "basicFvGeometryScheme::nonOrthDeltaCoeffs() : "
215 << "Constructing differencing factors array for face gradient"
216 << endl;
217 }
218
219 // Force the construction of the weighting factors
220 // needed to make sure deltaCoeffs are calculated for parallel runs.
221 weights();
222
223 auto tnonOrthDeltaCoeffs =
225 (
227 (
228 "nonOrthDeltaCoeffs",
229 mesh_.pointsInstance(),
230 mesh_,
233 false // Do not register
234 ),
235 mesh_,
237 );
238
239 auto& nonOrthDeltaCoeffs = tnonOrthDeltaCoeffs.ref();
240 nonOrthDeltaCoeffs.setOriented();
241
242
243 // Set local references to mesh data
244 const volVectorField& C = mesh_.C();
245 const labelUList& owner = mesh_.owner();
246 const labelUList& neighbour = mesh_.neighbour();
247 const surfaceVectorField& Sf = mesh_.Sf();
248 const surfaceScalarField& magSf = mesh_.magSf();
249
250 forAll(owner, facei)
251 {
252 vector delta = C[neighbour[facei]] - C[owner[facei]];
253 vector unitArea = Sf[facei]/magSf[facei];
254
255 // Standard cell-centre distance form
256 //NonOrthDeltaCoeffs[facei] = (unitArea & delta)/magSqr(delta);
257
258 // Slightly under-relaxed form
259 //NonOrthDeltaCoeffs[facei] = 1.0/mag(delta);
260
261 // More under-relaxed form
262 //NonOrthDeltaCoeffs[facei] = 1.0/(mag(unitArea & delta) + VSMALL);
263
264 // Stabilised form for bad meshes
265 nonOrthDeltaCoeffs[facei] = 1.0/max(unitArea & delta, 0.05*mag(delta));
266 }
267
268 auto& nonOrthDeltaCoeffsBf = nonOrthDeltaCoeffs.boundaryFieldRef();
269
270 forAll(nonOrthDeltaCoeffsBf, patchi)
271 {
272 fvsPatchScalarField& patchDeltaCoeffs = nonOrthDeltaCoeffsBf[patchi];
273
274 const fvPatch& p = patchDeltaCoeffs.patch();
275
276 const vectorField patchDeltas(mesh_.boundary()[patchi].delta());
277
278 forAll(p, patchFacei)
279 {
280 vector unitArea =
281 Sf.boundaryField()[patchi][patchFacei]
282 /magSf.boundaryField()[patchi][patchFacei];
283
284 const vector& delta = patchDeltas[patchFacei];
285
286 patchDeltaCoeffs[patchFacei] =
287 1.0/max(unitArea & delta, 0.05*mag(delta));
288 }
289
290 // Optionally correct
291 p.makeNonOrthoDeltaCoeffs(patchDeltaCoeffs);
292 }
293 return tnonOrthDeltaCoeffs;
294}
295
296
299{
300 if (debug)
301 {
302 Pout<< "surfaceInterpolation::makeNonOrthCorrectionVectors() : "
303 << "Constructing non-orthogonal correction vectors"
304 << endl;
305 }
306
307 auto tnonOrthCorrectionVectors =
309 (
311 (
312 "nonOrthCorrectionVectors",
313 mesh_.pointsInstance(),
314 mesh_,
317 false // Do not register
318 ),
319 mesh_,
320 dimless
321 );
322
323 auto& corrVecs = tnonOrthCorrectionVectors.ref();
324 corrVecs.setOriented();
325
326 // Set local references to mesh data
327 const volVectorField& C = mesh_.C();
328 const labelUList& owner = mesh_.owner();
329 const labelUList& neighbour = mesh_.neighbour();
330 const surfaceVectorField& Sf = mesh_.Sf();
331 const surfaceScalarField& magSf = mesh_.magSf();
332 tmp<surfaceScalarField> tNonOrthDeltaCoeffs(nonOrthDeltaCoeffs());
333 const surfaceScalarField& NonOrthDeltaCoeffs = tNonOrthDeltaCoeffs();
334
335 forAll(owner, facei)
336 {
337 vector unitArea(Sf[facei]/magSf[facei]);
338 vector delta(C[neighbour[facei]] - C[owner[facei]]);
339
340 corrVecs[facei] = unitArea - delta*NonOrthDeltaCoeffs[facei];
341 }
342
343 // Boundary correction vectors set to zero for boundary patches
344 // and calculated consistently with internal corrections for
345 // coupled patches
346
347 auto& corrVecsBf = corrVecs.boundaryFieldRef();
348
349 forAll(corrVecsBf, patchi)
350 {
351 fvsPatchVectorField& patchCorrVecs = corrVecsBf[patchi];
352
353 const fvPatch& p = patchCorrVecs.patch();
354
355 if (!patchCorrVecs.coupled())
356 {
357 patchCorrVecs = Zero;
358 }
359 else
360 {
361 const auto& patchNonOrthDeltaCoeffs =
362 NonOrthDeltaCoeffs.boundaryField()[patchi];
363
364 const vectorField patchDeltas(mesh_.boundary()[patchi].delta());
365
366 forAll(p, patchFacei)
367 {
368 vector unitArea =
369 Sf.boundaryField()[patchi][patchFacei]
370 /magSf.boundaryField()[patchi][patchFacei];
371
372 const vector& delta = patchDeltas[patchFacei];
373
374 patchCorrVecs[patchFacei] =
375 unitArea - delta*patchNonOrthDeltaCoeffs[patchFacei];
376 }
377 }
378
379 // Optionally correct
380 p.makeNonOrthoCorrVectors(patchCorrVecs);
381 }
382
383 if (debug)
384 {
385 Pout<< "surfaceInterpolation::makeNonOrthCorrectionVectors() : "
386 << "Finished constructing non-orthogonal correction vectors"
387 << endl;
388 }
389 return tnonOrthCorrectionVectors;
390}
391
392
393// ************************************************************************* //
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:53
C()
Construct null.
Definition: C.C:43
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Default geometry calculation scheme. Slight stabilisation for bad meshes.
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:126
Abstract base class for geometry calculation schemes.
virtual void movePoints()
Update basic geometric properties from provided points.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual bool coupled() const
Return true if this patch field is coupled.
const fvPatch & patch() const
Return patch.
virtual void updateGeom()
Update all geometric data.
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
dynamicFvMesh & mesh
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.