AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
MeshCreator.hpp
1 #pragma once
2 
3 #include <algorithm>
4 #include <array>
5 #include <memory>
6 #include <string>
7 
8 #include <dune/common/filledarray.hh>
9 #include <dune/common/fvector.hh>
10 #include <dune/common/typeutilities.hh>
11 
12 #if HAVE_ALBERTA
13 #include <dune/grid/albertagrid/albertareader.hh>
14 #endif
15 
16 #include <dune/grid/io/file/gmshreader.hh>
17 #include <dune/grid/utility/structuredgridfactory.hh>
18 
19 #if HAVE_DUNE_VTK
20 #include <dune/vtk/vtkreader.hh>
21 #include <dune/vtk/gridcreators/lagrangegridcreator.hh>
22 #endif
23 
24 #if HAVE_DUNE_GMSH4
25 #include <dune/gmsh4/gmsh4reader.hh>
26 #include <dune/gmsh4/gridcreators/lagrangegridcreator.hh>
27 #include <dune/gmsh4/utility/version.hh>
28 #endif
29 
30 #include <amdis/AdaptiveGrid.hpp>
31 #include <amdis/Initfile.hpp>
32 #include <amdis/Output.hpp>
33 #include <amdis/common/Concepts.hpp>
34 #include <amdis/common/Filesystem.hpp>
35 #include <amdis/common/TypeTraits.hpp>
36 #include <amdis/utility/MacroGridFactory.hpp>
37 
38 namespace Dune
39 {
40  // forward declarations
41  template <int dim, class Coordinates> class YaspGrid;
42  template <class ct, int dim> class EquidistantCoordinates;
43  template <class ct, int dim> class EquidistantOffsetCoordinates;
44  template <class ct, int dim> class TensorProductCoordinates;
45 }
46 
47 
48 namespace AMDiS
49 {
51  template <class G>
52  struct MeshCreator
53  {
54  enum { dimension = G::dimension };
55  enum { dimworld = G::dimensionworld };
56 
57  using Grid = AdaptiveGrid_t<G>;
58  using HostGrid = typename Grid::HostGrid;
59 
60  using ctype = typename Grid::ctype;
61 
63 
66  MeshCreator(std::string const& name)
67  : name_(name)
68  {}
69 
71  static std::shared_ptr<Grid> create(std::string name)
72  {
73  return MeshCreator{name}.create();
74  }
75 
77 
88  std::unique_ptr<Grid> create() const
89  {
90  auto filename = Parameters::get<std::string>(name_ + "->macro file name");
91  auto structured = Parameters::get<std::string>(name_ + "->structured");
92  std::unique_ptr<HostGrid> gridPtr;
93  if (filename) {
94  // read a macro file
95  gridPtr = create_unstructured_grid(filename.value());
96  } else if (structured) {
97  // create structured grids
98  if (structured.value() == "cube") {
99  gridPtr = create_cube_grid();
100  } else if (structured && structured.value() == "simplex") {
101  gridPtr = create_simplex_grid();
102  } else {
103  error_exit("Unknown structured grid type. Must be either `cube` or `simplex` in parameter [meshName]->structured.");
104  }
105  } else {
106  // decide by inspecting the grid type how to create the grid
107  gridPtr = create_by_gridtype<HostGrid>(Dune::PriorityTag<42>{});
108  }
109 
110  // Perform initial refinement and load balance if requested in the initfile
111  auto globalRefinements = Parameters::get<int>(name_ + "->global refinements");
112  if (globalRefinements)
113  gridPtr->globalRefine(globalRefinements.value());
114  auto loadBalance = Parameters::get<bool>(name_ + "->load balance");
115  if (loadBalance && loadBalance.value())
116  gridPtr->loadBalance();
117 
118  return construct(std::move(gridPtr));
119  }
120 
122  std::unique_ptr<HostGrid> create_cube_grid() const
123  {
124  return create_structured_grid([](auto&& lower, auto&& upper, auto&& numCells)
125  {
126  return Dune::StructuredGridFactory<HostGrid>::createCubeGrid(lower, upper, numCells);
127  });
128  }
129 
131  std::unique_ptr<HostGrid> create_simplex_grid() const
132  {
133  return create_structured_grid([](auto&& lower, auto&& upper, auto&& numCells)
134  {
135  return Dune::StructuredGridFactory<HostGrid>::createSimplexGrid(lower, upper, numCells);
136  });
137  }
138 
140  std::vector<int> const& boundaryIds() const
141  {
142  return boundaryIds_;
143  }
144 
146  std::vector<int> const& elementIds() const
147  {
148  return elementIds_;
149  }
150 
151  private:
152  static std::unique_ptr<Grid> construct(std::unique_ptr<Grid> hostGrid)
153  {
154  return std::move(hostGrid);
155  }
156 
157  template <class HG>
158  static std::unique_ptr<Grid> construct(std::unique_ptr<HG> hostGrid)
159  {
160  return std::make_unique<Grid>(std::move(hostGrid));
161  }
162 
163  // use the structured grid factory to create the grid
164  template <class Size = unsigned int, class Factory>
165  std::unique_ptr<HostGrid> create_structured_grid(Factory factory) const
166  {
167  if constexpr (int(Grid::dimension) == int(Grid::dimensionworld)) {
168  // Lower left corner of the domain
169  Dune::FieldVector<ctype,int(dimworld)> lower(0);
170  // Upper right corner of the domain
171  Dune::FieldVector<ctype,int(dimworld)> upper(1);
172  // number of blocks in each coordinate direction
173  auto numCells = Dune::filledArray<std::size_t(dimension),Size>(1);
174 
175  Parameters::get(name_ + "->min corner", lower);
176  Parameters::get(name_ + "->max corner", upper);
177  Parameters::get(name_ + "->num cells", numCells);
178  return factory(lower, upper, numCells);
179  } else {
180  error_exit("Structured grids can be created for dim == dow only.");
181  return nullptr;
182  }
183  }
184 
185  // read a filename from `[meshName]->macro file name` and determine from the extension the fileformat
186  std::unique_ptr<HostGrid> create_unstructured_grid(std::string const& filename) const
187  {
188  filesystem::path fn(filename);
189  auto ext = fn.extension();
190 
191  if (ext == ".msh") {
192 #if HAVE_DUNE_GMSH4
193  if (Dune::Gmsh4::fileVersion(filename)[0] >= 4)
194  return Dune::Gmsh4Reader<HostGrid, Dune::Gmsh4::LagrangeGridCreator<HostGrid>>::createGridFromFile(filename);
195  else
196 #endif
197  return read_gmsh_file<HostGrid>(filename, Dune::PriorityTag<42>{});
198  }
199 #if HAVE_DUNE_VTK
200  else if (ext == ".vtu") {
201  return Dune::VtkReader<HostGrid, Dune::Vtk::LagrangeGridCreator<HostGrid>>::createGridFromFile(filename);
202  }
203 #endif
204  else if (ext == ".1d" || ext == ".2d" || ext == ".3d" || ext == ".amc") {
205  return read_alberta_file<HostGrid>(filename, Dune::PriorityTag<42>{});
206  }
207  else {
208  error_exit("Cannot read grid file. Unknown file extension.");
209  return {};
210  }
211  }
212 
213  template <class GridType>
214  using SupportsGmshReader = decltype(std::declval<Dune::GridFactory<GridType>>().insertBoundarySegment(
215  std::declval<std::vector<unsigned int>>(),
216  std::declval<std::shared_ptr<Dune::BoundarySegment<GridType::dimension, GridType::dimensionworld> >>()) );
217 
218  template <class GridType, class LC>
219  using SupportsInsertionIndex = decltype(std::declval<Dune::GridFactory<GridType>>().insertionIndex(std::declval<LC>()));
220 
221  // use GmshReader if GridFactory supports insertBoundarySegments
222  template <class GridType = HostGrid,
223  REQUIRES(Dune::Std::is_detected<SupportsGmshReader, GridType>::value)>
224  std::unique_ptr<GridType> read_gmsh_file(std::string const& filename, Dune::PriorityTag<1>) const
225  {
226  Dune::GmshReader<GridType> reader;
227  Dune::GridFactory<GridType> factory;
228  std::vector<int> boundaryIds, elementIds;
229  reader.read(factory, filename, boundaryIds, elementIds);
230 
231  using HasInsertionIndexEntity
232  = Dune::Std::is_detected<SupportsInsertionIndex, GridType, typename GridType::template Codim<0>::Entity>;
233  using HasInsertionIndexIntersection
234  = Dune::Std::is_detected<SupportsInsertionIndex, GridType, typename GridType::LeafIntersection>;
235 
236  auto gridPtr = factory.createGrid();
237  if ((boundaryIds.empty() && elementIds.empty()) ||
238  (!HasInsertionIndexEntity::value && !HasInsertionIndexIntersection::value))
239  return std::unique_ptr<GridType>(std::move(gridPtr));
240 
241  // map boundaryIds and elementIds read from file to grid indexing.
242  if (!boundaryIds.empty() && HasInsertionIndexIntersection::value)
243  boundaryIds_.resize(gridPtr->numBoundarySegments());
244  if (!elementIds.empty() && HasInsertionIndexEntity::value)
245  elementIds_.resize(gridPtr->size(0));
246 
247  auto const& indexSet = gridPtr->leafIndexSet();
248  for (auto const& e : elements(gridPtr->leafGridView())) {
249  if constexpr(HasInsertionIndexEntity::value) {
250  if (!elementIds.empty())
251  elementIds_[indexSet.index(e)] = elementIds[factory.insertionIndex(e)];
252  }
253 
254  if (boundaryIds.empty() || !e.hasBoundaryIntersections())
255  continue;
256 
257  for (auto const& it : intersections(gridPtr->leafGridView(), e)) {
258  if constexpr(HasInsertionIndexIntersection::value) {
259  if (it.boundary())
260  boundaryIds_[it.boundarySegmentIndex()]
261  = boundaryIds[factory.insertionIndex(it)];
262  }
263  }
264  }
265 
266  return std::unique_ptr<GridType>(std::move(gridPtr));
267  }
268 
269  // fallback if GmshReader cannot be used
270  template <class GridType = HostGrid>
271  std::unique_ptr<GridType> read_gmsh_file(std::string const&, Dune::PriorityTag<0>) const
272  {
273  error_exit("Gmsh reader not supported for this grid type!");
274  return {};
275  }
276 
277  // read from Alberta file
278 
279 #if HAVE_ALBERTA
280  template <class GridType>
281  using IsAlbertaGrid = decltype(std::declval<GridType>().alberta2dune(0,0));
282 
283  // construct the albertagrid directly from a filename
284  template <class GridType = HostGrid,
285  REQUIRES(GridType::dimensionworld == DIM_OF_WORLD),
286  REQUIRES(Dune::Std::is_detected<IsAlbertaGrid, GridType>::value)>
287  std::unique_ptr<GridType> read_alberta_file(std::string const& filename, Dune::PriorityTag<3>) const
288  {
289  return std::make_unique<GridType>(filename);
290  }
291 
292  // use a gridfactory and the generic AlbertaReader
293  template <class GridType = HostGrid,
294  REQUIRES(GridType::dimensionworld == DIM_OF_WORLD)>
295  std::unique_ptr<GridType> read_alberta_file(std::string const& filename, Dune::PriorityTag<2>) const
296  {
297  Dune::GridFactory<GridType> factory;
298  if (Environment::mpiRank() == 0) {
299  Dune::AlbertaReader<GridType> reader;
300  reader.readGrid(filename, factory);
301  }
302  return std::unique_ptr<GridType>{factory.createGrid()};
303  }
304 
305  // error if WORLDDIM not the same as Grid::dimensionworld
306  template <class GridType = HostGrid,
307  REQUIRES(GridType::dimensionworld != DIM_OF_WORLD)>
308  std::unique_ptr<GridType> read_alberta_file(std::string const& filename, Dune::PriorityTag<1>) const
309  {
310  error_exit("AlbertaGrid (and AlbertaReader) require WORLDDIM == Grid::dimensionworld. Change the cmake parameters!");
311  return {};
312  }
313 #else
314  // fallback if no Alberta is available
315  template <class GridType>
316  std::unique_ptr<GridType> read_alberta_file(std::string const&, Dune::PriorityTag<0>) const
317  {
318  error_exit("AlbertaGrid (and AlbertaReader) not available. Set AlbertaFlags to your executable in cmake!");
319  return {};
320  }
321 #endif
322 
323 
324 #if HAVE_ALBERTA
325  // albertagrid -> simplex
326  template <class GridType = HostGrid,
327  REQUIRES(Dune::Std::is_detected<IsAlbertaGrid, GridType>::value)>
328  std::unique_ptr<GridType> create_by_gridtype(Dune::PriorityTag<3>) const
329  {
330  return create_structured_grid([](auto&& lower, auto&& upper, auto&& numCells)
331  {
332  return MacroGridFactory<GridType>::createSimplexGrid(lower, upper, numCells);
333  });
334  }
335 #endif
336 
337  // yasp grid -> cube
338  template <class GridType = HostGrid,
339  class = typename GridType::YGrid>
340  std::unique_ptr<GridType> create_by_gridtype(Dune::PriorityTag<2>) const
341  {
342  int overlap = 0;
343  Parameters::get(name_ + "->overlap", overlap);
344  std::string periodic(dimension, '0');
345  Parameters::get(name_ + "->periodic", periodic); // e.g. 000 01 111
346 
347  return create_yaspgrid(Types<GridType>{}, overlap, std::bitset<dimension>(periodic));
348  }
349 
350  template <int dim, class ct>
351  std::unique_ptr<HostGrid> create_yaspgrid(Types<Dune::YaspGrid<dim,Dune::EquidistantCoordinates<ct,dim>>>,
352  int overlap, std::bitset<dimension> const& per) const
353  {
354  return create_structured_grid<int>([&](auto&& /*lower*/, auto&& upper, std::array<int,dimension> const& numCells)
355  {
356  return std::make_unique<HostGrid>(upper, numCells, per, overlap);
357  });
358  }
359 
360  template <int dim, class ct>
361  std::unique_ptr<HostGrid> create_yaspgrid(Types<Dune::YaspGrid<dim,Dune::EquidistantOffsetCoordinates<ct,dim>>>,
362  int overlap, std::bitset<dimension> const& per) const
363  {
364  return create_structured_grid<int>([&](auto&& lower, auto&& upper, std::array<int,dimension> const& numCells)
365  {
366  return std::make_unique<HostGrid>(lower, upper, numCells, per, overlap);
367  });
368  }
369 
370  template <int dim, class ct>
371  std::unique_ptr<HostGrid> create_yaspgrid(Types<Dune::YaspGrid<dim,Dune::TensorProductCoordinates<ct,dim>>>,
372  int, std::bitset<dimension> const&) const
373  {
374  error_exit("MeshCreator cannot create YaspGrid with TensorProductCoordinates.");
375  return {};
376  }
377 
378 
379 #if HAVE_DUNE_SPGRID
380  // spgrid -> cube
381  template <class GridType = HostGrid,
382  class = typename GridType::ReferenceCube,
383  class = typename GridType::MultiIndex>
384  std::unique_ptr<GridType> create_by_gridtype(Dune::PriorityTag<1>) const
385  {
386  return create_structured_grid<int>([](auto&& lower, auto&& upper, std::array<int,dimension> const& numCells)
387  {
388  typename GridType::MultiIndex multiIndex(numCells);
389  return std::make_unique<GridType>(lower, upper, multiIndex);
390  });
391  }
392 #endif
393 
394  // final fallback
395  template <class GridType = HostGrid>
396  std::unique_ptr<GridType> create_by_gridtype(Dune::PriorityTag<0>) const
397  {
398  error_exit("Don't know how to create the grid.");
399  return {};
400  }
401 
402  private:
403  std::string name_;
404  mutable std::vector<int> boundaryIds_;
405  mutable std::vector<int> elementIds_;
406  };
407 
408 
409 } // end namespace AMDiS
constexpr bool MultiIndex
A multi-index type.
Definition: Concepts.hpp:150
path extension() const
Returns the file extension path component.
Definition: Filesystem.cpp:83
Definition: AdaptiveGrid.hpp:373
A variadic type list.
Definition: TypeTraits.hpp:128
Definition: AdaptBase.hpp:6
std::unique_ptr< HostGrid > create_simplex_grid() const
Create a structured simplex grid.
Definition: MeshCreator.hpp:131
std::vector< int > const & elementIds() const
Return the filled vector of element ids. NOTE: not all creators support reading this.
Definition: MeshCreator.hpp:146
A creator class for dune grids.
Definition: MeshCreator.hpp:52
static std::shared_ptr< Grid > create(std::string name)
Static create mthod. See create()
Definition: MeshCreator.hpp:71
std::unique_ptr< HostGrid > create_cube_grid() const
Create a structured cube grid.
Definition: MeshCreator.hpp:122
Definition: Filesystem.hpp:11
MeshCreator(std::string const &name)
Construct a new MeshCreator.
Definition: MeshCreator.hpp:66
std::unique_ptr< Grid > create() const
Create a new grid by inspecting the intifile parameter group [meshName]
Definition: MeshCreator.hpp:88
std::vector< int > const & boundaryIds() const
Return the filled vector of boundary ids. NOTE: not all creators support reading this.
Definition: MeshCreator.hpp:140
Definition: MacroGridFactory.hpp:97