fvMeshGeometry.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 | Copyright (C) 2011-2016 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("vs", 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  // If V00 is used then V0 should be stored for restart
267  V0Ptr_->writeOpt() = IOobject::AUTO_WRITE;
268  }
269 
270  return *V00Ptr_;
271 }
272 
273 
276 {
277  if (moving() && time().subCycling())
278  {
279  const TimeState& ts = time();
280  const TimeState& ts0 = time().prevTimeState();
281 
282  scalar tFrac =
283  (
284  ts.value() - (ts0.value() - ts0.deltaTValue())
285  )/ts0.deltaTValue();
286 
287  if (tFrac < (1 - SMALL))
288  {
289  return V0() + tFrac*(V() - V0());
290  }
291  else
292  {
293  return V();
294  }
295  }
296  else
297  {
298  return V();
299  }
300 }
301 
302 
305 {
306  if (moving() && time().subCycling())
307  {
308  const TimeState& ts = time();
309  const TimeState& ts0 = time().prevTimeState();
310 
311  scalar t0Frac =
312  (
313  (ts.value() - ts.deltaTValue())
314  - (ts0.value() - ts0.deltaTValue())
315  )/ts0.deltaTValue();
316 
317  if (t0Frac > SMALL)
318  {
319  return V0() + t0Frac*(V() - V0());
320  }
321  else
322  {
323  return V0();
324  }
325  }
326  else
327  {
328  return V0();
329  }
330 }
331 
332 
334 {
335  if (!SfPtr_)
336  {
337  makeSf();
338  }
339 
340  return *SfPtr_;
341 }
342 
343 
345 {
346  if (!magSfPtr_)
347  {
348  makeMagSf();
349  }
350 
351  return *magSfPtr_;
352 }
353 
354 
356 {
357  if (!CPtr_)
358  {
359  makeC();
360  }
361 
362  return *CPtr_;
363 }
364 
365 
367 {
368  if (!CfPtr_)
369  {
370  makeCf();
371  }
372 
373  return *CfPtr_;
374 }
375 
376 
378 {
379  if (debug)
380  {
381  InfoInFunction << "Calculating face deltas" << endl;
382  }
383 
385  (
387  (
388  IOobject
389  (
390  "delta",
391  pointsInstance(),
392  meshSubDir,
393  *this,
396  false
397  ),
398  *this,
399  dimLength
400  )
401  );
402  surfaceVectorField& delta = tdelta.ref();
403 
404  const volVectorField& C = this->C();
405  const labelUList& owner = this->owner();
406  const labelUList& neighbour = this->neighbour();
407 
408  forAll(owner, facei)
409  {
410  delta[facei] = C[neighbour[facei]] - C[owner[facei]];
411  }
412 
413  surfaceVectorField::Boundary& deltabf =
414  delta.boundaryFieldRef();
415 
416  forAll(deltabf, patchi)
417  {
418  deltabf[patchi] = boundary()[patchi].delta();
419  }
420 
421  return tdelta;
422 }
423 
424 
426 {
427  if (!phiPtr_)
428  {
430  << "mesh flux field does not exist, is the mesh actually moving?"
431  << abort(FatalError);
432  }
433 
434  // Set zero current time
435  // mesh motion fluxes if the time has been incremented
436  if (phiPtr_->timeIndex() != time().timeIndex())
437  {
438  (*phiPtr_) = dimensionedScalar("0", dimVolume/dimTime, 0.0);
439  }
440 
441  return *phiPtr_;
442 }
443 
444 
446 {
447  if (!phiPtr_)
448  {
450  << "mesh flux field does not exist, is the mesh actually moving?"
451  << abort(FatalError);
452  }
453 
454  return *phiPtr_;
455 }
456 
457 
458 // ************************************************************************* //
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
Foam::surfaceFields.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
Graphite solid properties.
Definition: C.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const surfaceVectorField & Sf() const
Return cell face area vectors.
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
static int debug
Debug switch.
Definition: fvSchemes.H:103
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
const surfaceScalarField & phi() const
Return cell face motion fluxes.
const vectorField & faceAreas() const
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
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:114
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
const Type & value() const
Return const reference to value.
const volVectorField & C() const
Return cell centres as volVectorField.
DimensionedField< scalar, volMesh > & setV0()
Return old-time cell volumes.
const vectorField & faceCentres() const
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
bool moving() const
Is mesh moving.
Definition: polyMesh.H:493
surfaceScalarField & setPhi()
Return cell face motion fluxes.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:288
SlicedGeometricField< vector, fvsPatchField, slicedFvsPatchField, surfaceMesh > slicedSurfaceVectorField
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
label timeIndex() const
Return the time index of the field.
const vectorField & cellCentres() const
word timeName
Definition: getTimeIndex.H:3
SlicedGeometricField< vector, fvPatchField, slicedFvPatchField, volMesh > slicedVolVectorField
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
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
const scalarField & cellVolumes() const
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycl old-time cell volumes.
The internalField of a SlicedGeometricField.
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 fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dimensioned< scalar > mag(const dimensioned< Type > &)
const TimeState & prevTimeState() const
Return previous TimeState if time is being sub-cycled.
Definition: Time.H:450
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:282
A class for managing temporary objects.
Definition: PtrList.H:54
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
#define InfoInFunction
Report an information message using Foam::Info.