fvMesh.H
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-2022 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 Class
25  Foam::fvMesh
26 
27 Description
28  Mesh data needed to do the Finite Volume discretisation.
29 
30  NOTE ON USAGE:
31  fvMesh contains all the topological and geometric information
32  related to the mesh. It is also responsible for keeping the data
33  up-to-date. This is done by deleting the cell volume, face area,
34  cell/face centre, addressing and other derived information as
35  required and recalculating it as necessary. The fvMesh therefore
36  reserves the right to delete the derived information upon every
37  topological (mesh refinement/morphing) or geometric change (mesh
38  motion). It is therefore unsafe to keep local references to the
39  derived data outside of the time loop.
40 
41 SourceFiles
42  fvMesh.C
43  fvMeshGeometry.C
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #ifndef fvMesh_H
48 #define fvMesh_H
49 
50 #include "polyMesh.H"
51 #include "lduMesh.H"
52 #include "primitiveMesh.H"
53 #include "fvBoundaryMesh.H"
54 #include "surfaceInterpolation.H"
55 #include "fvSchemes.H"
56 #include "fvSolution.H"
57 #include "data.H"
58 #include "DimensionedField.H"
59 #include "SlicedDimensionedField.H"
60 #include "volFieldsFwd.H"
61 #include "surfaceFieldsFwd.H"
62 #include "pointFieldsFwd.H"
63 #include "slicedVolFieldsFwd.H"
64 #include "slicedSurfaceFieldsFwd.H"
65 #include "className.H"
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 namespace Foam
70 {
71 
72 class fvMeshLduAddressing;
73 class fvMeshStitcher;
74 class fvMeshTopoChanger;
75 class fvMeshDistributor;
76 class fvMeshMover;
77 class volMesh;
78 class polyDistributionMap;
79 
80 template<class Type, template<class> class PatchField, class GeoMesh>
82 
83 template<class Type>
84 class UCompactListList;
85 
86 namespace fvMeshTopoChangers
87 {
88  class list;
89 }
90 
91 
92 /*---------------------------------------------------------------------------*\
93  Class fvMesh Declaration
94 \*---------------------------------------------------------------------------*/
95 
96 class fvMesh
97 :
98  public polyMesh,
99  public lduMesh,
100  public surfaceInterpolation,
101  public data
102 {
103 public:
104 
105  // Public data types
106 
107  //- Extent to which to stitch on read and readUpdate. By default, a
108  // full geometric stitch is performed. A non-geometric stitch can be
109  // done as an optimisation in situations when finite volume geometry
110  // is not needed (e.g., decomposition). Stitching can also be
111  // prevented altogether if that is appropriate (e.g., if the mesh is
112  // loaded for mapping in an un-stitched state).
113  enum class stitchType
114  {
115  none,
116  nonGeometric,
117  geometric
118  };
119 
120 
121 private:
122 
123  // Private Data
124 
125  //- Boundary mesh
126  fvBoundaryMesh boundary_;
127 
128  //- The stitcher function class
129  autoPtr<fvMeshStitcher> stitcher_;
130 
131  //- The topo-changer function class
132  autoPtr<fvMeshTopoChanger> topoChanger_;
133 
134  //- The distributor function class
135  autoPtr<fvMeshDistributor> distributor_;
136 
137  //- The mover function class
138  autoPtr<fvMeshMover> mover_;
139 
140 
141  // Demand-driven data
142 
143  //- Matrix addressing
144  mutable fvMeshLduAddressing* lduPtr_;
145 
146  //- Face-poly-face addressing
148  polyFacesBfPtr_;
149 
150  //- Offsets for poly-bFace-patch and patch-face addressing
151  mutable labelList* polyBFaceOffsetsPtr_;
152 
153  //- Poly-bFace-patch addressing
154  mutable labelList* polyBFaceOffsetPatchesPtr_;
155 
156  //- Poly-bFace-patch-face addressing
157  mutable labelList* polyBFaceOffsetPatchFacesPtr_;
158 
159  //- Poly-bFace-patch addressing
160  mutable UCompactListList<label>* polyBFacePatchesPtr_;
161 
162  //- Poly-bFace-patch-face addressing
163  mutable UCompactListList<label>* polyBFacePatchFacesPtr_;
164 
165  //- Current time index for cell volumes
166  // Note. The whole mechanism will be replaced once the
167  // dimensionedField is created and the dimensionedField
168  // will take care of the old-time levels.
169  mutable label curTimeIndex_;
170 
171  //- Cell volumes old time level
173 
174  //- Cell volumes old time level
175  mutable DimensionedField<scalar, volMesh>* V0Ptr_;
176 
177  //- Cell volumes old-old time level
178  mutable DimensionedField<scalar, volMesh>* V00Ptr_;
179 
180  //- Face area vectors
181  mutable slicedSurfaceVectorField* SfSlicePtr_;
182 
183  //- Face area vectors
184  mutable surfaceVectorField* SfPtr_;
185 
186  //- Mag face area vectors
187  mutable slicedSurfaceScalarField* magSfSlicePtr_;
188 
189  //- Mag face area vectors
190  mutable surfaceScalarField* magSfPtr_;
191 
192  //- Cell centres
193  mutable slicedVolVectorField* CSlicePtr_;
194 
195  //- Cell centres
196  mutable volVectorField* CPtr_;
197 
198  //- Face centres
199  mutable slicedSurfaceVectorField* CfSlicePtr_;
200 
201  //- Face centres
202  mutable surfaceVectorField* CfPtr_;
203 
204  //- Face motion fluxes
205  mutable surfaceScalarField* phiPtr_;
206 
207  //- fvSchemes created on demand
208  mutable autoPtr<fvSchemes> fvSchemes_;
209 
210  //- fvSolution created on demand
211  mutable autoPtr<fvSolution> fvSolution_;
212 
213 
214  // Private Member Functions
215 
216  // Storage management
217 
218  //- Clear geometry but not the old-time cell volumes
219  void clearGeomNotOldVol();
220 
221  //- Clear geometry like clearGeomNotOldVol but recreate any
222  // geometric demand-driven data that was set
223  void updateGeomNotOldVol();
224 
225  //- Clear geometry
226  void clearGeom();
227 
228  //- Clear addressing
229  void clearAddressing(const bool isMeshUpdate = false);
230 
231  //- Preserve old volume(s)
232  void storeOldVol(const scalarField&);
233 
234 
235  // Make topological data
236 
237  //- Get the patch types for the poly faces field
238  wordList polyFacesPatchTypes() const;
239 
240  //- Modify face-poly-face addressing
242  polyFacesBfRef();
243 
244 
245  // Make geometric data
246 
247  //- Make the face area vectors
248  void makeSf() const;
249 
250  //- Make the mag face area vectors
251  void makeMagSf() const;
252 
253  //- Make the cell centres
254  void makeC() const;
255 
256  //- Make the face centres
257  void makeCf() const;
258 
259  //- Modify old-time cell volumes
261 
262  //- Modify old-old-time cell volumes
264 
265  //- Modify cell face area vectors
266  surfaceVectorField& SfRef();
267 
268  //- Modify cell face area magnitudes
269  surfaceScalarField& magSfRef();
270 
271  //- Modify cell centres
272  volVectorField& CRef();
273 
274  //- Modify face centres
275  surfaceVectorField& CfRef();
276 
277  //- Modify cell face motion fluxes
278  surfaceScalarField& phiRef();
279 
280 
281 public:
282 
283  // Public Typedefs
285  typedef fvMesh Mesh;
287 
289 
290  // Declare name of the class and its debug switch
291  ClassName("fvMesh");
292 
293 
294  // Constructors
295 
296  //- Construct from IOobject
297  // with the option to not instantiate the mesh changers and the option
298  // to prevent some or all of the stitching
299  explicit fvMesh
300  (
301  const IOobject& io,
302  const bool changers = true,
303  const stitchType stitch = stitchType::geometric
304  );
305 
306  //- Construct from cellShapes with boundary.
307  fvMesh
308  (
309  const IOobject& io,
310  pointField&& points,
311  const cellShapeList& shapes,
312  const faceListList& boundaryFaces,
313  const wordList& boundaryPatchNames,
314  const PtrList<dictionary>& boundaryDicts,
315  const word& defaultBoundaryPatchName,
316  const word& defaultBoundaryPatchType,
317  const bool syncPar = true
318  );
319 
320  //- Construct from components without boundary.
321  // Boundary is added using addFvPatches() member function
322  fvMesh
323  (
324  const IOobject& io,
325  pointField&& points,
326  faceList&& faces,
327  labelList&& allOwner,
328  labelList&& allNeighbour,
329  const bool syncPar = true
330  );
331 
332  //- Construct without boundary from cells rather than owner/neighbour.
333  // Boundary is added using addPatches() member function
334  fvMesh
335  (
336  const IOobject& io,
337  pointField&& points,
338  faceList&& faces,
339  cellList&& cells,
340  const bool syncPar = true
341  );
342 
343  //- Disallow default bitwise copy construction
344  fvMesh(const fvMesh&);
345 
346 
347  //- Destructor
348  virtual ~fvMesh();
349 
350 
351  // Member Functions
352 
353  // Helpers
354 
355  //- Add boundary patches. Constructor helper
356  void addFvPatches
357  (
358  const List<polyPatch*>&,
359  const bool validBoundary = true
360  );
361 
362  //- Update the mesh based on the mesh files saved in time
363  // directories
364  readUpdateState readUpdate
365  (
366  const stitchType stitch = stitchType::geometric
367  );
368 
369 
370  // Access
371 
372  //- Return the top-level database
373  const Time& time() const
374  {
375  return polyMesh::time();
376  }
377 
378  //- Return the object registry - resolve conflict polyMesh/lduMesh
379  virtual const objectRegistry& thisDb() const
380  {
381  return polyMesh::thisDb();
382  }
383 
384  //- Return reference to name
385  // Note: name() is currently ambiguous due to derivation from
386  // surfaceInterpolation
387  const word& name() const
388  {
389  return polyMesh::name();
390  }
391 
392  //- Return reference to boundary mesh
393  const fvBoundaryMesh& boundary() const;
394 
395  //- Return reference to polyMesh
396  const polyMesh& operator()() const
397  {
398  return *this;
399  }
400 
401  //- Return ldu addressing
402  virtual const lduAddressing& lduAddr() const;
403 
404  //- Return a list of pointers for each patch
405  // with only those pointing to interfaces being set
406  virtual lduInterfacePtrsList interfaces() const
407  {
408  return boundary().interfaces();
409  }
410 
411  //- Return communicator used for parallel communication
412  virtual label comm() const
413  {
414  return polyMesh::comm();
415  }
416 
417  //- Internal face owner
418  const labelUList& owner() const
419  {
420  return lduAddr().lowerAddr();
421  }
422 
423  //- Internal face neighbour
424  const labelUList& neighbour() const
425  {
426  return lduAddr().upperAddr();
427  }
428 
429  //- Return whether the fvMesh is conformal with the polyMesh
430  bool conformal() const;
431 
432  //- Get the poly faces IO object
433  IOobject polyFacesBfIO(IOobject::readOption r) const;
434 
435  //- Return face-poly-face addressing
437  polyFacesBf() const;
438 
439  //- Return poly-bFace-patch addressing
440  const UCompactListList<label>& polyBFacePatches() const;
441 
442  //- Return poly-bFace-patch-face addressing
443  const UCompactListList<label>& polyBFacePatchFaces() const;
444 
445  //- Return the stitcher function class
446  const fvMeshStitcher& stitcher() const;
447 
448  //- Return the topo-changer function class
449  const fvMeshTopoChanger& topoChanger() const;
450 
451  //- Return the distributor function class
452  const fvMeshDistributor& distributor() const;
453 
454  //- Return the mover function class
455  const fvMeshMover& mover() const;
456 
457  //- Return cell volumes
458  const DimensionedField<scalar, volMesh>& V() const;
459 
460  //- Return old-time cell volumes
461  const DimensionedField<scalar, volMesh>& V0() const;
462 
463  //- Return old-old-time cell volumes
464  const DimensionedField<scalar, volMesh>& V00() const;
465 
466  //- Return sub-cycle cell volumes
468 
469  //- Return sub-cycle old-time cell volumes
471 
472  //- Return cell face area vectors
473  const surfaceVectorField& Sf() const;
474 
475  //- Return cell face area magnitudes
476  const surfaceScalarField& magSf() const;
477 
478  //- Return cell centres
479  const volVectorField& C() const;
480 
481  //- Return face centres
482  const surfaceVectorField& Cf() const;
483 
484  //- Return face deltas as surfaceVectorField
486 
487  //- Return cell face motion fluxes
488  const surfaceScalarField& phi() const;
489 
490  //- Return a labelType of valid component indicators
491  // 1 : valid (solved)
492  // -1 : invalid (not solved)
493  template<class Type>
494  typename pTraits<Type>::labelType validComponents() const;
495 
496  //- Return the fvSchemes
497  const fvSchemes& schemes() const;
498 
499  //- Return the fvSchemes
500  const fvSolution& solution() const;
501 
502 
503  // Edit
504 
505  //- Is this mesh dynamic?
506  bool dynamic() const;
507 
508  //- Update the mesh for topology change, mesh to mesh mapping
509  // or redistribution
510  bool update();
511 
512  //- Move the mesh
513  bool move();
514 
515  //- Clear all geometry and addressing
516  void clearOut();
517 
518  //- Update mesh corresponding to the given map
519  virtual void topoChange(const polyTopoChangeMap& map);
520 
521  //- Update from another mesh using the given map
522  virtual void mapMesh(const polyMeshMap&);
523 
524  //- Redistribute or update using the given distribution map
525  virtual void distribute(const polyDistributionMap& map);
526 
527  //- Reset the points
528  // without storing old points or returning swept volumes
529  virtual void setPoints(const pointField&);
530 
531  //- Move points, returns volumes swept by faces in motion
532  virtual tmp<scalarField> movePoints(const pointField&);
533 
534  //- Conform the fvMesh to the polyMesh
535  void conform
536  (
537  const surfaceScalarField& phi =
538  NullObjectRef<surfaceScalarField>()
539  );
540 
541  //- Unconform the fvMesh from the polyMesh
542  void unconform
543  (
545  polyFacesBf,
546  const surfaceVectorField& Sf,
547  const surfaceVectorField& Cf,
548  const surfaceScalarField& phi =
549  NullObjectRef<surfaceScalarField>(),
550  const bool sync = true
551  );
552 
553  //- Map all fields in time using given map.
554  void mapFields(const polyTopoChangeMap& map);
555 
556  //- Add/insert single patch. If validBoundary the new situation
557  // is consistent across processors.
558  virtual void addPatch
559  (
560  const label insertPatchi,
561  const polyPatch& patch,
562  const dictionary& patchFieldDict,
563  const word& defaultPatchFieldType,
564  const bool validBoundary
565  );
566 
567  //- Reorder and trim existing patches. If validBoundary the new
568  // situation is consistent across processors
569  virtual void reorderPatches
570  (
571  const labelUList& newToOld,
572  const bool validBoundary
573  );
574 
575  //- Remove boundary patches. Warning: fvPatchFields hold ref to
576  // these fvPatches.
577  void removeFvBoundary();
578 
579  //- Reset mesh
580  // For run-time mesh replacement and mesh to mesh mapping
581  void reset(const fvMesh&);
582 
583 
584  // Write
585 
586  //- Write the underlying polyMesh and other data
587  virtual bool writeObject
588  (
592  const bool write = true
593  ) const;
594 
595  //- Write mesh using IO settings from time
596  virtual bool write(const bool write = true) const;
597 
598 
599  // Member Operators
600 
601  //- Disallow default bitwise assignment
602  void operator=(const fvMesh&);
603 
604  bool operator!=(const fvMesh&) const;
605  bool operator==(const fvMesh&) const;
606 };
607 
608 
609 template<>
611 fvMesh::validComponents<sphericalTensor>() const;
612 
613 
614 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
615 
616 } // End namespace Foam
617 
618 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
619 
620 #ifdef NoRepository
621  #include "fvMeshTemplates.C"
622 #endif
623 
624 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
625 
626 #endif
627 
628 // ************************************************************************* //
scalar delta
Graphite solid properties.
Definition: C.H:48
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
const word & name() const
Return name.
Definition: IOobject.H:315
Specialisation of DimensionedField which holds a slice of a given complete field in such a form that ...
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Sequence of mesh topology changes applied in order.
Mesh manipulator that uses the intersection provided by the cyclic non-conformal poly patches to crea...
Abstract base class for fvMesh movers.
Definition: fvMeshMover.H:52
Traits class for primitives.
Definition: pTraits.H:50
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:59
Generic GeometricField class.
Cell to surface interpolation scheme. Included in fvMesh.
readOption
Enumeration defining the read options.
Definition: IOobject.H:116
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition: className.H:65
mesh schemes().setFluxRequired(pcorr.name())
stitchType
Extent to which to stitch on read and readUpdate. By default, a.
Definition: fvMesh.H:112
const cellShapeList & cells
const pointField & points
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
void mapMesh(const meshToMesh &interp, const HashSet< word > &selectedFields, const bool noLagrangian)
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
phi
Definition: correctPhi.H:3
faceListList boundary(nPatches)
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1453
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
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
Foam::fvMeshLduAddressing.
const Time & time() const
Return time.
Database for solution and other reduced data.
Definition: data.H:51
Abstract base class for fvMesh movers.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Generic GeometricBoundaryField class.
Definition: fvMesh.H:80
Unallocated base class of CompactListList.
Definition: fvMesh.H:83
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:47
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Foam::fvBoundaryMesh.
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:50
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Version number type.
Definition: IOstream.H:96
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:516
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
Specialisation of GeometricField which holds slices of given complete fields in a form that they act ...
bool operator!=(const particle &, const particle &)
Definition: particle.C:1282
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Abstract base class for fvMesh movers.
Namespace for OpenFOAM.