NEMoSys  0.63.0
A modular, extensible resource with robust automated mesh generation, mesh quality analysis, adaptive mesh refinement, and data transfer between arbitrary meshes.
AMRFoam.H
Go to the documentation of this file.
1 /*******************************************************************************
2 * Promesh *
3 * Copyright (C) 2022, IllinoisRocstar LLC. All rights reserved. *
4 * *
5 * Promesh is the property of IllinoisRocstar LLC. *
6 * *
7 * IllinoisRocstar LLC *
8 * Champaign, IL *
9 * www.illinoisrocstar.com *
10 * promesh@illinoisrocstar.com *
11 *******************************************************************************/
12 /*******************************************************************************
13 * This file is part of Promesh *
14 * *
15 * This version of Promesh is free software: you can redistribute it and/or *
16 * modify it under the terms of the GNU Lesser General Public License as *
17 * published by the Free Software Foundation, either version 3 of the License, *
18 * or (at your option) any later version. *
19 * *
20 * Promesh is distributed in the hope that it will be useful, but WITHOUT ANY *
21 * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS *
22 * FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more *
23 * details. *
24 * *
25 * You should have received a copy of the GNU Lesser General Public License *
26 * along with this program. If not, see <https://www.gnu.org/licenses/>. *
27 * *
28 *******************************************************************************/
29 #ifndef NEMOSYS_AMRFOAM_H_
30 #define NEMOSYS_AMRFOAM_H_
31 
32 #include <vector>
33 
34 #include <bitSet.H>
35 #include <Switch.H>
36 #include <dynamicFvMesh.H>
37 #include <hexRef8.H>
38 #include <motionSolver.H>
39 
40 namespace Foam {
41 
42 /** @brief AMRFoam class is a standalone adaptive mesh refinement application
43  that uses polyMesh as its base and provides useful interface for
44  use of different solvers and mesh formats. Initially the supported
45  refinement criterias would be size field defined over mesh as well
46  as a gradient field of any simulation quantity.
47 
48  Any object/pointer created for this class in C++ runtime will act
49  as finiteVolume mesh of OpenFOAM.
50 **/
51 // class motionSolver;
52 class AMRFoam : public dynamicFvMesh {
53  // Constructor and destructor
54  public:
55  /** @brief Standard constructor. Accepts polyMesh object for initialization.
56  It needs dynamicMeshDict present in constant folder. Method for
57  same is available in MeshManipulationFoam class as
58  "createDynamicMeshDict"
59  @param iomesh A polyMesh input object
60  **/
61  AMRFoam(const Foam::IOobject &iomesh);
62 
63  /** @brief Class destructor
64  **/
66  // Nothing
67  }
68 
69  public:
70  /** @brief updateAMR is a master method for this process. It drives the
71  refinement and unrefinement procedure between any two given
72  timesteps.
73  @param refineInterval Refinement interval/frequency
74  @param maxRefinement Maximum refinement allowed by user
75  @param vFld A scalar field of refinement criteria quantity
76  @param lowerRefineLevel Lower refinement level (Based on vFld)
77  @param upperRefineLevel Upper refinement level (Based on vFld)
78  @param unrefineAbove Criteria for unrefinement above this value
79  @param unrefineBelow Criteria for unrefinement below this value
80  @param nBufferLayers Buffer layers for slower than 2:1 refinement
81  @param maxCells Maximum number of cells allowed in mesh. Once this level
82  is reached, refinement will be stopped.
83  @return A boolean indicating mesh change
84  **/
85  bool updateAMR(const int &refineInterval, const int &maxRefinement,
86  volScalarField &vFld, const double &lowerRefineLevel,
87  const double &upperRefineLevel, const double &unrefineAbove,
88  const double &unrefineBelow, const int &nBufferLayers,
89  const int &maxCells);
90 
91  /** @brief New method for machine learning ML. Takes output from ML model
92  and performs refinement/unrefinement.
93  @param refineInterval Refinement interval/frequency
94  @param maxRefinement Maximum refinement allowed by user
95  @param nBufferLayers Buffer layers for slower than 2:1 refinement
96  @param maxCells Maximum number of cells allowed in mesh. Once this level
97  is reached, refinement will be stopped.
98  @param vFld Machine learning model output vector
99  @return A boolean indicating mesh change
100  **/
101  bool updateAMRML(const int &refineInterval, const int &maxRefinement,
102  const int &nBufferLayers, const int &maxCells,
103  volScalarField &vFld);
104 
105  /** @brief A mesh motion solver. This method allows solving of mesh motion
106  along with mesh refinment. Flux correction takes place after motion
107  and refinement operations are performed. Experimental method.
108  **/
109  void solveMotion();
110 
111  /** @brief Converts volume field to point field
112  @param vFld Volume Field.
113  @return Point Field.
114  **/
115  Foam::scalarField cell2Pt(const scalarField &vFld);
116 
117  /** @brief Writes refinement history in current time directory.
118  **/
119  void writeHistory();
120 
121  /** @brief Access to current mesh at every time step
122  @return A polyMesh object
123  **/
124  // const polyMesh& getMesh();
125  const polyMesh &getMesh();
126 
127  /** @brief Ability to write current mesh into current time directory
128  **/
129  void writeMesh();
130 
131  /** @brief Ability to explicitely disable motion
132  **/
133  void disableMotion();
134 
135  /** @brief Converts incoming size-field/solution-field into pointScalarField
136  defined over all points of mesh
137  @param inName
138  @return A pointScalarField with appended incoming values
139  **/
140  pointScalarField readIncomingPtField(const std::string &inName,
141  const std::string &outName);
142 
143  /** @brief Converts pointScalarField to volScalarField
144  @param pntF A point field defined on mesh points
145  @param outName Object name for output scalar
146  @return A Scalar Field defined on cells
147  **/
148  volScalarField pF2vF(const pointScalarField &pntF,
149  const std::string &outName);
150 
151  /** @brief Reads the refinement field and registers it.
152  @param fldName Name of field present in current runTime directory
153  **/
154  volScalarField readInitialField(const std::string &fldName);
155 
156  /** @brief Creates an empty pointScalarField based on current mesh
157  @param fldName Name of field to create in current runTime directory
158  **/
159  pointScalarField initPntSField(const std::string &fldName);
160 
161  /** @brief Creates an empty volScalarField based on current mesh
162  @param fldName Name of field to create in current runTime directory
163  **/
164  volScalarField initScalarField(const std::string &fldName);
165 
166  /** @brief Creates an empty volScalarField for gradient magnitude
167  @param fldName Name of field to create in current runTime directory
168  **/
169  volScalarField initGradientField(const std::string &fldName);
170 
171  /** @brief Converts vector to volScalarField for refinement/unrefinement
172  @param vec Incoming vector
173  @return volScalarField
174  **/
175  volScalarField assignToVolScalarField(const std::vector<int> &vec);
176 
177  /** @brief Access to updated field at every time step
178  **/
179  void enableUpdatedField();
180 
181  /** @brief Access to refinement history data
182  **/
183  void enableRefHistoryData();
184 
185  /** @brief Access to mesh data
186  **/
187  void enableMeshWriting();
188 
189  /** @brief reads incoming field and converts to volScalarField
190  @param inName Input file name
191  @param outName Object name for volScalarField
192  @return Scalar Field
193  **/
194  volScalarField readIncomingCellField(const std::string &inName,
195  const std::string &outName);
196 
197  /** @brief Transforms incoming scalar field into refinement levels understood
198  by OpenFOAM environement. Appends this information to lower, upper
199  and unrefine level vectors for AMR
200  **/
201  void transformField(volScalarField &inField);
202 
203  /** @brief Returns a gradient of scalar field over mesh
204  param fldGrd Scalar Field for Gradient Calculation
205  @return scalar field of gradient over mesh
206  **/
207  volScalarField getGradient(const volScalarField &fldGrd);
208 
209  /** @brief Vector to store lower refinement levels
210  **/
211  std::vector<double> lowerRefLvl;
212 
213  /** @brief Vector to store upper refinement levels
214  **/
215  std::vector<double> upperRefLvl;
216 
217  /** @brief Vector to store unrefine levels
218  **/
219  std::vector<double> unrefAboveLvl;
220 
221  /** @brief Vector to store unrefine levels
222  **/
223  std::vector<double> unrefBelowLvl;
224 
225  /** @brief Step size for one time refinement
226  **/
227  int stepSz = 1;
228 
229  /**@brief Boolean to check if mesh motion is needed
230  **/
231  bool enableMotion{false};
232 
233  // Private Attributes
234  private:
235  /** @brief Boolean to enable writing of updated solution field every timestep
236  **/
237  bool writeField{false};
238 
239  /** @brief Boolean to enable writing of refinement history
240  **/
241  bool writeRefHistory{false};
242 
243  /** @brief Boolean to enable writing of mesh
244  **/
245  bool writeMeshData{false};
246 
247  /** @brief Vector to store user defined vector
248  **/
249  std::vector<double> incomingField;
250 
251  protected:
252  //- Mesh cutting engine
253  hexRef8 meshCutter_;
254 
255  //- Dump cellLevel for postprocessing
257 
258  //- Fluxes to map
259  HashTable<word> correctFluxes_;
260 
261  //- Number of refinement/unrefinement steps done so far.
263 
264  //- Protected cells (usually since not hexes)
266 
267  // Protected Member Functions
268  /** @brief Checks for availability of mesh motion keywords in dictionary
269  **/
270  void checkForMotion();
271 
272  /** @brief Calculates and stores non-hex cells in protected cells catagory to
273  avoid their refinement
274  **/
275  void calcProtectedCells();
276  /** @brief Unrefines selected cells in a mesh.
277  @param splitPoints Selected hanging nodes that can be erased to unrefine
278  previously refined cells. Selection is performed using criteria and
279  protected cells are excluded.
280  @return Mapped Polymesh.
281  **/
282  virtual Foam::autoPtr<Foam::mapPolyMesh> refine(const labelList &);
283 
284  /** @brief Refines selected cells in a mesh.
285  @param cellsToRefine Cell list selected for refinement.
286  @return Mapped Polymesh.
287  **/
288  virtual Foam::autoPtr<Foam::mapPolyMesh> unrefine(const labelList &);
289 
290  //- Calculate cells that cannot be refined since would trigger
291  // refinement of protectedCell_ (since 2:1 refinement cascade)
292  void calculateProtectedCells(bitSet &unrefineableCell) const;
293 
294  //- Read the projection parameters from dictionary
295  void readDict();
296 
297  //- Calculates approximate value for refinement level so
298  // we don't go above maxCell
299  scalar getRefineLevel(const label maxCells, const label maxRefinement,
300  const scalar refineLevel, const scalarField &) const;
301 
302  //- Get per cell max of connected point
303  scalarField maxPointField(const scalarField &) const;
304 
305  //- Get point max of connected cell
306  scalarField maxCellField(const volScalarField &) const;
307 
308  scalarField cellToPoint(const scalarField &vFld) const;
309 
310  scalarField error(const scalarField &fld, const scalar minLevel,
311  const scalar maxLevel) const;
312 
313  //- Select candidate cells for refinement
314  virtual void selectRefineCandidates(const scalar lowerRefineLevel,
315  const scalar upperRefineLevel,
316  const scalarField &vFld,
317  bitSet &candidateCell) const;
318 
319  //- Subset candidate cells for refinement
320  virtual labelList selectRefineCells(const label maxCells,
321  const label maxRefinement,
322  const bitSet &candidateCell) const;
323 
324  /** @brief Marked the cells to unrefine based on user defined levels
325  @param unrefineAbove Unrefine cells above this thresold.
326  @param unrefineBelow Unrefine cells below this thresold.
327  @param markedCell Unrefinable cells list
328  @param pFld PointField at current runTime
329  @return List of cells/points to unrefine.
330  **/
331  virtual labelList selectUnrefinePoints(const scalar unrefineAbove,
332  const scalar unrefineBelow,
333  const bitSet &markedCell,
334  const scalarField &pFld) const;
335 
336  //- Extend markedCell with cell-face-cell.
337  void extendMarkedCells(bitSet &markedCell) const;
338 
339  //- Check all cells have 8 anchor points
340  void checkEightAnchorPoints(bitSet &protectedCell) const;
341 
342  //- Map single non-flux surface<Type>Field
343  // for new internal faces (e.g. AMR refine). This currently
344  // interpolates values from surrounding faces (faces on
345  // neighbouring cells) that do have a value.
346  template <class T>
347  void mapNewInternalFaces(const labelList &faceMap,
348  GeometricField<T, fvsPatchField, surfaceMesh> &);
349 
350  //- Correct surface fields for new faces
351  template <class T>
352  void mapNewInternalFaces(const labelList &faceMap);
353 
354  //- Correct surface fields for new faces. Converts any oriented
355  // fields into non-oriented (i.e. phi -> Uf) before mapping
356  template <class T>
357  void mapNewInternalFaces(const surfaceVectorField &Sf,
358  const surfaceScalarField &magSf,
359  const labelList &faceMap);
360 
361  public:
362  virtual bool update();
363 
364  // Pointer for motion solver
365  autoPtr<motionSolver> motionPtr_;
366 
367  // Member Functions
368 
369  //- Direct access to the refinement engine
370  const hexRef8 &meshCutter() const { return meshCutter_; }
371 
372  //- Cells which should not be refined/unrefined
373  const bitSet &protectedCell() const { return protectedCell_; }
374 
375  //- Cells which should not be refined/unrefined
376  bitSet &protectedCell() { return protectedCell_; }
377 
378  //- Write using given format, version and compression
379  bool writeObjectAMR(const bool valid) const;
380 
381  //- Map all fields in time using given map.
382  virtual void mapFields(const mapPolyMesh &mpm);
383 };
384 
385 } // namespace Foam
386 
387 #endif // NEMOSYS_AMRFOAM_H_
virtual Foam::autoPtr< Foam::mapPolyMesh > unrefine(const labelList &)
Refines selected cells in a mesh.
Definition: AMRFoam.C:1138
pointScalarField initPntSField(const std::string &fldName)
Creates an empty pointScalarField based on current mesh.
Definition: AMRFoam.C:814
std::vector< double > incomingField
Vector to store user defined vector.
Definition: AMRFoam.H:249
scalarField error(const scalarField &fld, const scalar minLevel, const scalar maxLevel) const
Definition: AMRFoam.C:1299
scalar getRefineLevel(const label maxCells, const label maxRefinement, const scalar refineLevel, const scalarField &) const
bitSet protectedCell_
Definition: AMRFoam.H:265
void enableMeshWriting()
Access to mesh data.
Definition: AMRFoam.C:832
virtual labelList selectUnrefinePoints(const scalar unrefineAbove, const scalar unrefineBelow, const bitSet &markedCell, const scalarField &pFld) const
Marked the cells to unrefine based on user defined levels.
Definition: AMRFoam.C:616
void enableUpdatedField()
Access to updated field at every time step.
Definition: AMRFoam.C:830
bool enableMotion
Boolean to check if mesh motion is needed.
Definition: AMRFoam.H:231
bool writeObjectAMR(const bool valid) const
Definition: AMRFoam.C:1425
void calculateProtectedCells(bitSet &unrefineableCell) const
Definition: AMRFoam.C:847
~AMRFoam()
Class destructor.
Definition: AMRFoam.H:65
virtual bool update()
Definition: AMRFoam.C:233
void checkForMotion()
Checks for availability of mesh motion keywords in dictionary.
Definition: AMRFoam.C:72
void readDict()
Definition: AMRFoam.C:920
Foam::scalarField cell2Pt(const scalarField &vFld)
Converts volume field to point field.
Definition: AMRFoam.C:680
void solveMotion()
A mesh motion solver.
Definition: AMRFoam.C:1449
volScalarField getGradient(const volScalarField &fldGrd)
Returns a gradient of scalar field over mesh param fldGrd Scalar Field for Gradient Calculation...
Definition: AMRFoam.C:834
AMRFoam class is a standalone adaptive mesh refinement application that uses polyMesh as its base and...
Definition: AMRFoam.H:52
void mapNewInternalFaces(const labelList &faceMap, GeometricField< T, fvsPatchField, surfaceMesh > &)
Definition: AMRFoam.C:1458
label nRefinementIterations_
Definition: AMRFoam.H:262
volScalarField assignToVolScalarField(const std::vector< int > &vec)
Converts vector to volScalarField for refinement/unrefinement.
Definition: AMRFoam.C:750
void calcProtectedCells()
Calculates and stores non-hex cells in protected cells catagory to avoid their refinement.
Definition: AMRFoam.C:86
const hexRef8 & meshCutter() const
Definition: AMRFoam.H:370
virtual void mapFields(const mapPolyMesh &mpm)
Definition: AMRFoam.C:937
void writeHistory()
Writes refinement history in current time directory.
Definition: AMRFoam.C:696
AMRFoam(const Foam::IOobject &iomesh)
Standard constructor.
Definition: AMRFoam.C:62
HashTable< word > correctFluxes_
Definition: AMRFoam.H:259
bool updateAMR(const int &refineInterval, const int &maxRefinement, volScalarField &vFld, const double &lowerRefineLevel, const double &upperRefineLevel, const double &unrefineAbove, const double &unrefineBelow, const int &nBufferLayers, const int &maxCells)
updateAMR is a master method for this process.
Definition: AMRFoam.C:386
bool writeRefHistory
Boolean to enable writing of refinement history.
Definition: AMRFoam.H:241
hexRef8 meshCutter_
Definition: AMRFoam.H:253
volScalarField initScalarField(const std::string &fldName)
Creates an empty volScalarField based on current mesh.
Definition: AMRFoam.C:798
volScalarField initGradientField(const std::string &fldName)
Creates an empty volScalarField for gradient magnitude.
Definition: AMRFoam.C:806
const polyMesh & getMesh()
Access to current mesh at every time step.
Definition: AMRFoam.C:702
void enableRefHistoryData()
Access to refinement history data.
Definition: AMRFoam.C:828
scalarField maxCellField(const volScalarField &) const
Definition: AMRFoam.C:1273
const bitSet & protectedCell() const
Definition: AMRFoam.H:373
bool updateAMRML(const int &refineInterval, const int &maxRefinement, const int &nBufferLayers, const int &maxCells, volScalarField &vFld)
New method for machine learning ML.
Definition: AMRFoam.C:510
void transformField(volScalarField &inField)
Transforms incoming scalar field into refinement levels understood by OpenFOAM environement.
Definition: AMRFoam.C:708
virtual void selectRefineCandidates(const scalar lowerRefineLevel, const scalar upperRefineLevel, const scalarField &vFld, bitSet &candidateCell) const
Definition: AMRFoam.C:1310
std::vector< double > unrefAboveLvl
Vector to store unrefine levels.
Definition: AMRFoam.H:219
pointScalarField readIncomingPtField(const std::string &inName, const std::string &outName)
Converts incoming size-field/solution-field into pointScalarField defined over all points of mesh...
Definition: AMRFoam.C:756
void checkEightAnchorPoints(bitSet &protectedCell) const
Definition: AMRFoam.C:1401
void disableMotion()
Ability to explicitely disable motion.
Definition: AMRFoam.C:84
void extendMarkedCells(bitSet &markedCell) const
Definition: AMRFoam.C:1379
void writeMesh()
Ability to write current mesh into current time directory.
Definition: AMRFoam.C:704
scalarField maxPointField(const scalarField &) const
Definition: AMRFoam.C:1260
int stepSz
Step size for one time refinement.
Definition: AMRFoam.H:227
std::vector< double > unrefBelowLvl
Vector to store unrefine levels.
Definition: AMRFoam.H:223
volScalarField readIncomingCellField(const std::string &inName, const std::string &outName)
reads incoming field and converts to volScalarField
Definition: AMRFoam.C:723
autoPtr< motionSolver > motionPtr_
Definition: AMRFoam.H:365
bool dumpLevel_
Definition: AMRFoam.H:256
scalarField cellToPoint(const scalarField &vFld) const
Definition: AMRFoam.C:1286
volScalarField pF2vF(const pointScalarField &pntF, const std::string &outName)
Converts pointScalarField to volScalarField.
Definition: AMRFoam.C:792
std::vector< double > upperRefLvl
Vector to store upper refinement levels.
Definition: AMRFoam.H:215
volScalarField readInitialField(const std::string &fldName)
Reads the refinement field and registers it.
Definition: AMRFoam.C:784
bitSet & protectedCell()
Definition: AMRFoam.H:376
bool writeMeshData
Boolean to enable writing of mesh.
Definition: AMRFoam.H:245
std::vector< double > lowerRefLvl
Vector to store lower refinement levels.
Definition: AMRFoam.H:211
bool writeField
Boolean to enable writing of updated solution field every timestep.
Definition: AMRFoam.H:237
virtual Foam::autoPtr< Foam::mapPolyMesh > refine(const labelList &)
Unrefines selected cells in a mesh.
Definition: AMRFoam.C:1082
virtual labelList selectRefineCells(const label maxCells, const label maxRefinement, const bitSet &candidateCell) const
Definition: AMRFoam.C:1325
Definition: AMRFoam.C:60