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"
80 #include "hexRef8.H"
81 #include "bitSet.H"
82 
83 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85 namespace Foam
86 {
87 
88 /*---------------------------------------------------------------------------*\
89  Class dynamicRefineFvMesh Declaration
90 \*---------------------------------------------------------------------------*/
91 
93 :
94  //public dynamicFvMesh
96 {
97 protected:
98 
99  //- Mesh cutting engine
101 
102  //- Fluxes to map
104 
105  //- Protected cells (usually since not hexes)
108  //- Number of refinement/unrefinement steps done so far.
110 
111  //- Dump cellLevel for post-processing
113 
114 
115  // Protected Member Functions
116 
117  //- Calculate cells that cannot be refined since would trigger
118  // refinement of protectedCell_ (since 2:1 refinement cascade)
119  void calculateProtectedCells(bitSet& unrefineableCell) const;
120 
121  //- Read the projection parameters from dictionary
122  void readDict();
123 
124 
125  //- Refine cells. Update mesh and fields.
126  virtual autoPtr<mapPolyMesh> refine(const labelList&);
127 
128  //- Unrefine cells. Gets passed in centre points of cells to combine.
129  virtual autoPtr<mapPolyMesh> unrefine(const labelList&);
130 
131 
132  // Selection of cells to un/refine
133 
134  //- Calculates approximate value for refinement level so
135  // we don't go above maxCell
136  scalar getRefineLevel
137  (
138  const label maxCells,
139  const label maxRefinement,
140  const scalar refineLevel,
141  const scalarField&
142  ) const;
143 
144  //- Get per cell max of connected point
145  scalarField maxPointField(const scalarField&) const;
146 
147  //- Get point max of connected cell
149 
150  scalarField cellToPoint(const scalarField& vFld) const;
151 
153  (
154  const scalarField& fld,
155  const scalar minLevel,
156  const scalar maxLevel
157  ) const;
158 
159  //- Select candidate cells for refinement
160  virtual void selectRefineCandidates
161  (
162  const scalar lowerRefineLevel,
163  const scalar upperRefineLevel,
164  const scalarField& vFld,
165  bitSet& candidateCell
166  ) const;
167 
168  //- Subset candidate cells for refinement
170  (
171  const label maxCells,
172  const label maxRefinement,
173  const bitSet& candidateCell
174  ) const;
175 
176  //- Select points that can be unrefined.
178  (
179  const scalar unrefineLevel,
180  const bitSet& markedCell,
181  const scalarField& pFld
182  ) const;
183 
184  //- Extend markedCell with cell-face-cell.
185  void extendMarkedCells(bitSet& markedCell) const;
186 
187  //- Check all cells have 8 anchor points
189 
190  //- Map single non-flux surface<Type>Field
191  // for new internal faces (e.g. AMR refine). This currently
192  // interpolates values from surrounding faces (faces on
193  // neighbouring cells) that do have a value.
194  template <class T>
196  (
197  const labelList& faceMap,
199  );
200 
201  //- Correct surface fields for new faces
202  template <class T>
204 
205  //- Correct surface fields for new faces. Converts any oriented
206  // fields into non-oriented (i.e. phi -> Uf) before mapping
207  template <class T>
209  (
210  const surfaceVectorField& Sf,
211  const surfaceScalarField& magSf,
212  const labelList& faceMap
213  );
214 
215  //- Update topology (refinement, unrefinement)
216  bool updateTopology();
217 
218 
219 private:
220 
221  //- No copy construct
222  dynamicRefineFvMesh(const dynamicRefineFvMesh&) = delete;
223 
224  //- No copy assignment
225  void operator=(const dynamicRefineFvMesh&) = delete;
226 
227 public:
228 
229  //- Runtime type information
230  TypeName("dynamicRefineFvMesh");
231 
232 
233  // Constructors
234 
235  //- Construct from IOobject
236  explicit dynamicRefineFvMesh
237  (
238  const IOobject& io,
239  const bool doInit=true
240  );
241 
242 
243  //- Destructor
244  virtual ~dynamicRefineFvMesh() = default;
245 
246 
247  // Member Functions
248 
249  //- Initialise all non-demand-driven data
250  virtual bool init(const bool doInit);
251 
252  //- Direct access to the refinement engine
253  const hexRef8& meshCutter() const
254  {
255  return meshCutter_;
256  }
257 
258  //- Cells which should not be refined/unrefined
259  const bitSet& protectedCell() const
260  {
261  return protectedCell_;
262  }
263 
264  //- Cells which should not be refined/unrefined
265  bitSet& protectedCell()
266  {
267  return protectedCell_;
268  }
269 
270  //- Update the mesh for both mesh motion and topology change
271  virtual bool update();
272 
273  //- Map all fields in time using given map.
274  virtual void mapFields(const mapPolyMesh& mpm);
275 
276 
277  // Writing
278 
279  //- Write using stream options
280  virtual bool writeObject
281  (
282  IOstreamOption streamOpt,
283  const bool writeOnProc
284  ) const;
285 };
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace Foam
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 #ifdef NoRepository
295 #endif
296 
297 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
298 
299 #endif
300 
301 // ************************************************************************* //
void calculateProtectedCells(bitSet &unrefineableCell) const
Calculate cells that cannot be refined since would trigger.
const surfaceVectorField & Sf() const
Return cell face area vectors.
bool dumpLevel_
Dump cellLevel for post-processing.
bool updateTopology()
Update topology (refinement, unrefinement)
scalarField maxPointField(const scalarField &) const
Get per cell max of connected point.
const hexRef8 & meshCutter() const
Direct access to the refinement engine.
scalar getRefineLevel(const label maxCells, const label maxRefinement, const scalar refineLevel, const scalarField &) const
Calculates approximate value for refinement level so.
void checkEightAnchorPoints(bitSet &protectedCell) const
Check all cells have 8 anchor points.
scalarField cellToPoint(const scalarField &vFld) const
bitSet protectedCell_
Protected cells (usually since not hexes)
HashTable< word > correctFluxes_
Fluxes to map.
TypeName("dynamicRefineFvMesh")
Runtime type information.
A fvMesh with built-in refinement.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
virtual labelList selectUnrefinePoints(const scalar unrefineLevel, const bitSet &markedCell, const scalarField &pFld) const
Select points that can be unrefined.
virtual ~dynamicRefineFvMesh()=default
Destructor.
virtual bool update()
Update the mesh for both mesh motion and topology change.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write using stream options.
void mapNewInternalFaces(const labelList &faceMap, GeometricField< T, fvsPatchField, surfaceMesh > &)
Map single non-flux surface<Type>Field.
virtual autoPtr< mapPolyMesh > unrefine(const labelList &)
Unrefine cells. Gets passed in centre points of cells to combine.
virtual void selectRefineCandidates(const scalar lowerRefineLevel, const scalar upperRefineLevel, const scalarField &vFld, bitSet &candidateCell) const
Select candidate cells for refinement.
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:62
hexRef8 meshCutter_
Mesh cutting engine.
virtual labelList selectRefineCells(const label maxCells, const label maxRefinement, const bitSet &candidateCell) const
Subset candidate cells for refinement.
virtual autoPtr< mapPolyMesh > refine(const labelList &)
Refine cells. Update mesh and fields.
A HashTable similar to std::unordered_map.
Definition: HashTable.H:108
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Dynamic mesh able to handle multiple motion solvers. NOTE: If the word entry "solvers" is not found i...
label nRefinementIterations_
Number of refinement/unrefinement steps done so far.
void extendMarkedCells(bitSet &markedCell) const
Extend markedCell with cell-face-cell.
const bitSet & protectedCell() const
Cells which should not be refined/unrefined.
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))
void readDict()
Read the projection parameters from dictionary.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
scalarField error(const scalarField &fld, const scalar minLevel, const scalar maxLevel) const
Namespace for OpenFOAM.
scalarField maxCellField(const volScalarField &) const
Get point max of connected cell.