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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017,2022 OpenFOAM Foundation
9 Copyright (C) 2020-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "fvMesh.H"
30#include "Time.H"
31#include "volFields.H"
32#include "surfaceFields.H"
33#include "slicedVolFields.H"
34#include "slicedSurfaceFields.H"
35#include "SubField.H"
36#include "cyclicFvPatchFields.H"
38
39// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40
42{
43 DebugInFunction << "Assembling face areas" << endl;
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
55 (
57 (
58 "S",
61 *this,
64 false
65 ),
66 *this,
67 dimArea,
68 faceAreas()
69 );
70
72}
73
74
76{
77 DebugInFunction << "Assembling mag face areas" << endl;
78
79 // It is an error to attempt to recalculate
80 // if the pointer is already set
81 if (magSfPtr_)
82 {
84 << "mag face areas already exist"
85 << abort(FatalError);
86 }
87
88 // Note: Added stabilisation for faces with exactly zero area.
89 // These should be caught on mesh checking but at least this stops
90 // the code from producing Nans.
91 magSfPtr_ = new surfaceScalarField
92 (
94 (
95 "magSf",
96 pointsInstance(),
97 meshSubDir,
98 *this,
101 false
102 ),
103 mag(Sf()) + dimensionedScalar("vs", dimArea, VSMALL)
104 );
105}
106
107
109{
110 DebugInFunction << "Assembling cell centres" << endl;
111
112 // It is an error to attempt to recalculate
113 // if the pointer is already set
114 if (CPtr_)
115 {
117 << "cell centres already exist"
118 << abort(FatalError);
119 }
120
121 // Construct as slices. Only preserve processor (not e.g. cyclic)
122
123 CPtr_ = new slicedVolVectorField
124 (
126 (
127 "C",
128 pointsInstance(),
129 meshSubDir,
130 *this,
133 false
134 ),
135 *this,
136 dimLength,
137 cellCentres(),
138 faceCentres(),
139 true, //preserveCouples
140 true //preserveProcOnly
141 );
142}
143
144
146{
147 DebugInFunction << "Assembling face centres" << endl;
148
149 // It is an error to attempt to recalculate
150 // if the pointer is already set
151 if (CfPtr_)
152 {
154 << "face centres already exist"
155 << abort(FatalError);
156 }
157
158 CfPtr_ = new slicedSurfaceVectorField
159 (
161 (
162 "Cf",
163 pointsInstance(),
164 meshSubDir,
165 *this,
168 false
169 ),
170 *this,
171 dimLength,
172 faceCentres()
173 );
174}
175
176
177// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
178
180{
181 if (!VPtr_)
182 {
184 << "Constructing from primitiveMesh::cellVolumes()" << endl;
185
187 (
189 (
190 "V",
191 time().timeName(),
192 *this,
195 false
196 ),
197 *this,
198 dimVolume,
199 cellVolumes()
200 );
201 }
202
203 return *VPtr_;
204}
205
206
208{
209 if (!V0Ptr_)
210 {
212 << "V0 is not available"
213 << abort(FatalError);
214 }
215
216 return *V0Ptr_;
217}
218
219
221{
222 if (!V0Ptr_)
223 {
225 << "V0 is not available"
226 << abort(FatalError);
227 }
228
229 return *V0Ptr_;
230}
231
232
234{
235 if (!V00Ptr_)
236 {
237 DebugInFunction << "Constructing from V0" << endl;
238
240 (
242 (
243 "V00",
244 time().timeName(),
245 *this,
248 true
249 ),
250 V0()
251 );
252
253
254 // If V00 is used then V0 should be stored for restart
255 V0Ptr_->writeOpt(IOobject::AUTO_WRITE);
256 }
257
258 return *V00Ptr_;
259}
260
261
263{
264 if (!steady() && moving() && time().subCycling())
265 {
266 const TimeState& ts = time();
267 const TimeState& ts0 = time().prevTimeState();
268
269 scalar tFrac =
270 (
271 ts.value() - (ts0.value() - ts0.deltaTValue())
272 )/ts0.deltaTValue();
273
274 if (tFrac < (1 - SMALL))
275 {
276 return V0() + tFrac*(V() - V0());
277 }
278 }
279
280 return V();
281}
282
283
285{
286 if (!steady() && moving() && time().subCycling())
287 {
288 const TimeState& ts = time();
289 const TimeState& ts0 = time().prevTimeState();
290
291 scalar t0Frac =
292 (
293 (ts.value() - ts.deltaTValue())
294 - (ts0.value() - ts0.deltaTValue())
295 )/ts0.deltaTValue();
296
297 if (t0Frac > SMALL)
298 {
299 return V0() + t0Frac*(V() - V0());
300 }
301 }
302
303 return V0();
304}
305
306
308{
309 if (!SfPtr_)
310 {
311 makeSf();
312 }
313
314 return *SfPtr_;
315}
316
317
319{
320 if (!magSfPtr_)
321 {
322 makeMagSf();
323 }
324
325 return *magSfPtr_;
326}
327
328
330{
331 if (!CPtr_)
332 {
333 makeC();
334 }
335
336 return *CPtr_;
337}
338
339
341{
342 if (!CfPtr_)
343 {
344 makeCf();
345 }
346
347 return *CfPtr_;
348}
349
350
352{
353 DebugInFunction << "Calculating face deltas" << endl;
354
356 (
358 (
360 (
361 "delta",
362 pointsInstance(),
363 meshSubDir,
364 *this,
367 false
368 ),
369 *this,
371 )
372 );
373 surfaceVectorField& delta = tdelta.ref();
374 delta.setOriented();
375
376 const volVectorField& C = this->C();
377 const labelUList& owner = this->owner();
378 const labelUList& neighbour = this->neighbour();
379
380 forAll(owner, facei)
381 {
382 delta[facei] = C[neighbour[facei]] - C[owner[facei]];
383 }
384
385 surfaceVectorField::Boundary& deltabf = delta.boundaryFieldRef();
386
387 forAll(deltabf, patchi)
388 {
389 deltabf[patchi] = boundary()[patchi].delta();
390 }
391
392 return tdelta;
393}
394
395
397{
398 if (!phiPtr_)
399 {
401 << "mesh flux field does not exist, is the mesh actually moving?"
402 << abort(FatalError);
403 }
404
405 // Set zero current time
406 // mesh motion fluxes if the time has been incremented
407 if (!time().subCycling() && phiPtr_->timeIndex() != time().timeIndex())
408 {
410 }
411
412 phiPtr_->setOriented();
413
414 return *phiPtr_;
415}
416
417
419{
420 if (!phiPtr_)
421 {
422 return nullptr;
423 }
424 else
425 {
426 // Return non-const reference
428 p.ref(*phiPtr_);
429 return p;
430 }
431}
432
433
434// ************************************************************************* //
scalar delta
Graphite solid properties.
Definition: C.H:53
void setOriented(const bool oriented=true) noexcept
Set the oriented flag.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Specialisation of DimensionedField that holds a slice of a given field so that it acts as a Dimension...
The time value with time-stepping information, user-defined remapping, etc.
Definition: TimeState.H:54
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:43
const Type & value() const
Return const reference to value.
slicedSurfaceVectorField * SfPtr_
Face area vectors.
Definition: fvMesh.H:119
const volVectorField & C() const
Return cell centres as volVectorField.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
DimensionedField< scalar, volMesh > & setV0()
Return old-time cell volumes.
void makeMagSf() const
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
void makeC() const
const surfaceScalarField & phi() const
Return cell face motion fluxes.
void makeSf() const
void makeCf() const
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
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.
refPtr< surfaceScalarField > setPhi()
Return cell face motion fluxes (or null)
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:860
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:324
const vectorField & faceAreas() const
A class for managing references or pointers (no reference counting)
Definition: refPtr.H:58
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
volScalarField & p
faceListList boundary
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
word timeName
Definition: getTimeIndex.H:3
#define DebugInFunction
Report an information message using Foam::Info.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
SlicedGeometricField< vector, fvsPatchField, slicedFvsPatchField, surfaceMesh > slicedSurfaceVectorField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
SlicedGeometricField< vector, fvPatchField, slicedFvPatchField, volMesh > slicedVolVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
label timeIndex
Definition: getTimeIndex.H:30
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.