dynamicRefineFvMesh.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 OpenFOAM Foundation
9  Copyright (C) 2017-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 Class
28  Foam::dynamicRefineFvMesh
29 
30 Description
31  A fvMesh with built-in refinement.
32 
33  Determines which cells to refine/unrefine and does all in update().
34 
35  \verbatim
36  {
37  // How often to refine
38  refineInterval 1;
39  // Field to be refinement on
40  field alpha.water;
41  // Refine field inbetween lower..upper
42  lowerRefineLevel 0.001;
43  upperRefineLevel 0.999;
44  // If value < unrefineLevel (default=GREAT) unrefine
45  //unrefineLevel 10;
46  // Have slower than 2:1 refinement
47  nBufferLayers 1;
48  // Refine cells only up to maxRefinement levels
49  maxRefinement 2;
50  // Stop refinement if maxCells reached
51  maxCells 200000;
52  // Flux field and corresponding velocity field. Fluxes on changed
53  // faces get recalculated by interpolating the velocity. Use 'none'
54  // on surfaceScalarFields that do not need to be reinterpolated, use
55  // NaN to detect use of mapped variable
56  correctFluxes
57  (
58  (phi none) //NaN) //none)
59  (nHatf none) //none)
60  (rho*phi none) //none)
61  (ghf none) //NaN) //none)
62  );
63 
64  // Write the refinement level as a volScalarField
65  dumpLevel true;
66  }
67  \endverbatim
68 
69 
70 SourceFiles
71  dynamicRefineFvMesh.C
72 
73 \*---------------------------------------------------------------------------*/
74 
75 #ifndef dynamicRefineFvMesh_H
76 #define dynamicRefineFvMesh_H
77 
78 #include "dynamicFvMesh.H"
79 #include "hexRef8.H"
80 #include "bitSet.H"
81 
82 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83 
84 namespace Foam
85 {
86 
87 /*---------------------------------------------------------------------------*\
88  Class dynamicRefineFvMesh Declaration
89 \*---------------------------------------------------------------------------*/
90 
92 :
93  public dynamicFvMesh
94 {
95 protected:
96 
97  //- Mesh cutting engine
99 
100  //- Fluxes to map
102 
103  //- Protected cells (usually since not hexes)
105 
106  //- Number of refinement/unrefinement steps done so far.
108 
109  //- Dump cellLevel for post-processing
110  bool dumpLevel_;
111 
112 
113  // Protected Member Functions
114 
115  //- Calculate cells that cannot be refined since would trigger
116  // refinement of protectedCell_ (since 2:1 refinement cascade)
117  void calculateProtectedCells(bitSet& unrefineableCell) const;
118 
119  //- Read the projection parameters from dictionary
120  void readDict();
121 
122 
123  //- Refine cells. Update mesh and fields.
124  virtual autoPtr<mapPolyMesh> refine(const labelList&);
125 
126  //- Unrefine cells. Gets passed in centre points of cells to combine.
127  virtual autoPtr<mapPolyMesh> unrefine(const labelList&);
128 
129 
130  // Selection of cells to un/refine
131 
132  //- Calculates approximate value for refinement level so
133  // we don't go above maxCell
134  scalar getRefineLevel
135  (
136  const label maxCells,
137  const label maxRefinement,
138  const scalar refineLevel,
139  const scalarField&
140  ) const;
141 
142  //- Get per cell max of connected point
143  scalarField maxPointField(const scalarField&) const;
144 
145  //- Get point max of connected cell
147 
148  scalarField cellToPoint(const scalarField& vFld) const;
149 
151  (
152  const scalarField& fld,
153  const scalar minLevel,
154  const scalar maxLevel
155  ) const;
156 
157  //- Select candidate cells for refinement
158  virtual void selectRefineCandidates
159  (
160  const scalar lowerRefineLevel,
161  const scalar upperRefineLevel,
162  const scalarField& vFld,
163  bitSet& candidateCell
164  ) const;
165 
166  //- Subset candidate cells for refinement
168  (
169  const label maxCells,
170  const label maxRefinement,
171  const bitSet& candidateCell
172  ) const;
173 
174  //- Select points that can be unrefined.
176  (
177  const scalar unrefineLevel,
178  const bitSet& markedCell,
179  const scalarField& pFld
180  ) const;
181 
182  //- Extend markedCell with cell-face-cell.
183  void extendMarkedCells(bitSet& markedCell) const;
184 
185  //- Check all cells have 8 anchor points
187 
188  //- Map single non-flux surface<Type>Field
189  // for new internal faces (e.g. AMR refine). This currently
190  // interpolates values from surrounding faces (faces on
191  // neighbouring cells) that do have a value.
192  template <class T>
194  (
195  const labelList& faceMap,
197  );
198 
199  //- Correct surface fields for new faces
200  template <class T>
202 
203  //- Correct surface fields for new faces. Converts any oriented
204  // fields into non-oriented (i.e. phi -> Uf) before mapping
205  template <class T>
207  (
208  const surfaceVectorField& Sf,
209  const surfaceScalarField& magSf,
210  const labelList& faceMap
211  );
212 
213 
214 private:
215 
216  //- No copy construct
217  dynamicRefineFvMesh(const dynamicRefineFvMesh&) = delete;
218 
219  //- No copy assignment
220  void operator=(const dynamicRefineFvMesh&) = delete;
221 
222 public:
223 
224  //- Runtime type information
225  TypeName("dynamicRefineFvMesh");
226 
227 
228  // Constructors
229 
230  //- Construct from IOobject
231  explicit dynamicRefineFvMesh
232  (
233  const IOobject& io,
234  const bool doInit=true
235  );
236 
237 
238  //- Destructor
239  virtual ~dynamicRefineFvMesh() = default;
240 
241 
242  // Member Functions
243 
244  //- Initialise all non-demand-driven data
245  virtual bool init(const bool doInit);
246 
247  //- Direct access to the refinement engine
248  const hexRef8& meshCutter() const
249  {
250  return meshCutter_;
251  }
252 
253  //- Cells which should not be refined/unrefined
254  const bitSet& protectedCell() const
255  {
256  return protectedCell_;
257  }
258 
259  //- Cells which should not be refined/unrefined
261  {
262  return protectedCell_;
263  }
264 
265  //- Update the mesh for both mesh motion and topology change
266  virtual bool update();
267 
268  //- Map all fields in time using given map.
269  virtual void mapFields(const mapPolyMesh& mpm);
270 
271 
272  // Writing
273 
274  //- Write using stream options
275  virtual bool writeObject
276  (
277  IOstreamOption streamOpt,
278  const bool valid
279  ) const;
280 };
281 
282 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
283 
284 } // End namespace Foam
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 #ifdef NoRepository
290 #endif
291 
292 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 
294 #endif
295 
296 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::dynamicRefineFvMesh::error
scalarField error(const scalarField &fld, const scalar minLevel, const scalar maxLevel) const
Definition: dynamicRefineFvMesh.C:715
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::dynamicRefineFvMesh::extendMarkedCells
void extendMarkedCells(bitSet &markedCell) const
Extend markedCell with cell-face-cell.
Definition: dynamicRefineFvMesh.C:946
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::dynamicRefineFvMesh::protectedCell_
bitSet protectedCell_
Protected cells (usually since not hexes)
Definition: dynamicRefineFvMesh.H:103
Foam::dynamicRefineFvMesh::TypeName
TypeName("dynamicRefineFvMesh")
Runtime type information.
Foam::dynamicRefineFvMesh::init
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: dynamicRefineFvMesh.C:1040
Foam::dynamicFvMesh
Abstract base class for geometry and/or topology changing fvMesh.
Definition: dynamicFvMesh.H:78
Foam::dynamicRefineFvMesh::meshCutter_
hexRef8 meshCutter_
Mesh cutting engine.
Definition: dynamicRefineFvMesh.H:97
Foam::dynamicRefineFvMesh::mapFields
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: dynamicRefineFvMesh.C:197
Foam::dynamicRefineFvMesh::update
virtual bool update()
Update the mesh for both mesh motion and topology change.
Definition: dynamicRefineFvMesh.C:1211
Foam::dynamicRefineFvMesh::mapNewInternalFaces
void mapNewInternalFaces(const labelList &faceMap, GeometricField< T, fvsPatchField, surfaceMesh > &)
Map single non-flux surface<Type>Field.
Definition: dynamicRefineFvMeshTemplates.C:34
Foam::dynamicRefineFvMesh::selectRefineCells
virtual labelList selectRefineCells(const label maxCells, const label maxRefinement, const bitSet &candidateCell) const
Subset candidate cells for refinement.
Definition: dynamicRefineFvMesh.C:771
Foam::dynamicRefineFvMesh::correctFluxes_
HashTable< word > correctFluxes_
Fluxes to map.
Definition: dynamicRefineFvMesh.H:100
Foam::dynamicRefineFvMesh::getRefineLevel
scalar getRefineLevel(const label maxCells, const label maxRefinement, const scalar refineLevel, const scalarField &) const
Calculates approximate value for refinement level so.
bitSet.H
Foam::dynamicRefineFvMesh::readDict
void readDict()
Read the projection parameters from dictionary.
Definition: dynamicRefineFvMesh.C:166
dynamicRefineFvMeshTemplates.C
Foam::dynamicRefineFvMesh::writeObject
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
Definition: dynamicRefineFvMesh.C:1417
Foam::dynamicRefineFvMesh::checkEightAnchorPoints
void checkEightAnchorPoints(bitSet &protectedCell) const
Check all cells have 8 anchor points.
Definition: dynamicRefineFvMesh.C:980
Foam::Field< scalar >
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:330
hexRef8.H
Foam::dynamicRefineFvMesh::selectRefineCandidates
virtual void selectRefineCandidates(const scalar lowerRefineLevel, const scalar upperRefineLevel, const scalarField &vFld, bitSet &candidateCell) const
Select candidate cells for refinement.
Definition: dynamicRefineFvMesh.C:737
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:63
Foam::dynamicRefineFvMesh::maxCellField
scalarField maxCellField(const volScalarField &) const
Get point max of connected cell.
Definition: dynamicRefineFvMesh.C:677
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::dynamicRefineFvMesh::refine
virtual autoPtr< mapPolyMesh > refine(const labelList &)
Refine cells. Update mesh and fields.
Definition: dynamicRefineFvMesh.C:393
Foam::dynamicRefineFvMesh::protectedCell
const bitSet & protectedCell() const
Cells which should not be refined/unrefined.
Definition: dynamicRefineFvMesh.H:253
Foam::dynamicRefineFvMesh::protectedCell
bitSet & protectedCell()
Cells which should not be refined/unrefined.
Definition: dynamicRefineFvMesh.H:259
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dynamicRefineFvMesh::selectUnrefinePoints
virtual labelList selectUnrefinePoints(const scalar unrefineLevel, const bitSet &markedCell, const scalarField &pFld) const
Select points that can be unrefined.
Definition: dynamicRefineFvMesh.C:851
Foam::dynamicRefineFvMesh::~dynamicRefineFvMesh
virtual ~dynamicRefineFvMesh()=default
Destructor.
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:105
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::dynamicRefineFvMesh::calculateProtectedCells
void calculateProtectedCells(bitSet &unrefineableCell) const
Calculate cells that cannot be refined since would trigger.
Definition: dynamicRefineFvMesh.C:53
Foam::hexRef8
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:67
Foam::List< label >
Foam::dynamicRefineFvMesh
A fvMesh with built-in refinement.
Definition: dynamicRefineFvMesh.H:90
Foam::dynamicRefineFvMesh::dumpLevel_
bool dumpLevel_
Dump cellLevel for post-processing.
Definition: dynamicRefineFvMesh.H:109
dynamicFvMesh.H
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::dynamicRefineFvMesh::unrefine
virtual autoPtr< mapPolyMesh > unrefine(const labelList &)
Unrefine cells. Gets passed in centre points of cells to combine.
Definition: dynamicRefineFvMesh.C:482
Foam::dynamicRefineFvMesh::cellToPoint
scalarField cellToPoint(const scalarField &vFld) const
Definition: dynamicRefineFvMesh.C:695
Foam::dynamicRefineFvMesh::meshCutter
const hexRef8 & meshCutter() const
Direct access to the refinement engine.
Definition: dynamicRefineFvMesh.H:247
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::dynamicRefineFvMesh::maxPointField
scalarField maxPointField(const scalarField &) const
Get per cell max of connected point.
Definition: dynamicRefineFvMesh.C:659
Foam::dynamicRefineFvMesh::nRefinementIterations_
label nRefinementIterations_
Number of refinement/unrefinement steps done so far.
Definition: dynamicRefineFvMesh.H:106
Foam::fvMesh::Sf
const surfaceVectorField & Sf() const
Return cell face area vectors.
Definition: fvMeshGeometry.C:319