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

Merge branch 'feature/spgrid-initialize-simple' into 'master'

[spgrid] Add gridmanager.init(...) with lowerLeft/upperRight/cells

See merge request !3372
parents a8bf33c3 0bfb46ee
No related branches found
No related tags found
1 merge request!3372[spgrid] Add gridmanager.init(...) with lowerLeft/upperRight/cells
Pipeline #26447 passed
......@@ -83,15 +83,7 @@ public:
cells = getParamFromGroup<IntArray>(paramGroup, "Grid.Cells", cells);
const auto periodic = getParamFromGroup<std::bitset<dim>>(paramGroup, "Grid.Periodic", std::bitset<dim>{});
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();
init(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
}
else
{
......@@ -99,6 +91,26 @@ public:
DUNE_THROW(ParameterException, "Please supply a grid file in " << prefix << "Grid.File or " << prefix << "Grid.UpperRight/Cells.");
}
}
void init(const Dune::FieldVector<ct, dim>& lowerLeft,
const Dune::FieldVector<ct, dim>& upperRight,
const std::array<int, dim>& cells,
const std::string& paramGroup = "",
const int overlap = 1,
const std::bitset<dim> periodic = std::bitset<dim>{})
{
if (overlap == 0)
DUNE_THROW(Dune::NotImplemented, "dune-spgrid does currently not support zero overlap!");
using IntArray = std::array<int, dim>;
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();
}
};
#endif // HAVE_DUNE_SPGRID
......
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