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-2019 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 (SfPtr_)
48  {
50  << "face areas already exist"
51  << abort(FatalError);
52  }
53 
54  SfPtr_ = new slicedSurfaceVectorField
55  (
56  IOobject
57  (
58  "S",
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 (magSfPtr_)
83  {
85  << "mag face areas already exist"
86  << abort(FatalError);
87  }
88 
89  // Note: Added stabilisation for faces with exactly zero area.
90  // These should be caught on mesh checking but at least this stops
91  // the code from producing Nans.
92  magSfPtr_ = new surfaceScalarField
93  (
94  IOobject
95  (
96  "magSf",
98  meshSubDir,
99  *this,
102  false
103  ),
104  mag(Sf()) + dimensionedScalar(dimArea, vSmall)
105  );
106 }
107 
108 
109 void Foam::fvMesh::makeC() const
110 {
111  if (debug)
112  {
113  InfoInFunction << "Assembling cell centres" << endl;
114  }
115 
116  // It is an error to attempt to recalculate
117  // if the pointer is already set
118  if (CPtr_)
119  {
121  << "cell centres already exist"
122  << abort(FatalError);
123  }
124 
125  // Construct as slices. Only preserve processor (not e.g. cyclic)
126 
127  CPtr_ = new slicedVolVectorField
128  (
129  IOobject
130  (
131  "C",
132  pointsInstance(),
133  meshSubDir,
134  *this,
137  false
138  ),
139  *this,
140  dimLength,
141  cellCentres(),
142  faceCentres(),
143  true, // preserveCouples
144  true // preserveProcOnly
145  );
146 }
147 
148 
149 void Foam::fvMesh::makeCf() const
150 {
151  if (debug)
152  {
153  InfoInFunction << "Assembling face centres" << endl;
154  }
155 
156  // It is an error to attempt to recalculate
157  // if the pointer is already set
158  if (CfPtr_)
159  {
161  << "face centres already exist"
162  << abort(FatalError);
163  }
164 
165  CfPtr_ = new slicedSurfaceVectorField
166  (
167  IOobject
168  (
169  "Cf",
170  pointsInstance(),
171  meshSubDir,
172  *this,
175  false
176  ),
177  *this,
178  dimLength,
179  faceCentres()
180  );
181 }
182 
183 
184 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
185 
187 {
188  if (!VPtr_)
189  {
190  if (debug)
191  {
193  << "Constructing from primitiveMesh::cellVolumes()" << endl;
194  }
195 
197  (
198  IOobject
199  (
200  "V",
201  time().timeName(),
202  *this,
205  false
206  ),
207  *this,
208  dimVolume,
209  cellVolumes()
210  );
211  }
212 
213  return *static_cast<slicedVolScalarField::Internal*>(VPtr_);
214 }
215 
216 
218 {
219  if (!V0Ptr_)
220  {
222  << "V0 is not available"
223  << abort(FatalError);
224  }
225 
226  return *V0Ptr_;
227 }
228 
229 
231 {
232  if (!V0Ptr_)
233  {
235  << "V0 is not available"
236  << abort(FatalError);
237  }
238 
239  return *V0Ptr_;
240 }
241 
242 
244 {
245  if (!V00Ptr_)
246  {
247  if (debug)
248  {
249  InfoInFunction << "Constructing from V0" << endl;
250  }
251 
253  (
254  IOobject
255  (
256  "V00",
257  time().timeName(),
258  *this,
261  false
262  ),
263  V0()
264  );
265  }
266 
267  return *V00Ptr_;
268 }
269 
270 
273 {
274  if (moving() && time().subCycling())
275  {
276  const TimeState& ts = time();
277  const TimeState& ts0 = time().prevTimeState();
278 
279  scalar tFrac =
280  (
281  ts.value() - (ts0.value() - ts0.deltaTValue())
282  )/ts0.deltaTValue();
283 
284  if (tFrac < (1 - small))
285  {
286  return V0() + tFrac*(V() - V0());
287  }
288  else
289  {
290  return V();
291  }
292  }
293  else
294  {
295  return V();
296  }
297 }
298 
299 
302 {
303  if (moving() && time().subCycling())
304  {
305  const TimeState& ts = time();
306  const TimeState& ts0 = time().prevTimeState();
307 
308  scalar t0Frac =
309  (
310  (ts.value() - ts.deltaTValue())
311  - (ts0.value() - ts0.deltaTValue())
312  )/ts0.deltaTValue();
313 
314  if (t0Frac > small)
315  {
316  return V0() + t0Frac*(V() - V0());
317  }
318  else
319  {
320  return V0();
321  }
322  }
323  else
324  {
325  return V0();
326  }
327 }
328 
329 
331 {
332  if (!SfPtr_)
333  {
334  makeSf();
335  }
336 
337  return *SfPtr_;
338 }
339 
340 
342 {
343  if (!magSfPtr_)
344  {
345  makeMagSf();
346  }
347 
348  return *magSfPtr_;
349 }
350 
351 
353 {
354  if (!CPtr_)
355  {
356  makeC();
357  }
358 
359  return *CPtr_;
360 }
361 
362 
364 {
365  if (!CfPtr_)
366  {
367  makeCf();
368  }
369 
370  return *CfPtr_;
371 }
372 
373 
375 {
376  if (debug)
377  {
378  InfoInFunction << "Calculating face deltas" << endl;
379  }
380 
382  (
384  (
385  "delta",
386  *this,
387  dimLength
388  )
389  );
390  surfaceVectorField& delta = tdelta.ref();
391 
392  const volVectorField& C = this->C();
393  const labelUList& owner = this->owner();
394  const labelUList& neighbour = this->neighbour();
395 
396  forAll(owner, facei)
397  {
398  delta[facei] = C[neighbour[facei]] - C[owner[facei]];
399  }
400 
401  surfaceVectorField::Boundary& deltabf =
402  delta.boundaryFieldRef();
403 
404  forAll(deltabf, patchi)
405  {
406  deltabf[patchi] = boundary()[patchi].delta();
407  }
408 
409  return tdelta;
410 }
411 
412 
414 {
415  if (!phiPtr_)
416  {
418  << "mesh flux field does not exist, is the mesh actually moving?"
419  << abort(FatalError);
420  }
421 
422  // Set zero current time
423  // mesh motion fluxes if the time has been incremented
424  if (!time().subCycling() && phiPtr_->timeIndex() != time().timeIndex())
425  {
426  (*phiPtr_) = dimensionedScalar(dimVolume/dimTime, 0);
427  }
428 
429  return *phiPtr_;
430 }
431 
432 
434 {
435  if (!phiPtr_)
436  {
438  << "mesh flux field does not exist, is the mesh actually moving?"
439  << abort(FatalError);
440  }
441 
442  return *phiPtr_;
443 }
444 
445 
446 // ************************************************************************* //
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.
bool moving() const
Is mesh moving.
Definition: polyMesh.H:512
The time value with time-stepping information, user-defined remapping, etc.
Definition: TimeState.H:48
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
static int debug
Debug switch.
Definition: fvSchemes.H:97
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
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:209
const surfaceScalarField & phi() const
Return cell face motion fluxes.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
DimensionedField< scalar, volMesh > & setV0()
Return old-time cell volumes.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
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:434
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:284
static tmp< GeometricField< vector, fvsPatchField, surfaceMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvsPatchField< vector >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
surfaceScalarField & setPhi()
Return cell face motion fluxes.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
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:802
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:41
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-cycl 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:278
label timeIndex() const
Return the time index of the field.
The internalField of a SlicedGeometricField.
const vectorField & faceCentres() const
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
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 dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dimensioned< scalar > mag(const dimensioned< Type > &)
const volVectorField & C() const
Return cell centres as volVectorField.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
const scalarField & cellVolumes() const
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:540
#define InfoInFunction
Report an information message using Foam::Info.