surfaceInterpolation.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25  Cell to face interpolation scheme. Included in fvMesh.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fvMesh.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "demandDrivenData.H"
33 #include "coupledFvPatch.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
40 }
41 
42 
43 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
44 
46 {
47  deleteDemandDrivenData(weights_);
48  deleteDemandDrivenData(deltaCoeffs_);
49  deleteDemandDrivenData(nonOrthDeltaCoeffs_);
50  deleteDemandDrivenData(nonOrthCorrectionVectors_);
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 :
58  mesh_(fvm),
59  weights_(nullptr),
60  deltaCoeffs_(nullptr),
61  nonOrthDeltaCoeffs_(nullptr),
62  nonOrthCorrectionVectors_(nullptr)
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {
70  clearOut();
71 }
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
78 {
79  if (!weights_)
80  {
81  makeWeights();
82  }
83 
84  return (*weights_);
85 }
86 
87 
90 {
91  if (!deltaCoeffs_)
92  {
93  makeDeltaCoeffs();
94  }
95 
96  return (*deltaCoeffs_);
97 }
98 
99 
102 {
103  if (!nonOrthDeltaCoeffs_)
104  {
105  makeNonOrthDeltaCoeffs();
106  }
107 
108  return (*nonOrthDeltaCoeffs_);
109 }
110 
111 
114 {
115  if (!nonOrthCorrectionVectors_)
116  {
117  makeNonOrthCorrectionVectors();
118  }
119 
120  return (*nonOrthCorrectionVectors_);
121 }
122 
123 
125 {
126  deleteDemandDrivenData(weights_);
127  deleteDemandDrivenData(deltaCoeffs_);
128  deleteDemandDrivenData(nonOrthDeltaCoeffs_);
129  deleteDemandDrivenData(nonOrthCorrectionVectors_);
130 
131  return true;
132 }
133 
134 
135 void Foam::surfaceInterpolation::makeWeights() const
136 {
137  if (debug)
138  {
139  Pout<< "surfaceInterpolation::makeWeights() : "
140  << "Constructing weighting factors for face interpolation"
141  << endl;
142  }
143 
144  weights_ = new surfaceScalarField
145  (
146  IOobject
147  (
148  "weights",
149  mesh_.pointsInstance(),
150  mesh_,
153  false // Do not register
154  ),
155  mesh_,
156  dimless
157  );
158  surfaceScalarField& weights = *weights_;
159 
160  // Set local references to mesh data
161  const labelUList& owner = mesh_.owner();
162  const labelUList& neighbour = mesh_.neighbour();
163  const surfaceVectorField& Sf = mesh_.Sf();
164  const surfaceVectorField& Cf = mesh_.Cf();
165  const volVectorField& C = mesh_.C();
166 
167  // ... and reference to the internal field of the weighting factors
168  scalarField& w = weights.primitiveFieldRef();
169 
170  forAll(owner, facei)
171  {
172  // Note: mag in the dot-product.
173  // For all valid meshes, the non-orthogonality will be less that
174  // 90 deg and the dot-product will be positive. For invalid
175  // meshes (d & s <= 0), this will stabilise the calculation
176  // but the result will be poor.
177  const scalar SfdOwn = mag(Sf[facei]&(Cf[facei] - C[owner[facei]]));
178  const scalar SfdNei = mag(Sf[facei]&(C[neighbour[facei]] - Cf[facei]));
179  const scalar SfdOwnNei = SfdOwn + SfdNei;
180 
181  if (SfdNei/vGreat < SfdOwnNei)
182  {
183  w[facei] = SfdNei/SfdOwnNei;
184  }
185  else
186  {
187  const scalar dOwn = mag(Cf[facei] - C[owner[facei]]);
188  const scalar dNei = mag(C[neighbour[facei]] - Cf[facei]);
189  const scalar dOwnNei = dOwn + dNei;
190 
191  w[facei] = dNei/dOwnNei;
192  }
193  }
194 
196  weights.boundaryFieldRef();
197 
198  forAll(mesh_.boundary(), patchi)
199  {
200  mesh_.boundary()[patchi].makeWeights(wBf[patchi]);
201  }
202 
203  if (debug)
204  {
205  Pout<< "surfaceInterpolation::makeWeights() : "
206  << "Finished constructing weighting factors for face interpolation"
207  << endl;
208  }
209 }
210 
211 
212 void Foam::surfaceInterpolation::makeDeltaCoeffs() const
213 {
214  if (debug)
215  {
216  Pout<< "surfaceInterpolation::makeDeltaCoeffs() : "
217  << "Constructing interpolation factors array for face gradient"
218  << endl;
219  }
220 
221  // Force the construction of the weighting factors
222  // needed to make sure deltaCoeffs are calculated for parallel runs.
223  weights();
224 
225  deltaCoeffs_ = new surfaceScalarField
226  (
227  IOobject
228  (
229  "deltaCoeffs",
230  mesh_.pointsInstance(),
231  mesh_,
234  false // Do not register
235  ),
236  mesh_,
238  );
239  surfaceScalarField& deltaCoeffs = *deltaCoeffs_;
240 
241 
242  // Set local references to mesh data
243  const volVectorField& C = mesh_.C();
244  const labelUList& owner = mesh_.owner();
245  const labelUList& neighbour = mesh_.neighbour();
246 
247  forAll(owner, facei)
248  {
249  deltaCoeffs[facei] = 1.0/mag(C[neighbour[facei]] - C[owner[facei]]);
250  }
251 
252  surfaceScalarField::Boundary& deltaCoeffsBf =
253  deltaCoeffs.boundaryFieldRef();
254 
255  forAll(deltaCoeffsBf, patchi)
256  {
257  deltaCoeffsBf[patchi] = 1.0/mag(mesh_.boundary()[patchi].delta());
258  }
259 }
260 
261 
262 void Foam::surfaceInterpolation::makeNonOrthDeltaCoeffs() const
263 {
264  if (debug)
265  {
266  Pout<< "surfaceInterpolation::makeNonOrthDeltaCoeffs() : "
267  << "Constructing interpolation factors array for face gradient"
268  << endl;
269  }
270 
271  // Force the construction of the weighting factors
272  // needed to make sure deltaCoeffs are calculated for parallel runs.
273  weights();
274 
275  nonOrthDeltaCoeffs_ = new surfaceScalarField
276  (
277  IOobject
278  (
279  "nonOrthDeltaCoeffs",
280  mesh_.pointsInstance(),
281  mesh_,
284  false // Do not register
285  ),
286  mesh_,
288  );
289  surfaceScalarField& nonOrthDeltaCoeffs = *nonOrthDeltaCoeffs_;
290 
291 
292  // Set local references to mesh data
293  const volVectorField& C = mesh_.C();
294  const labelUList& owner = mesh_.owner();
295  const labelUList& neighbour = mesh_.neighbour();
296  const surfaceVectorField& Sf = mesh_.Sf();
297  const surfaceScalarField& magSf = mesh_.magSf();
298 
299  forAll(owner, facei)
300  {
301  vector delta = C[neighbour[facei]] - C[owner[facei]];
302  vector unitArea = Sf[facei]/magSf[facei];
303 
304  // Standard cell-centre distance form
305  // NonOrthDeltaCoeffs[facei] = (unitArea & delta)/magSqr(delta);
306 
307  // Slightly under-relaxed form
308  // NonOrthDeltaCoeffs[facei] = 1.0/mag(delta);
309 
310  // More under-relaxed form
311  // NonOrthDeltaCoeffs[facei] = 1.0/(mag(unitArea & delta) + vSmall);
312 
313  // Stabilised form for bad meshes
314  nonOrthDeltaCoeffs[facei] = 1.0/max(unitArea & delta, 0.05*mag(delta));
315  }
316 
317  surfaceScalarField::Boundary& nonOrthDeltaCoeffsBf =
318  nonOrthDeltaCoeffs.boundaryFieldRef();
319 
320  forAll(nonOrthDeltaCoeffsBf, patchi)
321  {
322  vectorField delta(mesh_.boundary()[patchi].delta());
323 
324  nonOrthDeltaCoeffsBf[patchi] =
325  1.0/max(mesh_.boundary()[patchi].nf() & delta, 0.05*mag(delta));
326  }
327 }
328 
329 
330 void Foam::surfaceInterpolation::makeNonOrthCorrectionVectors() const
331 {
332  if (debug)
333  {
334  Pout<< "surfaceInterpolation::makeNonOrthCorrectionVectors() : "
335  << "Constructing non-orthogonal correction vectors"
336  << endl;
337  }
338 
339  nonOrthCorrectionVectors_ = new surfaceVectorField
340  (
341  IOobject
342  (
343  "nonOrthCorrectionVectors",
344  mesh_.pointsInstance(),
345  mesh_,
348  false // Do not register
349  ),
350  mesh_,
351  dimless
352  );
353  surfaceVectorField& corrVecs = *nonOrthCorrectionVectors_;
354 
355  // Set local references to mesh data
356  const volVectorField& C = mesh_.C();
357  const labelUList& owner = mesh_.owner();
358  const labelUList& neighbour = mesh_.neighbour();
359  const surfaceVectorField& Sf = mesh_.Sf();
360  const surfaceScalarField& magSf = mesh_.magSf();
361  const surfaceScalarField& NonOrthDeltaCoeffs = nonOrthDeltaCoeffs();
362 
363  forAll(owner, facei)
364  {
365  vector unitArea = Sf[facei]/magSf[facei];
366  vector delta = C[neighbour[facei]] - C[owner[facei]];
367 
368  corrVecs[facei] = unitArea - delta*NonOrthDeltaCoeffs[facei];
369  }
370 
371  // Boundary correction vectors set to zero for boundary patches
372  // and calculated consistently with internal corrections for
373  // coupled patches
374 
375  surfaceVectorField::Boundary& corrVecsBf =
376  corrVecs.boundaryFieldRef();
377 
378  forAll(corrVecsBf, patchi)
379  {
380  fvsPatchVectorField& patchCorrVecs = corrVecsBf[patchi];
381 
382  if (!patchCorrVecs.coupled())
383  {
384  patchCorrVecs = Zero;
385  }
386  else
387  {
388  const fvsPatchScalarField& patchNonOrthDeltaCoeffs
389  = NonOrthDeltaCoeffs.boundaryField()[patchi];
390 
391  const fvPatch& p = patchCorrVecs.patch();
392 
393  const vectorField patchDeltas(mesh_.boundary()[patchi].delta());
394 
395  forAll(p, patchFacei)
396  {
397  vector unitArea =
398  Sf.boundaryField()[patchi][patchFacei]
399  /magSf.boundaryField()[patchi][patchFacei];
400 
401  const vector& delta = patchDeltas[patchFacei];
402 
403  patchCorrVecs[patchFacei] =
404  unitArea - delta*patchNonOrthDeltaCoeffs[patchFacei];
405  }
406  }
407  }
408 
409  if (debug)
410  {
411  Pout<< "surfaceInterpolation::makeNonOrthCorrectionVectors() : "
412  << "Finished constructing non-orthogonal correction vectors"
413  << endl;
414  }
415 }
416 
417 
418 // ************************************************************************* //
static const Foam::dimensionedScalar C("C", Foam::dimTemperature, 234.5)
scalar delta
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Cell to surface interpolation scheme. Included in fvMesh.
bool movePoints()
Do what is necessary if the mesh has moved.
const surfaceVectorField & nonOrthCorrectionVectors() const
Return reference to non-orthogonality correction vectors.
const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
surfaceInterpolation(const fvMesh &)
Construct given an fvMesh.
const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
void clearOut()
Clear all geometry and addressing.
Template functions to aid in the implementation of demand driven data.
label patchi
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
VolField< vector > volVectorField
Definition: volFieldsFwd.H:65
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
void deleteDemandDrivenData(DataType *&dataPtr)
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
SurfaceField< vector > surfaceVectorField
UList< label > labelUList
Definition: UList.H:65
fvsPatchField< vector > fvsPatchVectorField
fvsPatchField< scalar > fvsPatchScalarField
volScalarField & p
Foam::surfaceFields.