29 #ifndef NEMOSYS_AMRFOAM_H_ 30 #define NEMOSYS_AMRFOAM_H_ 36 #include <dynamicFvMesh.H> 38 #include <motionSolver.H> 61 AMRFoam(
const Foam::IOobject &iomesh);
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,
101 bool updateAMRML(
const int &refineInterval,
const int &maxRefinement,
102 const int &nBufferLayers,
const int &maxCells,
103 volScalarField &vFld);
115 Foam::scalarField
cell2Pt(
const scalarField &vFld);
141 const std::string &outName);
148 volScalarField
pF2vF(
const pointScalarField &pntF,
149 const std::string &outName);
195 const std::string &outName);
207 volScalarField
getGradient(
const volScalarField &fldGrd);
282 virtual Foam::autoPtr<Foam::mapPolyMesh>
refine(
const labelList &);
288 virtual Foam::autoPtr<Foam::mapPolyMesh>
unrefine(
const labelList &);
299 scalar
getRefineLevel(
const label maxCells,
const label maxRefinement,
300 const scalar refineLevel,
const scalarField &)
const;
308 scalarField
cellToPoint(
const scalarField &vFld)
const;
310 scalarField
error(
const scalarField &fld,
const scalar minLevel,
311 const scalar maxLevel)
const;
315 const scalar upperRefineLevel,
316 const scalarField &vFld,
317 bitSet &candidateCell)
const;
321 const label maxRefinement,
322 const bitSet &candidateCell)
const;
332 const scalar unrefineBelow,
333 const bitSet &markedCell,
334 const scalarField &pFld)
const;
348 GeometricField<T, fvsPatchField, surfaceMesh> &);
358 const surfaceScalarField &magSf,
359 const labelList &faceMap);
382 virtual void mapFields(
const mapPolyMesh &mpm);
387 #endif // NEMOSYS_AMRFOAM_H_ virtual Foam::autoPtr< Foam::mapPolyMesh > unrefine(const labelList &)
Refines selected cells in a mesh.
pointScalarField initPntSField(const std::string &fldName)
Creates an empty pointScalarField based on current mesh.
std::vector< double > incomingField
Vector to store user defined vector.
scalarField error(const scalarField &fld, const scalar minLevel, const scalar maxLevel) const
scalar getRefineLevel(const label maxCells, const label maxRefinement, const scalar refineLevel, const scalarField &) const
void enableMeshWriting()
Access to mesh data.
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.
void enableUpdatedField()
Access to updated field at every time step.
bool enableMotion
Boolean to check if mesh motion is needed.
bool writeObjectAMR(const bool valid) const
void calculateProtectedCells(bitSet &unrefineableCell) const
~AMRFoam()
Class destructor.
void checkForMotion()
Checks for availability of mesh motion keywords in dictionary.
Foam::scalarField cell2Pt(const scalarField &vFld)
Converts volume field to point field.
void solveMotion()
A mesh motion solver.
volScalarField getGradient(const volScalarField &fldGrd)
Returns a gradient of scalar field over mesh param fldGrd Scalar Field for Gradient Calculation...
AMRFoam class is a standalone adaptive mesh refinement application that uses polyMesh as its base and...
void mapNewInternalFaces(const labelList &faceMap, GeometricField< T, fvsPatchField, surfaceMesh > &)
label nRefinementIterations_
volScalarField assignToVolScalarField(const std::vector< int > &vec)
Converts vector to volScalarField for refinement/unrefinement.
void calcProtectedCells()
Calculates and stores non-hex cells in protected cells catagory to avoid their refinement.
const hexRef8 & meshCutter() const
virtual void mapFields(const mapPolyMesh &mpm)
void writeHistory()
Writes refinement history in current time directory.
AMRFoam(const Foam::IOobject &iomesh)
Standard constructor.
HashTable< word > correctFluxes_
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.
bool writeRefHistory
Boolean to enable writing of refinement history.
volScalarField initScalarField(const std::string &fldName)
Creates an empty volScalarField based on current mesh.
volScalarField initGradientField(const std::string &fldName)
Creates an empty volScalarField for gradient magnitude.
const polyMesh & getMesh()
Access to current mesh at every time step.
void enableRefHistoryData()
Access to refinement history data.
scalarField maxCellField(const volScalarField &) const
const bitSet & protectedCell() const
bool updateAMRML(const int &refineInterval, const int &maxRefinement, const int &nBufferLayers, const int &maxCells, volScalarField &vFld)
New method for machine learning ML.
void transformField(volScalarField &inField)
Transforms incoming scalar field into refinement levels understood by OpenFOAM environement.
virtual void selectRefineCandidates(const scalar lowerRefineLevel, const scalar upperRefineLevel, const scalarField &vFld, bitSet &candidateCell) const
std::vector< double > unrefAboveLvl
Vector to store unrefine levels.
pointScalarField readIncomingPtField(const std::string &inName, const std::string &outName)
Converts incoming size-field/solution-field into pointScalarField defined over all points of mesh...
void checkEightAnchorPoints(bitSet &protectedCell) const
void disableMotion()
Ability to explicitely disable motion.
void extendMarkedCells(bitSet &markedCell) const
void writeMesh()
Ability to write current mesh into current time directory.
scalarField maxPointField(const scalarField &) const
int stepSz
Step size for one time refinement.
std::vector< double > unrefBelowLvl
Vector to store unrefine levels.
volScalarField readIncomingCellField(const std::string &inName, const std::string &outName)
reads incoming field and converts to volScalarField
autoPtr< motionSolver > motionPtr_
scalarField cellToPoint(const scalarField &vFld) const
volScalarField pF2vF(const pointScalarField &pntF, const std::string &outName)
Converts pointScalarField to volScalarField.
std::vector< double > upperRefLvl
Vector to store upper refinement levels.
volScalarField readInitialField(const std::string &fldName)
Reads the refinement field and registers it.
bool writeMeshData
Boolean to enable writing of mesh.
std::vector< double > lowerRefLvl
Vector to store lower refinement levels.
bool writeField
Boolean to enable writing of updated solution field every timestep.
virtual Foam::autoPtr< Foam::mapPolyMesh > refine(const labelList &)
Unrefines selected cells in a mesh.
virtual labelList selectRefineCells(const label maxCells, const label maxRefinement, const bitSet &candidateCell) const