Skip to content
Snippets Groups Projects
Commit 88703a1a authored by Timo Koch's avatar Timo Koch
Browse files

[io][gridmanager] Add support for dune-spgrid

parent 179cbd53
No related branches found
No related tags found
1 merge request!1097Feature/periodic boundaries
...@@ -67,6 +67,12 @@ ...@@ -67,6 +67,12 @@
#include <dune/foamgrid/dgffoam.hh> #include <dune/foamgrid/dgffoam.hh>
#endif #endif
// SPGrid specific includes
#if HAVE_DUNE_SPGRID
#include <dune/grid/spgrid.hh>
#include <dune/grid/spgrid/dgfparser.hh>
#endif
#include <dumux/common/parameters.hh> #include <dumux/common/parameters.hh>
#include <dumux/discretization/methods.hh> #include <dumux/discretization/methods.hh>
...@@ -1281,6 +1287,69 @@ public: ...@@ -1281,6 +1287,69 @@ public:
#endif // HAVE_DUNE_FOAMGRID #endif // HAVE_DUNE_FOAMGRID
#if HAVE_DUNE_SPGRID
/*!
* \brief Provides a grid manager for SPGrid
*
* The following keys are recognized:
* - File : A DGF or gmsh file to load from, type detection by file extension
*
*/
template<class ct, int dim, template< int > class Ref, class Comm>
class GridManager<Dune::SPGrid<ct, dim, Ref, Comm>>
: public GridManagerBase<Dune::SPGrid<ct, dim, Ref, Comm>>
{
public:
using Grid = Dune::SPGrid<ct, dim, Ref, Comm>;
using ParentType = GridManagerBase<Grid>;
/*!
* \brief Make the grid. This is implemented by specializations of this method.
*/
void init(const std::string& paramGroup = "")
{
// try to create it from file
if (hasParamInGroup(paramGroup, "Grid.File"))
{
ParentType::makeGridFromDgfFile(getParamFromGroup<std::string>(paramGroup, "Grid.File"));
ParentType::maybeRefineGrid(paramGroup);
ParentType::loadBalance();
return;
}
// Didn't find a way to construct the grid
else if (hasParamInGroup(paramGroup, "Grid.UpperRight"))
{
using GlobalPosition = Dune::FieldVector<ct, dim>;
const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.LowerLeft", GlobalPosition(0.0));
const auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.UpperRight");
using IntArray = std::array<int, dim>;
IntArray cells; cells.fill(1);
cells = getParamFromGroup<IntArray>(paramGroup, "Grid.Cells", cells);
const auto periodic = getParamFromGroup<std::bitset<dim>>(paramGroup, "Grid.Periodic", std::bitset<dim>{});
const auto overlap = getParamFromGroup<int>(paramGroup, "Grid.Overlap", 1);
IntArray spOverlap; spOverlap.fill(overlap);
using Domain = typename Grid::Domain;
std::vector< typename Domain::Cube > cubes;
cubes.push_back( typename Domain::Cube( lowerLeft, upperRight ) );
Domain domain( cubes, typename Domain::Topology( static_cast<unsigned int>(periodic.to_ulong()) ) );
ParentType::gridPtr() = std::make_shared<Grid>( domain, cells, spOverlap );
ParentType::maybeRefineGrid(paramGroup);
ParentType::loadBalance();
}
else
{
const auto prefix = paramGroup == "" ? paramGroup : paramGroup + ".";
DUNE_THROW(ParameterException, "Please supply a grid file in " << prefix << "Grid.File or " << prefix << "Grid.UpperRight/Cells.");
}
}
};
#endif // HAVE_DUNE_SPGRID
} // end namespace Dumux } // end namespace Dumux
#endif #endif
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment