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