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  "C",
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::volScalarField::Internal& Foam::fvMesh::V0Ref()
183 {
184  if (!V0Ptr_)
185  {
187  << "V0 is not available"
188  << abort(FatalError);
189  }
190 
191  return *V0Ptr_;
192 }
193 
194 
195 Foam::volScalarField::Internal& Foam::fvMesh::V00Ref()
196 {
197  if (!V00Ptr_)
198  {
199  V00();
200  }
201 
202  return *V00Ptr_;
203 }
204 
205 
206 Foam::surfaceVectorField& Foam::fvMesh::SfRef()
207 {
208  if (!SfPtr_)
209  {
210  SfPtr_ = Sf().cloneUnSliced().ptr();
211 
212  deleteDemandDrivenData(SfSlicePtr_);
213  }
214 
215  return *SfPtr_;
216 }
217 
218 
219 Foam::surfaceScalarField& Foam::fvMesh::magSfRef()
220 {
221  if (!magSfPtr_)
222  {
223  magSfPtr_ = magSf().cloneUnSliced().ptr();
224 
225  deleteDemandDrivenData(magSfSlicePtr_);
226  }
227 
228  return *magSfPtr_;
229 }
230 
231 
232 Foam::volVectorField& Foam::fvMesh::CRef()
233 {
234  if (!CPtr_)
235  {
236  CPtr_ = C().cloneUnSliced().ptr();
237 
238  deleteDemandDrivenData(CSlicePtr_);
239  }
240 
241  return *CPtr_;
242 }
243 
244 
245 Foam::surfaceVectorField& Foam::fvMesh::CfRef()
246 {
247  if (!CfPtr_)
248  {
249  CfPtr_ = Cf().cloneUnSliced().ptr();
250 
251  deleteDemandDrivenData(CfSlicePtr_);
252  }
253 
254  return *CfPtr_;
255 }
256 
257 
258 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
259 
261 {
262  if (!VPtr_)
263  {
264  if (debug)
265  {
267  << "Constructing from primitiveMesh::cellVolumes()" << endl;
268  }
269 
271  (
272  IOobject
273  (
274  "V",
275  time().name(),
276  *this,
279  true
280  ),
281  *this,
282  dimVolume,
283  cellVolumes()
284  );
285  }
286 
287  return *VPtr_;
288 }
289 
290 
292 {
293  if (!V0Ptr_)
294  {
296  << "V0 is not available"
297  << abort(FatalError);
298  }
299 
300  return *V0Ptr_;
301 }
302 
303 
305 {
306  if (!V00Ptr_)
307  {
308  if (debug)
309  {
310  InfoInFunction << "Constructing from V0" << endl;
311  }
312 
314  (
315  IOobject
316  (
317  "V00",
318  time().name(),
319  *this,
322  false
323  ),
324  V0()
325  );
326  }
327 
328  return *V00Ptr_;
329 }
330 
331 
334 {
335  if (moving() && time().subCycling())
336  {
337  const TimeState& ts = time();
338  const TimeState& ts0 = time().prevTimeState();
339 
340  scalar tFrac =
341  (
342  ts.value() - (ts0.value() - ts0.deltaTValue())
343  )/ts0.deltaTValue();
344 
345  if (tFrac < (1 - small))
346  {
347  return V0() + tFrac*(V() - V0());
348  }
349  else
350  {
351  return V();
352  }
353  }
354  else
355  {
356  return V();
357  }
358 }
359 
360 
363 {
364  if (moving() && time().subCycling())
365  {
366  const TimeState& ts = time();
367  const TimeState& ts0 = time().prevTimeState();
368 
369  scalar t0Frac =
370  (
371  (ts.value() - ts.deltaTValue())
372  - (ts0.value() - ts0.deltaTValue())
373  )/ts0.deltaTValue();
374 
375  if (t0Frac > small)
376  {
377  return V0() + t0Frac*(V() - V0());
378  }
379  else
380  {
381  return V0();
382  }
383  }
384  else
385  {
386  return V0();
387  }
388 }
389 
390 
392 {
393  if (SfPtr_)
394  {
395  return *SfPtr_;
396  }
397 
398  if (!SfSlicePtr_)
399  {
400  makeSf();
401  }
402 
403  return *SfSlicePtr_;
404 }
405 
406 
408 {
409  if (magSfPtr_)
410  {
411  return *magSfPtr_;
412  }
413 
414  if (!magSfSlicePtr_)
415  {
416  makeMagSf();
417  }
418 
419  return *magSfSlicePtr_;
420 }
421 
422 
424 {
425  if (CPtr_)
426  {
427  return *CPtr_;
428  }
429 
430  if (!CSlicePtr_)
431  {
432  makeC();
433  }
434 
435  return *CSlicePtr_;
436 }
437 
438 
440 {
441  if (CfPtr_)
442  {
443  return *CfPtr_;
444  }
445 
446  if (!CfSlicePtr_)
447  {
448  makeCf();
449  }
450 
451  return *CfSlicePtr_;
452 }
453 
454 
456 {
457  if (debug)
458  {
459  InfoInFunction << "Calculating face deltas" << endl;
460  }
461 
463  (
465  (
466  "delta",
467  *this,
468  dimLength
469  )
470  );
471  surfaceVectorField& delta = tdelta.ref();
472 
473  const volVectorField& C = this->C();
474  const labelUList& owner = this->owner();
475  const labelUList& neighbour = this->neighbour();
476 
477  forAll(owner, facei)
478  {
479  delta[facei] = C[neighbour[facei]] - C[owner[facei]];
480  }
481 
483  delta.boundaryFieldRef();
484 
485  forAll(deltabf, patchi)
486  {
487  deltabf[patchi] = boundary()[patchi].delta();
488  }
489 
490  return tdelta;
491 }
492 
493 
495 {
496  if (!phiPtr_)
497  {
499  << "mesh flux field does not exist, is the mesh actually moving?"
500  << abort(FatalError);
501  }
502 
503  // Set zero current time
504  // mesh motion fluxes if the time has been incremented
505  if (!time().subCycling() && phiPtr_->timeIndex() != time().timeIndex())
506  {
507  (*phiPtr_) = dimensionedScalar(dimVolume/dimTime, 0);
508  }
509 
510  return *phiPtr_;
511 }
512 
513 
514 Foam::surfaceScalarField& Foam::fvMesh::phiRef()
515 {
516  if (!phiPtr_)
517  {
519  << "mesh flux field does not exist, is the mesh actually moving?"
520  << abort(FatalError);
521  }
522 
523  return *phiPtr_;
524 }
525 
526 
527 // ************************************************************************* //
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 >> &)
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.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:269
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:1020
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:306
label patchi
#define InfoInFunction
Report an information message using Foam::Info.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
SlicedGeometricField< vector, fvsPatchField, slicedFvsPatchField, surfaceMesh > slicedSurfaceVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
SlicedGeometricField< scalar, fvsPatchField, slicedFvsPatchField, surfaceMesh > slicedSurfaceScalarField
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimLength
void deleteDemandDrivenData(DataPtr &dataPtr)
SlicedGeometricField< vector, fvPatchField, slicedFvPatchField, volMesh > slicedVolVectorField
const dimensionSet dimTime
const dimensionSet dimVolume
error FatalError
const dimensionSet dimArea
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label timeIndex
Definition: getTimeIndex.H:4
faceListList boundary(nPatches)
Foam::surfaceFields.