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.
FoamRefineDriver.C
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 *******************************************************************************/
30 
31 #include <pointFieldsFwd.H>
32 #include <interpolatePointToCell.H>
33 #include "getDicts.H"
34 #include "Refinement/AMRFoam.H"
35 #include "AuxiliaryFunctions.H"
36 #include "Mesh/foamGeoMesh.H"
37 #include "Mesh/geoMeshFactory.H"
38 
39 #ifdef MLAMR
40 # include <Mesh/meshBase.H>
41 # include <fdeep/fdeep.hpp>
42 #endif
43 
44 // helpers for execute
45 namespace {
46 std::unique_ptr<Foam::AMRFoam> initializeRefine(
47  const NEM::DRV::RefineDriver::Files &files,
48  const NEM::DRV::FoamRefineOptsBase &opts, const Foam::Time &runTime) {
49  // Reading mesh from vtk file.
50  auto vgm = NEM::MSH::Read(files.inputMeshFile);
51  auto mesh = new NEM::MSH::foamGeoMesh();
52  mesh->takeGeoMesh(vgm);
53  mesh->write("");
54  mesh->Delete();
55  vgm->Delete();
56 
57  Foam::polyMesh mesh1(Foam::IOobject(Foam::polyMesh::defaultRegion,
58  runTime.timeName(), runTime,
59  Foam::IOobject::MUST_READ));
60 
61  // Initialized AMR class
62  auto amr = std::unique_ptr<Foam::AMRFoam>{new Foam::AMRFoam(mesh1)};
63 
64  if (opts.writeFieldData)
65  if (!opts.writeMesh) amr->enableUpdatedField();
66 
67  if (opts.writeRefHistory)
68  if (!opts.writeMesh) amr->enableRefHistoryData();
69 
70  if (opts.writeMesh) amr->enableMeshWriting();
71 
72  return amr;
73 }
74 
75 void finishRefine(Foam::AMRFoam *amr, const std::string &outputMeshFile,
76  const Foam::Time &runTime) {
77  amr->writeMesh();
78 
79  // Export final mesh to vtu format
80  auto fm = new Foam::fvMesh(Foam::IOobject("", runTime.timeName(), runTime,
81  Foam::IOobject::MUST_READ));
82  auto fgm_ = new NEM::MSH::foamGeoMesh(fm);
83  auto mesh = NEM::MSH::New(outputMeshFile);
84  mesh->takeGeoMesh(fgm_);
85  mesh->write(outputMeshFile);
86  mesh->Delete();
87  fgm_->Delete();
88 }
89 } // namespace
90 
91 namespace NEM {
92 namespace DRV {
93 
94 FoamRefineDriver::Opts::Opts(Criteria refCriteria, std::string inputFieldFile)
96  inputFieldFile(std::move(inputFieldFile)),
97  refCriteria(refCriteria) {}
98 
100  : RefineDriver(std::move(files)), opts_(std::move(opts)) {}
101 
103 
105  return opts_;
106 }
107 
108 void FoamRefineDriver::setOpts(Opts opts) { this->opts_ = std::move(opts); }
109 
111  // Initializes AMR workflow
112  std::unique_ptr<getDicts> initFoam;
113  initFoam = std::unique_ptr<getDicts>(new getDicts());
114  initFoam->writeBasicDicts("system/", this->opts_.startTime,
115  this->opts_.timeStep, this->opts_.endTime);
116  initFoam->createDynamicMeshDict(true);
117  auto controlDict_ = initFoam->createControlDict(true);
118  Foam::Time runTime(controlDict_.get(), ".", ".");
119 
120  auto amr = initializeRefine(this->files_, this->opts_, runTime);
121 
122  int nSteps = static_cast<int>(std::round(opts_.endTime / opts_.timeStep));
123 
124  // Creates scalar field from incoming text/CSV files
125  Foam::volScalarField meshFieldXY =
126  amr->readIncomingCellField(this->opts_.inputFieldFile, "meshFieldXY");
127 
129  for (int i = 0; i < nSteps; i++) {
130  amr->updateAMR(this->opts_.refineInterval, this->opts_.maxRefinement,
131  meshFieldXY, this->opts_.lowerRefineLevel,
132  this->opts_.upperRefineLevel, this->opts_.unrefineAbove,
133  this->opts_.unrefineBelow, this->opts_.nBufferLayers,
134  this->opts_.maxCells);
135 
136  runTime++;
137  }
138  meshFieldXY.write();
139  } else { // opts_.refCriteria == Opts::Criteria::GRADIENT
140  Foam::volScalarField magGrad = amr->getGradient(meshFieldXY);
141  for (int i = 0; i < nSteps; i++) {
142  amr->updateAMR(this->opts_.refineInterval, this->opts_.maxRefinement,
143  magGrad, this->opts_.lowerRefineLevel,
144  this->opts_.upperRefineLevel, this->opts_.unrefineAbove,
145  this->opts_.unrefineBelow, this->opts_.nBufferLayers,
146  this->opts_.maxCells);
147 
148  runTime++;
149  }
150  magGrad.write();
151  }
152 
153  finishRefine(amr.get(), this->files_.outputMeshFile, runTime);
154 
155  std::cout << "End!" << std::endl;
156 }
157 
158 #ifdef MLAMR
159 FoamMLRefineDriver::FoamMLRefineDriver(Files files, Opts opts)
160  : RefineDriver(std::move(files)), opts_(std::move(opts)) {}
161 
162 FoamMLRefineDriver::FoamMLRefineDriver() : FoamMLRefineDriver({{}, {}}, {}) {}
163 
164 const FoamMLRefineDriver::Opts &FoamMLRefineDriver::getOpts() const {
165  return opts_;
166 }
167 
168 void FoamMLRefineDriver::setOpts(Opts opts) { this->opts_ = std::move(opts); }
169 
170 void FoamMLRefineDriver::execute() const {
171  // Initializes AMR workflow
172  std::unique_ptr<getDicts> initFoam;
173  initFoam = std::unique_ptr<getDicts>(new getDicts());
174  initFoam->writeBasicDicts("system/", this->opts_.startTime,
175  this->opts_.timeStep, this->opts_.endTime);
176  initFoam->createDynamicMeshDict(true);
177  auto controlDict_ = initFoam->createControlDict(true);
178  Foam::Time runTime(controlDict_.get(), ".", ".");
179 
180  std::shared_ptr<meshBase> mesh =
182  auto amr = initializeRefine(this->files_, this->opts_, runTime);
183 
184  int nSteps = opts_.endTime / opts_.timeStep;
185 
186  std::vector<double> nonDimUGrad;
187  std::vector<double> X;
188  std::vector<double> Y;
189  std::vector<double> Z;
190  mesh->getCellDataArray("nonDimUGrad", nonDimUGrad);
191  mesh->getCellDataArray("X", X);
192  mesh->getCellDataArray("Y", Y);
193  mesh->getCellDataArray("Z", Z);
194 
195  // Loading ML Model
196  const auto model = fdeep::load_model(this->opts_.mlModel);
197 
198  // Making predictions
199  std::vector<int> refinementVec;
200  for (int i = 0; i < nonDimUGrad.size(); i++) { refinementVec.push_back(0); }
201  for (int i = 0; i < nonDimUGrad.size(); i++) {
202  const auto result = model.predict(
203  {fdeep::tensor(fdeep::tensor_shape(static_cast<double>(4)),
204  {X[i], Y[i], Z[i], nonDimUGrad[i]})});
205 
206  if (result[0].get(0, 0, 0, 0, 0) >= 0.5) refinementVec[i] = 1;
207  }
208  Foam::volScalarField meshField = amr->assignToVolScalarField(refinementVec);
209 
210  for (int i = 0; i < nSteps; i++) {
211  amr->updateAMRML(this->opts_.refineInterval, this->opts_.maxRefinement,
212  this->opts_.nBufferLayers, this->opts_.maxCells,
213  meshField);
214 
215  runTime++;
216  }
217 
218  finishRefine(amr.get(), this->files_.outputMeshFile, runTime);
219  std::cout << "End!" << std::endl;
220 }
221 #endif
222 
223 } // namespace DRV
224 } // namespace NEM
std::string outputMeshFile
Definition: NemDriver.H:91
void enableUpdatedField()
Access to updated field at every time step.
Definition: AMRFoam.C:830
geoMeshBase * Read(const std::string &fileName)
Read a mesh from file.
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
RefineDriver(Files files)
Definition: RefineDriver.C:34
STL namespace.
AMRFoam class is a standalone adaptive mesh refinement application that uses polyMesh as its base and...
Definition: AMRFoam.H:52
volScalarField assignToVolScalarField(const std::vector< int > &vec)
Converts vector to volScalarField for refinement/unrefinement.
Definition: AMRFoam.C:750
const Opts & getOpts() const
std::string inputMeshFile
Definition: NemDriver.H:90
void execute() const override
Run the workflow represented by the driver.
std::shared_ptr< meshBase > mesh
static std::shared_ptr< meshBase > CreateShared(const std::string &fname)
Create shared ptr from fname.
Definition: meshBase.C:171
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 writeMesh()
Ability to write current mesh into current time directory.
Definition: AMRFoam.C:704
A concrete implementation of geoMeshBase representing a mesh in a fvMesh.
Definition: foamGeoMesh.H:52