fvMeshGeometry.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 \*---------------------------------------------------------------------------*/
25 
26 #include "fvMesh.H"
27 #include "Time.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "slicedVolFields.H"
31 #include "slicedSurfaceFields.H"
32 #include "SubField.H"
33 #include "cyclicFvPatchFields.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 void Foam::fvMesh::makeSf() const
38 {
39  if (debug)
40  {
41  InfoInFunction << "Assembling face areas" << endl;
42  }
43 
44  // It is an error to attempt to recalculate
45  // if the pointer is already set
46  if (SfSlicePtr_ || SfPtr_)
47  {
49  << "face areas already exist"
50  << abort(FatalError);
51  }
52 
53  SfSlicePtr_ = new slicedSurfaceVectorField
54  (
55  IOobject
56  (
57  "Sf",
59  meshSubDir,
60  *this,
63  true
64  ),
65  *this,
66  dimArea,
67  faceAreas()
68  );
69 }
70 
71 
72 void Foam::fvMesh::makeMagSf() const
73 {
74  if (debug)
75  {
76  InfoInFunction << "Assembling mag face areas" << endl;
77  }
78 
79  // It is an error to attempt to recalculate
80  // if the pointer is already set
81  if (magSfSlicePtr_ || magSfPtr_)
82  {
84  << "mag face areas already exist"
85  << abort(FatalError);
86  }
87 
88  magSfSlicePtr_ = new slicedSurfaceScalarField
89  (
90  IOobject
91  (
92  "magSf",
93  pointsInstance(),
94  meshSubDir,
95  *this,
98  true
99  ),
100  *this,
101  dimArea,
102  magFaceAreas()
103  );
104 }
105 
106 
107 void Foam::fvMesh::makeC() const
108 {
109  if (debug)
110  {
111  InfoInFunction << "Assembling cell centres" << endl;
112  }
113 
114  // It is an error to attempt to recalculate
115  // if the pointer is already set
116  if (CSlicePtr_ || CPtr_)
117  {
119  << "cell centres already exist"
120  << abort(FatalError);
121  }
122 
123  // Construct as slices. Only preserve processor (not e.g. cyclic)
124 
125  CSlicePtr_ = new slicedVolVectorField
126  (
127  IOobject
128  (
129  "Cc",
130  pointsInstance(),
131  meshSubDir,
132  *this,
135  true
136  ),
137  *this,
138  dimLength,
139  cellCentres(),
140  faceCentres(),
141  true, // preserveCouples
142  true // preserveProcOnly
143  );
144 }
145 
146 
147 void Foam::fvMesh::makeCf() const
148 {
149  if (debug)
150  {
151  InfoInFunction << "Assembling face centres" << endl;
152  }
153 
154  // It is an error to attempt to recalculate
155  // if the pointer is already set
156  if (CfSlicePtr_ || CfPtr_)
157  {
159  << "face centres already exist"
160  << abort(FatalError);
161  }
162 
163  CfSlicePtr_ = new slicedSurfaceVectorField
164  (
165  IOobject
166  (
167  "Cf",
168  pointsInstance(),
169  meshSubDir,
170  *this,
173  true
174  ),
175  *this,
176  dimLength,
177  faceCentres()
178  );
179 }
180 
181 
182 Foam::surfaceVectorField& Foam::fvMesh::SfRef()
183 {
184  if (!SfPtr_)
185  {
186  SfPtr_ = Sf().cloneUnSliced().ptr();
187 
188  deleteDemandDrivenData(SfSlicePtr_);
189  }
190 
191  return *SfPtr_;
192 }
193 
194 
195 Foam::surfaceScalarField& Foam::fvMesh::magSfRef()
196 {
197  if (!magSfPtr_)
198  {
199  magSfPtr_ = magSf().cloneUnSliced().ptr();
200 
201  deleteDemandDrivenData(magSfSlicePtr_);
202  }
203 
204  return *magSfPtr_;
205 }
206 
207 
208 Foam::volVectorField& Foam::fvMesh::CRef()
209 {
210  if (!CPtr_)
211  {
212  CPtr_ = C().cloneUnSliced().ptr();
213 
214  deleteDemandDrivenData(CSlicePtr_);
215  }
216 
217  return *CPtr_;
218 }
219 
220 
221 Foam::surfaceVectorField& Foam::fvMesh::CfRef()
222 {
223  if (!CfPtr_)
224  {
225  CfPtr_ = Cf().cloneUnSliced().ptr();
226 
227  deleteDemandDrivenData(CfSlicePtr_);
228  }
229 
230  return *CfPtr_;
231 }
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
237 {
238  if (!VPtr_)
239  {
240  if (debug)
241  {
243  << "Constructing from primitiveMesh::cellVolumes()" << endl;
244  }
245 
247  (
248  IOobject
249  (
250  "Vc",
251  time().name(),
252  *this,
255  true
256  ),
257  *this,
258  dimVolume,
259  cellVolumes()
260  );
261  }
262 
263  return *VPtr_;
264 }
265 
266 
268 {
269  if (!V0Ptr_ || Foam::isNull(V0Ptr_))
270  {
272  << "Vc0 is not available"
273  << abort(FatalError);
274  }
275 
276  return *V0Ptr_;
277 }
278 
279 
281 {
282  if (!V00Ptr_)
283  {
284  V00Ptr_ = NullObjectPtr<DimensionedField<scalar, volMesh>>();
285  }
286 
287  if (Foam::isNull(V00Ptr_))
288  {
289  return V0();
290  }
291 
292  return *V00Ptr_;
293 }
294 
295 
297 {
298  if (moving() && time().subCycling())
299  {
300  const TimeState& ts = time();
301  const TimeState& ts0 = time().prevTimeState();
302 
303  scalar tFrac =
304  (
305  ts.value() - (ts0.value() - ts0.deltaTValue())
306  )/ts0.deltaTValue();
307 
308  if (tFrac < (1 - small))
309  {
310  return V0() + tFrac*(V() - V0());
311  }
312  else
313  {
314  return V();
315  }
316  }
317  else
318  {
319  return V();
320  }
321 }
322 
323 
325 {
326  if (moving() && time().subCycling())
327  {
328  const TimeState& ts = time();
329  const TimeState& ts0 = time().prevTimeState();
330 
331  scalar t0Frac =
332  (
333  (ts.value() - ts.deltaTValue())
334  - (ts0.value() - ts0.deltaTValue())
335  )/ts0.deltaTValue();
336 
337  if (t0Frac > small)
338  {
339  return V0() + t0Frac*(V() - V0());
340  }
341  else
342  {
343  return V0();
344  }
345  }
346  else
347  {
348  return V0();
349  }
350 }
351 
352 
354 {
355  if (SfPtr_)
356  {
357  return *SfPtr_;
358  }
359 
360  if (!SfSlicePtr_)
361  {
362  makeSf();
363  }
364 
365  return *SfSlicePtr_;
366 }
367 
368 
370 {
371  if (magSfPtr_)
372  {
373  return *magSfPtr_;
374  }
375 
376  if (!magSfSlicePtr_)
377  {
378  makeMagSf();
379  }
380 
381  return *magSfSlicePtr_;
382 }
383 
384 
386 {
387  if (CPtr_)
388  {
389  return *CPtr_;
390  }
391 
392  if (!CSlicePtr_)
393  {
394  makeC();
395  }
396 
397  return *CSlicePtr_;
398 }
399 
400 
402 {
403  if (CfPtr_)
404  {
405  return *CfPtr_;
406  }
407 
408  if (!CfSlicePtr_)
409  {
410  makeCf();
411  }
412 
413  return *CfSlicePtr_;
414 }
415 
416 
418 {
419  if (debug)
420  {
421  InfoInFunction << "Calculating face deltas" << endl;
422  }
423 
425  (
427  (
428  "delta",
429  *this,
430  dimLength
431  )
432  );
433  surfaceVectorField& delta = tdelta.ref();
434 
435  const volVectorField& C = this->C();
436  const labelUList& owner = this->owner();
437  const labelUList& neighbour = this->neighbour();
438 
439  forAll(owner, facei)
440  {
441  delta[facei] = C[neighbour[facei]] - C[owner[facei]];
442  }
443 
445  delta.boundaryFieldRef();
446 
447  forAll(deltabf, patchi)
448  {
449  deltabf[patchi] = boundary()[patchi].delta();
450  }
451 
452  return tdelta;
453 }
454 
455 
457 {
458  if (!phiPtr_)
459  {
461  << "mesh flux field does not exist, is the mesh actually moving?"
462  << abort(FatalError);
463  }
464 
465  // Set zero current time
466  // mesh motion fluxes if the time has been incremented
467  if (!time().subCycling() && phiPtr_->timeIndex() != time().timeIndex())
468  {
469  (*phiPtr_) = dimensionedScalar(dimVolume/dimTime, 0);
470  }
471 
472  return *phiPtr_;
473 }
474 
475 
476 Foam::surfaceScalarField& Foam::fvMesh::phiRef()
477 {
478  if (!phiPtr_)
479  {
481  << "mesh flux field does not exist, is the mesh actually moving?"
482  << abort(FatalError);
483  }
484 
485  return *phiPtr_;
486 }
487 
488 
489 // ************************************************************************* //
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
Graphite solid properties.
Definition: C.H:51
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricBoundaryField class.
Generic GeometricField class.
tmp< GeometricField< Type, PatchField, GeoMesh > > cloneUnSliced() const
Clone un-sliced.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:167
Specialisation of DimensionedField which holds a slice of a given complete field in such a form that ...
The time value with time-stepping information, user-defined remapping, etc.
Definition: TimeState.H:55
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
const Type & value() const
Return const reference to value.
const volVectorField & C() const
Return cell centres.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
const surfaceScalarField & phi() const
Return cell face motion fluxes.
const surfaceVectorField & Cf() const
Return face centres.
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
const surfaceVectorField & Sf() const
Return cell face area vectors.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:965
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:270
const vectorField & faceAreas() const
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
#define InfoInFunction
Report an information message using Foam::Info.
SlicedGeometricField< vector, fvsPatchField, slicedFvsPatchField, surfaceMesh > slicedSurfaceVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
SlicedGeometricField< scalar, fvsPatchField, slicedFvsPatchField, surfaceMesh > slicedSurfaceScalarField
void deleteDemandDrivenData(DataType *&dataPtr)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimLength
SlicedGeometricField< vector, fvPatchField, slicedFvPatchField, volMesh > slicedVolVectorField
const dimensionSet dimTime
const dimensionSet dimVolume
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
error FatalError
const dimensionSet dimArea
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label timeIndex
Definition: getTimeIndex.H:4
faceListList boundary(nPatches)
Foam::surfaceFields.