fvMesh.H
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) 2016-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
27Class
28 Foam::fvMesh
29
30Description
31 Mesh data needed to do the Finite Volume discretisation.
32
33 NOTE ON USAGE:
34 fvMesh contains all the topological and geometric information
35 related to the mesh. It is also responsible for keeping the data
36 up-to-date. This is done by deleting the cell volume, face area,
37 cell/face centre, addressing and other derived information as
38 required and recalculating it as necessary. The fvMesh therefore
39 reserves the right to delete the derived information upon every
40 topological (mesh refinement/morphing) or geometric change (mesh
41 motion). It is therefore unsafe to keep local references to the
42 derived data outside of the time loop.
43
44SourceFiles
45 fvMesh.C
46 fvMeshGeometry.C
47
48\*---------------------------------------------------------------------------*/
49
50#ifndef Foam_fvMesh_H
51#define Foam_fvMesh_H
52
53#include "polyMesh.H"
54#include "lduMesh.H"
55#include "primitiveMesh.H"
56#include "fvBoundaryMesh.H"
58#include "fvSchemes.H"
59#include "fvSolution.H"
60#include "data.H"
61#include "volFieldsFwd.H"
62#include "surfaceFieldsFwd.H"
63#include "pointFieldsFwd.H"
65#include "slicedVolFieldsFwd.H"
67#include "className.H"
68#include "SolverPerformance.H"
69
70// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71
72namespace Foam
73{
74
75// Forward Declarations
76class fvMeshLduAddressing;
77class volMesh;
78template<class Type> class fvMatrix;
79
80/*---------------------------------------------------------------------------*\
81 Class fvMesh Declaration
82\*---------------------------------------------------------------------------*/
84class fvMesh
85:
86 public polyMesh,
87 public lduMesh,
88 public fvSchemes,
89 public surfaceInterpolation, // needs input from fvSchemes
90 public fvSolution,
91 public data
92{
93protected:
94
95 // Private data
96
97 //- Boundary mesh
99
100 // Demand-driven data
103
104 //- Current time index for cell volumes
105 // Note. The whole mechanism will be replaced once the
106 // dimensionedField is created and the dimensionedField
107 // will take care of the old-time levels.
108 mutable label curTimeIndex_;
109
110 //- Cell volumes
112
113 //- Cell volumes old time level
115
116 //- Cell volumes old-old time level
118
119 //- Face area vectors
121
122 //- Mag face area vectors
124
125 //- Cell centres
127
128 //- Face centres
130
131 //- Face motion fluxes
133
134
135 // Private Member Functions
136
137 // Storage management
138
139 //- Clear geometry but not the old-time cell volumes
140 void clearGeomNotOldVol();
141
142 //- Clear geometry like clearGeomNotOldVol but recreate any
143 // geometric demand-driven data that was set
144 void updateGeomNotOldVol();
145
146 //- Clear local geometry
147 void clearGeom();
148
149 //- Clear local addressing
150 void clearAddressing(const bool isMeshUpdate = false);
151
152 //- Clear local-only storage (geometry, addressing etc)
153 void clearOutLocal();
154
155 //- Preserve old volume(s)
156 void storeOldVol(const scalarField&);
157
158
159 // Make geometric data
160
161 void makeSf() const;
162 void makeMagSf() const;
163
164 void makeC() const;
165 void makeCf() const;
166
167
168 //- No copy construct
169 fvMesh(const fvMesh&) = delete;
170
171 //- No copy assignment
172 void operator=(const fvMesh&) = delete;
173
174
175public:
176
177 // Public Typedefs
178
179 //- The mesh type
180 typedef fvMesh Mesh;
181
182 //- The boundary type associated with the mesh
184
185
186 // Declare name of the class and its debug switch
187 ClassName("fvMesh");
188
189
190 // Constructors
191
192 //- Construct from IOobject
193 explicit fvMesh(const IOobject& io, const bool doInit=true);
194
195 //- Construct from IOobject or as zero-sized mesh
196 // Boundary is added using addFvPatches() member function
197 fvMesh(const IOobject& io, const Foam::zero, bool syncPar=true);
198
199 //- Construct as copy (for dictionaries) and zero-sized components.
200 // Boundary is added using addFvPatches() member function
201 fvMesh
202 (
203 const IOobject& io,
204 const fvMesh& baseMesh,
205 const Foam::zero,
206 const bool syncPar = true
207 );
208
209 //- Construct from components without boundary.
210 // Boundary is added using addFvPatches() member function
211 fvMesh
212 (
213 const IOobject& io,
215 faceList&& faces,
216 labelList&& allOwner,
217 labelList&& allNeighbour,
218 const bool syncPar = true
219 );
220
221 //- Construct without boundary from cells rather than owner/neighbour.
222 // Boundary is added using addFvPatches() member function
223 fvMesh
224 (
225 const IOobject& io,
227 faceList&& faces,
228 cellList&& cells,
229 const bool syncPar = true
230 );
231
232 //- Copy construct (for dictionaries) with components, without boundary.
233 // Boundary is added using addFvPatches() member function
234 fvMesh
235 (
236 const IOobject& io,
237 const fvMesh& baseMesh,
239 faceList&& faces,
240 labelList&& allOwner,
241 labelList&& allNeighbour,
242 const bool syncPar = true
243 );
244
245 //- Copy construct (for dictionaries) with cells, without boundary.
246 // Boundary is added using addFvPatches() member function
247 fvMesh
248 (
249 const IOobject& io,
250 const fvMesh& baseMesh,
252 faceList&& faces,
253 cellList&& cells,
254 const bool syncPar = true
255 );
256
257
258 //- Destructor
259 virtual ~fvMesh();
260
261
262 // Member Functions
263
264 // Helpers
265
266 //- Initialise all non-demand-driven data
267 virtual bool init(const bool doInit);
268
269 //- Add boundary patches. Constructor helper
270 void addFvPatches
271 (
272 polyPatchList& plist,
273 const bool validBoundary = true
274 );
275
276 //- Add boundary patches. Constructor helper
277 void addFvPatches
278 (
279 const List<polyPatch*>& p,
280 const bool validBoundary = true
281 );
282
283 //- Update the mesh based on the mesh files saved in time
284 // directories
285 virtual readUpdateState readUpdate();
286
287
288 // Access
289
290 //- Return the top-level database
291 const Time& time() const
292 {
293 return polyMesh::time();
294 }
295
296 //- Return true if thisDb() is a valid DB
297 virtual bool hasDb() const
298 {
299 return true;
300 }
301
302 //- Return the object registry - resolve conflict polyMesh/lduMesh
303 virtual const objectRegistry& thisDb() const
304 {
305 return polyMesh::thisDb();
306 }
307
308 //- Return reference to name
309 // Note: name() is currently ambiguous due to derivation from
310 // surfaceInterpolation
311 const word& name() const
312 {
313 return polyMesh::name();
314 }
315
316 //- Return reference to boundary mesh
317 const fvBoundaryMesh& boundary() const;
318
319 //- Return ldu addressing
320 virtual const lduAddressing& lduAddr() const;
321
322 //- Return a list of pointers for each patch
323 // with only those pointing to interfaces being set
324 virtual lduInterfacePtrsList interfaces() const;
325
326 //- Return communicator used for parallel communication
327 virtual label comm() const
328 {
329 return polyMesh::comm();
330 }
331
332
333 // Overlap
334
335 //- Interpolate interpolationCells only
336 virtual void interpolate(volScalarField&) const
337 {}
338
339 //- Interpolate interpolationCells only
340 virtual void interpolate(volVectorField&) const
341 {}
342
343 //- Interpolate interpolationCells only
344 virtual void interpolate(volSphericalTensorField&) const
345 {}
346
347 //- Interpolate interpolationCells only
348 virtual void interpolate(volSymmTensorField&) const
349 {}
350
351 //- Interpolate interpolationCells only
352 virtual void interpolate(volTensorField&) const
353 {}
354
355 //- Interpolate interpolationCells only. No bcs.
356 virtual void interpolate(scalarField&) const
357 {}
358
359 //- Interpolate interpolationCells only. No bcs.
360 virtual void interpolate(vectorField&) const
361 {}
362
363 //- Interpolate interpolationCells only. No bcs.
364 virtual void interpolate(sphericalTensorField&) const
365 {}
366
367 //- Interpolate interpolationCells only. No bcs.
368 virtual void interpolate(symmTensorField&) const
369 {}
370
371 //- Interpolate interpolationCells only. No bcs.
372 virtual void interpolate(tensorField&) const
373 {}
374
375 //- Solve returning the solution statistics given convergence
376 //- tolerance. Use the given solver controls
378 (
380 const dictionary&
381 ) const;
382
383 //- Solve returning the solution statistics given convergence
384 //- tolerance. Use the given solver controls
386 (
388 const dictionary&
389 ) const;
390
391 //- Solve returning the solution statistics given convergence
392 //- tolerance. Use the given solver controls
394 (
396 const dictionary&
397 ) const;
398
399 //- Solve returning the solution statistics given convergence
400 //- tolerance. Use the given solver controls
402 (
404 const dictionary&
405 ) const;
406
407 //- Solve returning the solution statistics given convergence
408 //- tolerance. Use the given solver controls
410 (
412 const dictionary&
413 ) const;
414
415
416 //- Internal face owner. Note bypassing virtual mechanism so
417 // e.g. relaxation always gets done using original addressing
418 const labelUList& owner() const
419 {
420 return fvMesh::lduAddr().lowerAddr();
421 }
422
423 //- Internal face neighbour
424 const labelUList& neighbour() const
425 {
426 return fvMesh::lduAddr().upperAddr();
427 }
428
429 //- Return cell volumes
431
432 //- Return old-time cell volumes
434
435 //- Return old-old-time cell volumes
437
438 //- Return sub-cycle cell volumes
440
441 //- Return sub-cycle old-time cell volumes
443
444 //- Return cell face area vectors
445 const surfaceVectorField& Sf() const;
446
447 //- Return cell face area magnitudes
448 const surfaceScalarField& magSf() const;
449
450 //- Return cell face motion fluxes
451 const surfaceScalarField& phi() const;
452
453 //- Return cell centres as volVectorField
454 const volVectorField& C() const;
455
456 //- Return face centres as surfaceVectorField
457 const surfaceVectorField& Cf() const;
458
459 //- Return face deltas as surfaceVectorField
461
462 //- Return a labelType of valid component indicators
463 // 1 : valid (solved)
464 // -1 : invalid (not solved)
465 template<class Type>
466 typename pTraits<Type>::labelType validComponents() const;
467
468
469 // Edit
470
471 //- Clear all geometry and addressing
472 void clearOut();
473
474 //- Update mesh corresponding to the given map
475 virtual void updateMesh(const mapPolyMesh& mpm);
476
477 //- Avoid masking surfaceInterpolation method
479
480 //- Move points, returns volumes swept by faces in motion
481 virtual void movePoints(const pointField&);
482
483 //- Update all geometric data. This gets redirected up from
484 //- primitiveMesh level
485 virtual void updateGeom();
486
487 //- Map all fields in time using given map.
488 virtual void mapFields(const mapPolyMesh& mpm);
489
490 //- Remove boundary patches. Warning: fvPatchFields hold ref to
491 //- these fvPatches.
492 void removeFvBoundary();
493
494 //- Return cell face motion fluxes (or null)
496
497 //- Return old-time cell volumes
499
500
501 // Write
502
503 //- Write the underlying polyMesh and other data
504 virtual bool writeObject
505 (
506 IOstreamOption streamOpt,
507 const bool valid
508 ) const;
509
510 //- Write mesh using IO settings from time
511 virtual bool write(const bool valid = true) const;
512
513
514 // Member Operators
515
516 //- Compares addresses
517 bool operator!=(const fvMesh& rhs) const;
518
519 //- Compares addresses
520 bool operator==(const fvMesh& rhs) const;
521};
522
523
524template<>
525typename pTraits<sphericalTensor>::labelType
526fvMesh::validComponents<sphericalTensor>() const;
527
528
529// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
530
531} // End namespace Foam
532
533// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
534
535#ifdef NoRepository
536 #include "fvMeshTemplates.C"
537 #include "fvPatchFvMeshTemplates.C"
538#endif
539
540// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
541
542#endif
543
544// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
The IOstreamOption is a simple container for options an IOstream can normally have.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Specialisation of DimensionedField that holds a slice of a given field so that it acts as a Dimension...
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
Database for solution data, solver performance and other reduced data.
Definition: data.H:58
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
scalar time
Injection time - set at collection [s].
Foam::fvBoundaryMesh.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Foam::fvMeshLduAddressing.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
virtual bool movePoints()
Avoid masking surfaceInterpolation method.
pTraits< Type >::labelType validComponents() const
Return a labelType of valid component indicators.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:302
virtual void interpolate(volVectorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:339
slicedSurfaceVectorField * SfPtr_
Face area vectors.
Definition: fvMesh.H:119
const volVectorField & C() const
Return cell centres as volVectorField.
virtual void interpolate(volSphericalTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:343
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
DimensionedField< scalar, volMesh > & setV0()
Return old-time cell volumes.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
DimensionedField< scalar, volMesh > * V00Ptr_
Cell volumes old-old time level.
Definition: fvMesh.H:116
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.C:735
void makeMagSf() const
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
void addFvPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:632
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
void makeC() const
void operator=(const fvMesh &)=delete
No copy assignment.
const surfaceScalarField & phi() const
Return cell face motion fluxes.
slicedSurfaceVectorField * CfPtr_
Face centres.
Definition: fvMesh.H:128
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:718
surfaceScalarField * phiPtr_
Face motion fluxes.
Definition: fvMesh.H:131
virtual void interpolate(volTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:351
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: fvMesh.C:276
fvBoundaryMesh boundary_
Boundary mesh.
Definition: fvMesh.H:97
void clearGeom()
Clear local geometry.
Definition: fvMesh.C:116
fvMesh Mesh
The mesh type.
Definition: fvMesh.H:179
SlicedDimensionedField< scalar, volMesh > * VPtr_
Cell volumes.
Definition: fvMesh.H:110
virtual void interpolate(volScalarField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:335
fvMesh(const fvMesh &)=delete
No copy construct.
const word & name() const
Return reference to name.
Definition: fvMesh.H:310
void clearOutLocal()
Clear local-only storage (geometry, addressing etc)
Definition: fvMesh.C:227
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:568
void makeSf() const
virtual void interpolate(vectorField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:359
void updateGeomNotOldVol()
Clear geometry like clearGeomNotOldVol but recreate any.
Definition: fvMesh.C:78
void clearGeomNotOldVol()
Clear geometry but not the old-time cell volumes.
Definition: fvMesh.C:54
virtual void interpolate(volSymmTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:347
void makeCf() const
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:417
DimensionedField< scalar, volMesh > * V0Ptr_
Cell volumes old time level.
Definition: fvMesh.H:113
void storeOldVol(const scalarField &)
Preserve old volume(s)
Definition: fvMesh.C:164
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
bool operator!=(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:1095
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:326
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:971
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:423
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:675
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:1056
const surfaceVectorField & Sf() const
Return cell face area vectors.
virtual void interpolate(sphericalTensorField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:363
virtual bool hasDb() const
Return true if thisDb() is a valid DB.
Definition: fvMesh.H:296
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
bool operator==(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:1101
fvMeshLduAddressing * lduPtr_
Definition: fvMesh.H:101
void removeFvBoundary()
Definition: fvMesh.C:662
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
fvBoundaryMesh BoundaryMesh
The boundary type associated with the mesh.
Definition: fvMesh.H:182
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: fvMesh.C:741
virtual void interpolate(symmTensorField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:367
ClassName("fvMesh")
label curTimeIndex_
Current time index for cell volumes.
Definition: fvMesh.H:107
virtual void interpolate(scalarField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:355
virtual void interpolate(tensorField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:371
slicedVolVectorField * CPtr_
Cell centres.
Definition: fvMesh.H:125
refPtr< surfaceScalarField > setPhi()
Return cell face motion fluxes (or null)
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
surfaceScalarField * magSfPtr_
Mag face area vectors.
Definition: fvMesh.H:122
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:239
virtual void updateGeom()
Definition: fvMesh.C:961
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:60
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:60
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:63
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
label comm() const noexcept
Return communicator used for parallel communication.
Definition: polyMesh.C:1328
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
const objectRegistry & thisDb() const noexcept
Return the object registry.
Definition: polyMesh.H:519
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
void clearAddressing()
Clear topological data.
const cellList & cells() const
A class for managing references or pointers (no reference counting)
Definition: refPtr.H:58
Cell to surface interpolation scheme. Included in fvMesh.
virtual bool movePoints()
Do what is necessary if the mesh has moved.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
Macro definitions for declaring ClassName(), NamespaceName(), etc.
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition: className.H:67
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
Namespace for OpenFOAM.
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
runTime write()
Forwards and collection of common point field types.
CEqn solve()