dynamicRefineFvMesh.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-2021 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::dynamicRefineFvMesh
26 
27 Description
28  Dynamic mesh refinement/unrefinement based on volScalarField values.
29 
30  Refinement can optionally be specified in a cellZone or in multiple
31  regions, each controlled by a different volScalarField.
32 
33 Usage
34  Example of single field based refinement in all cells:
35  \verbatim
36  dynamicFvMesh dynamicRefineFvMesh;
37 
38  // How often to refine
39  refineInterval 1;
40 
41  // Field to be refinement on
42  field alpha.water;
43 
44  // Refine field in between lower..upper
45  lowerRefineLevel 0.001;
46  upperRefineLevel 0.999;
47 
48  // If value < unrefineLevel unrefine
49  unrefineLevel 10;
50 
51  // If value < unrefineLevel (default=great) unrefine
52  // unrefineLevel 10;
53  // Have slower than 2:1 refinement
54  nBufferLayers 1;
55 
56  // Refine cells only up to maxRefinement levels
57  maxRefinement 2;
58 
59  // Stop refinement if maxCells reached
60  maxCells 200000;
61 
62  // Flux field and corresponding velocity field. Fluxes on changed
63  // faces get recalculated by interpolating the velocity. Use 'none'
64  // on surfaceScalarFields that do not need to be reinterpolated.
65  correctFluxes
66  (
67  (phi none)
68  (nHatf none)
69  (rhoPhi none)
70  (alphaPhi0.water none)
71  (ghf none)
72  );
73 
74  // Write the refinement level as a volScalarField
75  dumpLevel true;
76  \endverbatim
77 
78  Example of single field based refinement in two regions:
79  \verbatim
80  dynamicFvMesh dynamicRefineFvMesh;
81 
82  // How often to refine
83  refineInterval 1;
84 
85  refinementRegions
86  {
87  region1
88  {
89  cellZone refinementRegion1;
90 
91  // Field to be refinement on
92  field alpha.water;
93 
94  // Refine field in between lower..upper
95  lowerRefineLevel 0.001;
96  upperRefineLevel 0.999;
97 
98  // Refine cells only up to maxRefinement levels
99  maxRefinement 1;
100 
101  // If value < unrefineLevel unrefine
102  unrefineLevel 10;
103  }
104 
105  region2
106  {
107  cellZone refinementRegion2;
108 
109  // Field to be refinement on
110  field alpha.water;
111 
112  // Refine field in between lower..upper
113  lowerRefineLevel 0.001;
114  upperRefineLevel 0.999;
115 
116  // Refine cells only up to maxRefinement levels
117  maxRefinement 2;
118 
119  // If value < unrefineLevel unrefine
120  unrefineLevel 10;
121  }
122  }
123 
124  // If value < unrefineLevel (default=great) unrefine
125  // unrefineLevel 10;
126  // Have slower than 2:1 refinement
127  nBufferLayers 1;
128 
129  // Stop refinement if maxCells reached
130  maxCells 200000;
131 
132  // Flux field and corresponding velocity field. Fluxes on changed
133  // faces get recalculated by interpolating the velocity. Use 'none'
134  // on surfaceScalarFields that do not need to be reinterpolated.
135  correctFluxes
136  (
137  (phi none)
138  (nHatf none)
139  (rhoPhi none)
140  (alphaPhi0.water none)
141  (ghf none)
142  );
143 
144  // Write the refinement level as a volScalarField
145  dumpLevel true;
146  \endverbatim
147 
148 
149 SourceFiles
150  dynamicRefineFvMesh.C
151 
152 \*---------------------------------------------------------------------------*/
153 
154 #ifndef dynamicRefineFvMesh_H
155 #define dynamicRefineFvMesh_H
156 
157 #include "dynamicFvMesh.H"
158 #include "hexRef8.H"
159 #include "PackedBoolList.H"
160 #include "Switch.H"
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 namespace Foam
165 {
166 
167 /*---------------------------------------------------------------------------*\
168  Class dynamicRefineFvMesh Declaration
169 \*---------------------------------------------------------------------------*/
172 :
173  public dynamicFvMesh
174 {
175  // Private Member Data
176 
177  //- Mesh cutting engine
178  hexRef8 meshCutter_;
179 
180  //- Dump cellLevel for postprocessing
181  Switch dumpLevel_;
182 
183  //- Fluxes to map
184  HashTable<word> correctFluxes_;
185 
186  //- Number of refinement/unrefinement steps done so far.
187  label nRefinementIterations_;
188 
189  //- Protected cells (usually since not hexes)
190  PackedBoolList protectedCells_;
191 
192 
193  // Private Member Functions
194 
195  //- Count set/unset elements in packedlist.
196  static label count(const PackedBoolList&, const unsigned int);
197 
198  //- Calculate cells that cannot be refined since would trigger
199  // refinement of protectedCells_ (since 2:1 refinement cascade)
200  void calculateProtectedCells(PackedBoolList& unrefineableCells) const;
201 
202  //- Read the projection parameters from dictionary
203  void readDict();
204 
205 
206  //- Refine cells. Update mesh and fields.
207  autoPtr<mapPolyMesh> refine(const labelList&);
208 
209  //- Unrefine cells. Gets passed in centre points of cells to combine.
210  autoPtr<mapPolyMesh> unrefine(const labelList&);
211 
212 
213  // Selection of cells to un/refine
214 
215  const cellZone& findCellZone
216  (
217  const word& cellZoneName
218  ) const;
219 
220  scalarField cellToPoint(const scalarField& vFld) const;
221 
223  (
224  const scalarField& fld,
225  const scalar minLevel,
226  const scalar maxLevel
227  ) const;
228 
229  scalarField error
230  (
231  const scalarField& fld,
232  const labelList& cells,
233  const scalar minLevel,
234  const scalar maxLevel
235  ) const;
236 
237  //- Select candidate cells for refinement
238  virtual void selectRefineCandidates
239  (
240  PackedBoolList& candidateCell,
241  const scalar lowerRefineLevel,
242  const scalar upperRefineLevel,
243  const scalar maxRefinement,
244  const scalarField& vFld
245  ) const;
246 
247  //- Select candidate cells for refinement
248  virtual void selectRefineCandidates
249  (
250  PackedBoolList& candidateCell,
251  const scalar lowerRefineLevel,
252  const scalar upperRefineLevel,
253  const scalar maxRefinement,
254  const scalarField& vFld,
255  const labelList& cells
256  ) const;
257 
258  //- Select candidate cells for refinement
259  virtual scalar selectRefineCandidates
260  (
261  PackedBoolList& candidateCell,
262  const dictionary& refineDict
263  ) const;
264 
265  //- Subset candidate cells for refinement
266  virtual labelList selectRefineCells
267  (
268  const label maxCells,
269  const label maxRefinement,
270  const PackedBoolList& candidateCell
271  ) const;
272 
273  //- Get point max of connected cell
274  void selectUnrefineCandidates
275  (
276  boolList& unrefineCandidates,
277  const volScalarField&,
278  const scalar unrefineLevel
279  ) const;
280 
281  //- Get point max of connected cell
282  void selectUnrefineCandidates
283  (
284  boolList& unrefineCandidates,
285  const volScalarField&,
286  const cellZone&,
287  const scalar unrefineLevel
288  ) const;
289 
290  //- Get point max of connected cell
291  void selectUnrefineCandidates
292  (
293  boolList& unrefineCandidates,
294  const dictionary& refineDict
295  ) const;
296 
297  //- Select points that can be unrefined.
298  virtual labelList selectUnrefinePoints
299  (
300  const PackedBoolList& markedCell,
301  const boolList& unrefineCandidates
302  ) const;
303 
304  //- Extend markedCell with cell-face-cell.
305  void extendMarkedCells(PackedBoolList& markedCell) const;
306 
307  //- Check all cells have 8 anchor points
308  void checkEightAnchorPoints
309  (
311  label& nProtected
312  ) const;
313 
314 
315 public:
316 
317  //- Runtime type information
318  TypeName("dynamicRefineFvMesh");
319 
320 
321  // Constructors
322 
323  //- Construct from IOobject
324  explicit dynamicRefineFvMesh(const IOobject& io);
325 
326  //- Disallow default bitwise copy construction
327  dynamicRefineFvMesh(const dynamicRefineFvMesh&) = delete;
328 
329 
330  //- Destructor
331  virtual ~dynamicRefineFvMesh();
332 
333 
334  // Member Functions
335 
336  //- Direct access to the refinement engine
337  const hexRef8& meshCutter() const
338  {
339  return meshCutter_;
340  }
341 
342  //- Cells which should not be refined/unrefined
343  const PackedBoolList& protectedCell() const
344  {
345  return protectedCells_;
346  }
347 
348  //- Cells which should not be refined/unrefined
350  {
351  return protectedCells_;
352  }
353 
354  //- Update the mesh for both mesh motion and topology change
355  virtual bool update();
356 
357 
358  // Writing
359 
360  //- Write using given format, version and compression
361  virtual bool writeObject
362  (
366  const bool write = true
367  ) const;
368 
369 
370  // Member Operators
371 
372  //- Disallow default bitwise assignment
373  void operator=(const dynamicRefineFvMesh&) = delete;
374 };
375 
376 
377 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
378 
379 } // End namespace Foam
380 
381 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
382 
383 #endif
384 
385 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const hexRef8 & meshCutter() const
Direct access to the refinement engine.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write using given format, version and compression.
const cellList & cells() const
TypeName("dynamicRefineFvMesh")
Runtime type information.
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1093
Dynamic mesh refinement/unrefinement based on volScalarField values.
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:66
virtual bool update()
Update the mesh for both mesh motion and topology change.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void operator=(const dynamicRefineFvMesh &)=delete
Disallow default bitwise assignment.
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:64
A class for handling words, derived from string.
Definition: word.H:59
const PackedBoolList & protectedCell() const
Cells which should not be refined/unrefined.
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
An STL-conforming hash table.
Definition: HashTable.H:61
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
A topoSetSource to select points based on usage in cells.
Definition: cellToPoint.H:49
A subset of mesh cells.
Definition: cellZone.H:61
A bit-packed bool list.
dynamicRefineFvMesh(const IOobject &io)
Construct from IOobject.
Abstract base class for geometry and/or topology changing fvMesh.
Definition: dynamicFvMesh.H:49
Version number type.
Definition: IOstream.H:96
virtual ~dynamicRefineFvMesh()
Destructor.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.