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