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-2022 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 #include "cyclicAMIFvPatchFields.H"
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
38 void Foam::fvMesh::makeSf() const
39 {
40  if (debug)
41  {
42  InfoInFunction << "Assembling face areas" << endl;
43  }
44 
45  // It is an error to attempt to recalculate
46  // if the pointer is already set
47  if (SfSlicePtr_ || SfPtr_)
48  {
50  << "face areas already exist"
51  << abort(FatalError);
52  }
53 
54  SfSlicePtr_ = new slicedSurfaceVectorField
55  (
56  IOobject
57  (
58  "Sf",
60  meshSubDir,
61  *this,
64  false
65  ),
66  *this,
67  dimArea,
68  faceAreas()
69  );
70 }
71 
72 
73 void Foam::fvMesh::makeMagSf() const
74 {
75  if (debug)
76  {
77  InfoInFunction << "Assembling mag face areas" << endl;
78  }
79 
80  // It is an error to attempt to recalculate
81  // if the pointer is already set
82  if (magSfSlicePtr_ || magSfPtr_)
83  {
85  << "mag face areas already exist"
86  << abort(FatalError);
87  }
88 
89  magSfSlicePtr_ = new slicedSurfaceScalarField
90  (
91  IOobject
92  (
93  "magSf",
95  meshSubDir,
96  *this,
99  false
100  ),
101  *this,
102  dimArea,
103  magFaceAreas()
104  );
105 }
106 
107 
108 void Foam::fvMesh::makeC() const
109 {
110  if (debug)
111  {
112  InfoInFunction << "Assembling cell centres" << endl;
113  }
114 
115  // It is an error to attempt to recalculate
116  // if the pointer is already set
117  if (CSlicePtr_ || CPtr_)
118  {
120  << "cell centres already exist"
121  << abort(FatalError);
122  }
123 
124  // Construct as slices. Only preserve processor (not e.g. cyclic)
125 
126  CSlicePtr_ = new slicedVolVectorField
127  (
128  IOobject
129  (
130  "C",
131  pointsInstance(),
132  meshSubDir,
133  *this,
136  false
137  ),
138  *this,
139  dimLength,
140  cellCentres(),
141  faceCentres(),
142  true, // preserveCouples
143  true // preserveProcOnly
144  );
145 }
146 
147 
148 void Foam::fvMesh::makeCf() const
149 {
150  if (debug)
151  {
152  InfoInFunction << "Assembling face centres" << endl;
153  }
154 
155  // It is an error to attempt to recalculate
156  // if the pointer is already set
157  if (CfSlicePtr_ || CfPtr_)
158  {
160  << "face centres already exist"
161  << abort(FatalError);
162  }
163 
164  CfSlicePtr_ = new slicedSurfaceVectorField
165  (
166  IOobject
167  (
168  "Cf",
169  pointsInstance(),
170  meshSubDir,
171  *this,
174  false
175  ),
176  *this,
177  dimLength,
178  faceCentres()
179  );
180 }
181 
182 
183 Foam::volScalarField::Internal& Foam::fvMesh::V0Ref()
184 {
185  if (!V0Ptr_)
186  {
188  << "V0 is not available"
189  << abort(FatalError);
190  }
191 
192  return *V0Ptr_;
193 }
194 
195 
196 Foam::volScalarField::Internal& Foam::fvMesh::V00Ref()
197 {
198  if (!V00Ptr_)
199  {
200  V00();
201  }
202 
203  return *V00Ptr_;
204 }
205 
206 
207 Foam::surfaceVectorField& Foam::fvMesh::SfRef()
208 {
209  if (!SfPtr_)
210  {
211  SfPtr_ = Sf().cloneUnSliced().ptr();
212 
213  deleteDemandDrivenData(SfSlicePtr_);
214  }
215 
216  return *SfPtr_;
217 }
218 
219 
220 Foam::surfaceScalarField& Foam::fvMesh::magSfRef()
221 {
222  if (!magSfPtr_)
223  {
224  magSfPtr_ = magSf().cloneUnSliced().ptr();
225 
226  deleteDemandDrivenData(magSfSlicePtr_);
227  }
228 
229  return *magSfPtr_;
230 }
231 
232 
233 Foam::volVectorField& Foam::fvMesh::CRef()
234 {
235  if (!CPtr_)
236  {
237  CPtr_ = C().cloneUnSliced().ptr();
238 
239  deleteDemandDrivenData(CSlicePtr_);
240  }
241 
242  return *CPtr_;
243 }
244 
245 
246 Foam::surfaceVectorField& Foam::fvMesh::CfRef()
247 {
248  if (!CfPtr_)
249  {
250  CfPtr_ = Cf().cloneUnSliced().ptr();
251 
252  deleteDemandDrivenData(CfSlicePtr_);
253  }
254 
255  return *CfPtr_;
256 }
257 
258 
259 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
260 
262 {
263  if (!VPtr_)
264  {
265  if (debug)
266  {
268  << "Constructing from primitiveMesh::cellVolumes()" << endl;
269  }
270 
272  (
273  IOobject
274  (
275  "V",
276  time().timeName(),
277  *this,
280  false
281  ),
282  *this,
283  dimVolume,
284  cellVolumes()
285  );
286  }
287 
288  return *VPtr_;
289 }
290 
291 
293 {
294  if (!V0Ptr_)
295  {
297  << "V0 is not available"
298  << abort(FatalError);
299  }
300 
301  return *V0Ptr_;
302 }
303 
304 
306 {
307  if (!V00Ptr_)
308  {
309  if (debug)
310  {
311  InfoInFunction << "Constructing from V0" << endl;
312  }
313 
315  (
316  IOobject
317  (
318  "V00",
319  time().timeName(),
320  *this,
323  false
324  ),
325  V0()
326  );
327  }
328 
329  return *V00Ptr_;
330 }
331 
332 
335 {
336  if (moving() && time().subCycling())
337  {
338  const TimeState& ts = time();
339  const TimeState& ts0 = time().prevTimeState();
340 
341  scalar tFrac =
342  (
343  ts.value() - (ts0.value() - ts0.deltaTValue())
344  )/ts0.deltaTValue();
345 
346  if (tFrac < (1 - small))
347  {
348  return V0() + tFrac*(V() - V0());
349  }
350  else
351  {
352  return V();
353  }
354  }
355  else
356  {
357  return V();
358  }
359 }
360 
361 
364 {
365  if (moving() && time().subCycling())
366  {
367  const TimeState& ts = time();
368  const TimeState& ts0 = time().prevTimeState();
369 
370  scalar t0Frac =
371  (
372  (ts.value() - ts.deltaTValue())
373  - (ts0.value() - ts0.deltaTValue())
374  )/ts0.deltaTValue();
375 
376  if (t0Frac > small)
377  {
378  return V0() + t0Frac*(V() - V0());
379  }
380  else
381  {
382  return V0();
383  }
384  }
385  else
386  {
387  return V0();
388  }
389 }
390 
391 
393 {
394  if (SfPtr_)
395  {
396  return *SfPtr_;
397  }
398 
399  if (!SfSlicePtr_)
400  {
401  makeSf();
402  }
403 
404  return *SfSlicePtr_;
405 }
406 
407 
409 {
410  if (magSfPtr_)
411  {
412  return *magSfPtr_;
413  }
414 
415  if (!magSfSlicePtr_)
416  {
417  makeMagSf();
418  }
419 
420  return *magSfSlicePtr_;
421 }
422 
423 
425 {
426  if (CPtr_)
427  {
428  return *CPtr_;
429  }
430 
431  if (!CSlicePtr_)
432  {
433  makeC();
434  }
435 
436  return *CSlicePtr_;
437 }
438 
439 
441 {
442  if (CfPtr_)
443  {
444  return *CfPtr_;
445  }
446 
447  if (!CfSlicePtr_)
448  {
449  makeCf();
450  }
451 
452  return *CfSlicePtr_;
453 }
454 
455 
457 {
458  if (debug)
459  {
460  InfoInFunction << "Calculating face deltas" << endl;
461  }
462 
464  (
466  (
467  "delta",
468  *this,
469  dimLength
470  )
471  );
472  surfaceVectorField& delta = tdelta.ref();
473 
474  const volVectorField& C = this->C();
475  const labelUList& owner = this->owner();
476  const labelUList& neighbour = this->neighbour();
477 
478  forAll(owner, facei)
479  {
480  delta[facei] = C[neighbour[facei]] - C[owner[facei]];
481  }
482 
484  delta.boundaryFieldRef();
485 
486  forAll(deltabf, patchi)
487  {
488  deltabf[patchi] = boundary()[patchi].delta();
489  }
490 
491  return tdelta;
492 }
493 
494 
496 {
497  if (!phiPtr_)
498  {
500  << "mesh flux field does not exist, is the mesh actually moving?"
501  << abort(FatalError);
502  }
503 
504  // Set zero current time
505  // mesh motion fluxes if the time has been incremented
506  if (!time().subCycling() && phiPtr_->timeIndex() != time().timeIndex())
507  {
508  (*phiPtr_) = dimensionedScalar(dimVolume/dimTime, 0);
509  }
510 
511  return *phiPtr_;
512 }
513 
514 
515 Foam::surfaceScalarField& Foam::fvMesh::phiRef()
516 {
517  if (!phiPtr_)
518  {
520  << "mesh flux field does not exist, is the mesh actually moving?"
521  << abort(FatalError);
522  }
523 
524  return *phiPtr_;
525 }
526 
527 
528 // ************************************************************************* //
Foam::surfaceFields.
Graphite solid properties.
Definition: C.H:48
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const surfaceVectorField & Sf() const
Return cell face area vectors.
const dimensionSet dimArea
bool moving() const
Is mesh moving.
Definition: polyMesh.H:525
The time value with time-stepping information, user-defined remapping, etc.
Definition: TimeState.H:52
Specialisation of DimensionedField which holds a slice of a given complete field in such a form that ...
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const surfaceVectorField & Cf() const
Return face centres.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
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
static tmp< GeometricField< vector, fvsPatchField, surfaceMesh > > New(const word &name, const Internal &, const PtrList< fvsPatchField< vector >> &)
Return a temporary field constructed from name,.
const surfaceScalarField & phi() const
Return cell face motion fluxes.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const TimeState & prevTimeState() const
Return previous TimeState if time is being sub-cycled.
Definition: Time.H:451
const dimensionSet dimLength
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:423
const scalarField & magFaceAreas() const
const dimensionSet dimTime
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:882
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
SlicedGeometricField< vector, fvsPatchField, slicedFvsPatchField, surfaceMesh > slicedSurfaceVectorField
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
const Type & value() const
Return const reference to value.
word timeName
Definition: getTimeIndex.H:3
SlicedGeometricField< vector, fvPatchField, slicedFvPatchField, volMesh > slicedVolVectorField
const vectorField & cellCentres() const
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:417
label timeIndex() const
Return the time index of the field.
const vectorField & faceCentres() const
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
tmp< GeometricField< Type, PatchField, GeoMesh > > cloneUnSliced() const
Clone un-sliced.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const vectorField & faceAreas() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimVolume
const volVectorField & C() const
Return cell centres.
A class for managing temporary objects.
Definition: PtrList.H:53
void deleteDemandDrivenData(DataPtr &dataPtr)
SlicedGeometricField< scalar, fvsPatchField, slicedFvsPatchField, surfaceMesh > slicedSurfaceScalarField
const scalarField & cellVolumes() const
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800
#define InfoInFunction
Report an information message using Foam::Info.