Loading...
Searching...
No Matches
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-2024 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"
35#include "SubField.H"
36#include "cyclicFvPatchFields.H"
38
39// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40
41void Foam::fvMesh::makeSf() const
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
54 SfPtr_ = std::make_unique<slicedSurfaceVectorField>
55 (
57 (
58 "S",
61 *this,
65 ),
66 *this,
67 dimArea,
69 );
70
71 SfPtr_->setOriented();
72}
73
74
75void Foam::fvMesh::makeMagSf() const
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_ = std::make_unique<surfaceScalarField>
92 (
93 IOobject
94 (
95 "magSf",
96 pointsInstance(),
97 meshSubDir,
98 *this,
102 ),
103 mag(Sf()) + dimensionedScalar("vs", dimArea, VSMALL)
104 );
105}
106
107
108void Foam::fvMesh::makeC() const
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_ = std::make_unique<slicedVolVectorField>
124 (
125 IOobject
126 (
127 "C",
128 pointsInstance(),
129 meshSubDir,
130 *this,
134 ),
135 *this,
136 dimLength,
137 cellCentres(),
139 true, //preserveCouples
140 true //preserveProcOnly
141 );
142}
143
144
145void Foam::fvMesh::makeCf() const
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_ = std::make_unique<slicedSurfaceVectorField>
159 (
160 IOobject
161 (
162 "Cf",
163 pointsInstance(),
164 meshSubDir,
165 *this,
169 ),
170 *this,
171 dimLength,
173 );
174}
175
176
177// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
178
180{
181 if (!VPtr_)
182 {
184 << "Constructing from primitiveMesh::cellVolumes()" << endl;
185
186 VPtr_ = std::make_unique<SlicedDimensionedField<scalar, volMesh>>
187 (
189 (
190 "V",
191 time().timeName(),
192 *this,
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"
214 }
215
216 return *V0Ptr_;
217}
218
219
221{
222 if (!V0Ptr_)
223 {
225 << "V0 is not available"
227 }
228
229 return *V0Ptr_;
230}
231
232
234{
235 if (!V00Ptr_)
236 {
237 DebugInFunction << "Constructing from V0" << endl;
238
239 V00Ptr_ = std::make_unique<DimensionedField<scalar, volMesh>>
240 (
241 IOobject
242 (
243 "V00",
244 time().timeName(),
245 *this,
249 ),
250 V0()
251 );
252
253
254 // If V00 is used then V0 should be stored for restart
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());
278 }
279
280 return V();
281}
282
283
284Foam::tmp<Foam::volScalarField::Internal> Foam::fvMesh::Vsc0() const
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());
301 }
302
303 return V0();
304}
305
306
308{
309 if (!SfPtr_)
310 {
312 }
313
314 return *SfPtr_;
315}
316
317
319{
320 if (!magSfPtr_)
321 {
323 }
324
325 return *magSfPtr_;
326}
327
328
330{
331 auto tunitVectors = tmp<surfaceVectorField>::New
332 (
333 IOobject
334 (
335 "unit(Sf)",
336 pointsInstance(),
337 meshSubDir,
338 *this,
342 ),
343 *this,
344 dimless,
345 (this->Sf() / this->magSf())
346 );
347
348 tunitVectors.ref().oriented() = this->Sf().oriented();
349 return tunitVectors;
350}
351
352
354{
355 if (!CPtr_)
356 {
358 }
359
360 return *CPtr_;
361}
362
363
365{
366 if (!CfPtr_)
367 {
369 }
370
371 return *CfPtr_;
372}
373
374
376{
377 DebugInFunction << "Calculating face deltas" << endl;
378
379 auto tdelta = tmp<surfaceVectorField>::New
380 (
381 IOobject
382 (
383 "delta",
384 pointsInstance(),
385 meshSubDir,
386 *this,
390 ),
391 *this,
393 );
394 auto& delta = tdelta.ref();
395 delta.setOriented();
396
397 const volVectorField& C = this->C();
398 const labelUList& owner = this->owner();
399 const labelUList& neighbour = this->neighbour();
400
401 forAll(owner, facei)
402 {
403 delta[facei] = C[neighbour[facei]] - C[owner[facei]];
404 }
405
406 auto& deltabf = delta.boundaryFieldRef();
407
408 forAll(deltabf, patchi)
409 {
410 deltabf[patchi] = boundary()[patchi].delta();
411 }
412
413 return tdelta;
414}
415
416
418{
419 if (!phiPtr_)
420 {
422 << "mesh flux field does not exist, is the mesh actually moving?"
423 << abort(FatalError);
424 }
425
426 // Set zero current time
427 // mesh motion fluxes if the time has been incremented
428 if (!time().subCycling() && phiPtr_->timeIndex() != time().timeIndex())
429 {
430 (*phiPtr_) = Zero;
431 }
433 phiPtr_->setOriented();
434
435 return *phiPtr_;
436}
437
438
440{
442
443 // Return non-const reference, or nullptr if not available
444 phiref.ref(phiPtr_.get());
445
446 return phiref;
447}
448
449
450// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
scalar delta
DimensionedField< scalar, volMesh > Internal
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
The time value with time-stepping information, user-defined remapping, etc.
Definition TimeState.H:50
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
const Type & value() const noexcept
Return const reference to value.
std::unique_ptr< slicedSurfaceVectorField > SfPtr_
Face area vectors.
Definition fvMesh.H:127
const volVectorField & C() const
Return cell centres as volVectorField.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
std::unique_ptr< surfaceScalarField > magSfPtr_
Mag face area vectors.
Definition fvMesh.H:132
DimensionedField< scalar, volMesh > & setV0()
Return old-time cell volumes.
void makeMagSf() const
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
void makeC() const
const surfaceScalarField & phi() const
Return cell face motion fluxes.
std::unique_ptr< DimensionedField< scalar, volMesh > > V00Ptr_
Cell volumes old-old time level.
Definition fvMesh.H:122
std::unique_ptr< DimensionedField< scalar, volMesh > > V0Ptr_
Cell volumes old time level.
Definition fvMesh.H:117
void makeSf() const
void makeCf() const
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition fvMesh.H:572
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
tmp< surfaceVectorField > unitSf() const
Return cell face unit normals.
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
const labelUList & neighbour() const
Internal face neighbour.
Definition fvMesh.H:580
const surfaceVectorField & Sf() const
Return cell face area vectors.
std::unique_ptr< slicedVolVectorField > CPtr_
Cell centres.
Definition fvMesh.H:137
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
std::unique_ptr< slicedSurfaceVectorField > CfPtr_
Face centres.
Definition fvMesh.H:142
refPtr< surfaceScalarField > setPhi()
Return cell face motion fluxes, if any (can be nullptr).
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
std::unique_ptr< surfaceScalarField > phiPtr_
Face motion fluxes.
Definition fvMesh.H:147
std::unique_ptr< SlicedDimensionedField< scalar, volMesh > > VPtr_
Cell volumes.
Definition fvMesh.H:112
bool moving() const noexcept
Is mesh moving.
Definition polyMesh.H:732
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition polyMesh.C:838
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
const vectorField & faceCentres() const
const scalarField & cellVolumes() const
const vectorField & cellCentres() const
const vectorField & faceAreas() const
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition refPtrI.H:230
bool steady() const noexcept
True if default ddt scheme is steady-state.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
faceListList boundary
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
word timeName
Definition getTimeIndex.H:3
#define DebugInFunction
Report an information message using Foam::Info.
const expr V(m.psi().mesh().V())
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimArea(sqr(dimLength))
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
label timeIndex
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.