Skip to content
Snippets Groups Projects
Commit 8c58db99 authored by Kilian Weishaupt's avatar Kilian Weishaupt
Browse files

[freeflow][iofields] Add createStaggeredFVNameFunction as free function

* takes a multitypeblockvector and an index constant to automatically
  figure out the pv names for cellcenter and face privars
parent 4aa41755
No related branches found
No related tags found
1 merge request!1212Feature/iofields
......@@ -26,6 +26,8 @@
#include <dune/common/fvector.hh>
#include <dune/common/deprecated.hh>
#include <dune/common/indices.hh>
#include <dune/istl/multitypeblockvector.hh> // TODO: needed? or is forward declare enough?
#include <dumux/common/parameters.hh>
#include <dumux/discretization/methods.hh>
......@@ -34,6 +36,44 @@
namespace Dumux
{
/*!
* \ingroup InputOutput
* \ingroup NavierStokesModel
* \brief helper function to determine the names of cell-centered primary variables of a model with staggered grid discretization
* \note use this as input for the load solution function
*/
template<class ModelTraits, class IOFields, class FluidSystem, std::size_t cellCenterIdx = 0, class ...Args, std::size_t i>
std::function<std::string(int,int)> createStaggeredPVNameFunction(const Dune::MultiTypeBlockVector<Args...>&,
Dune::index_constant<i> idx,
const std::string& paramGroup = "")
{
using CellCenterPriVarsType = typename std::tuple_element_t<cellCenterIdx, std::tuple<Args...>>::block_type;
static constexpr auto offset = ModelTraits::numEq() - CellCenterPriVarsType::dimension;
if (i == cellCenterIdx) // cell center
{
if (hasParamInGroup(paramGroup, "LoadSolution.CellCenterPriVarNames"))
{
const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, "LoadSolution.CellCenterPriVarNames");
return [n = std::move(pvName)](int pvIdx, int state = 0){ return n[pvIdx]; };
}
else
// add an offset to the pvIdx so that the velocities are skipped
return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits ,FluidSystem>(pvIdx + offset, state); };
}
else // face
{
if (hasParamInGroup(paramGroup, "LoadSolution.FacePriVarNames"))
{
const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, "LoadSolution.FacePriVarNames");
return [n = std::move(pvName)](int pvIdx, int state = 0){ return n[pvIdx]; };
}
else
// there is only one scalar-type primary variable called "v" living on the faces
return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem>(pvIdx , state); };
}
}
/*!
* \ingroup NavierStokesModel
* \brief Adds I/O fields for the Navier-Stokes model
......
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