24#ifndef DUMUX_STAGGERED_FF_VELOCITYOUTPUT_HH
25#define DUMUX_STAGGERED_FF_VELOCITYOUTPUT_HH
36template<
class Gr
idVariables,
class SolutionVector>
40 using GridGeometry =
typename GridVariables::GridGeometry;
41 using Scalar =
typename GridVariables::Scalar;
42 using FVElementGeometry =
typename GridGeometry::LocalView;
43 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
44 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
45 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
46 using ElementFluxVarsCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
47 using VolumeVariables =
typename GridVariables::VolumeVariables;
48 using FluidSystem =
typename VolumeVariables::FluidSystem;
49 using GridView =
typename GridGeometry::GridView;
51 using Problem =
typename std::decay_t<decltype(std::declval<GridVolumeVariables>().problem())>;
52 using Element =
typename GridView::template Codim<0>::Entity;
53 using CoordScalar =
typename GridView::ctype;
65 : gridVariables_(gridVariables)
70 enableOutput_ = getParamFromGroup<bool>(gridVariables.curGridVolVars().problem().paramGroup(),
"Vtk.AddVelocity",
true);
77 std::string
phaseName(
int phaseIdx)
const override {
return FluidSystem::phaseName(phaseIdx); }
80 int numFluidPhases()
const override {
return VolumeVariables::numFluidPhases(); }
85 const Element& element,
86 const FVElementGeometry& fvGeometry,
87 const ElementVolumeVariables& elemVolVars,
88 const ElementFluxVarsCache& elemFluxVarsCache,
89 int phaseIdx)
const override
91 auto elemFaceVars =
localView(gridVariables_.curGridFaceVars());
92 elemFaceVars.bindElement(element, fvGeometry, sol_);
93 for (
auto&& scv : scvs(fvGeometry))
95 auto dofIdxGlobal = scv.dofIndex();
97 for (
auto&& scvf : scvfs(fvGeometry))
99 auto dirIdx = scvf.directionIndex();
100 velocity[dofIdxGlobal][dirIdx] += 0.5*elemFaceVars[scvf].velocitySelf();
106 const GridVariables& gridVariables_;
107 const SolutionVector& sol_;
109 bool enableOutput_ =
true;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Velocity output for staggered free-flow models.
Definition: discretization/staggered/freeflow/velocityoutput.hh:38
void calculateVelocity(VelocityVector &velocity, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVarsCache &elemFluxVarsCache, int phaseIdx) const override
Definition: discretization/staggered/freeflow/velocityoutput.hh:84
StaggeredFreeFlowVelocityOutput(const GridVariables &gridVariables, const SolutionVector &sol)
Constructor initializes the static data with the initial solution.
Definition: discretization/staggered/freeflow/velocityoutput.hh:64
bool enableOutput() const override
Returns whether to enable the velocity output or not.
Definition: discretization/staggered/freeflow/velocityoutput.hh:74
std::string phaseName(int phaseIdx) const override
returns the phase name of a given phase index
Definition: discretization/staggered/freeflow/velocityoutput.hh:77
int numFluidPhases() const override
returns the number of phases
Definition: discretization/staggered/freeflow/velocityoutput.hh:80
typename ParentType::VelocityVector VelocityVector
Definition: discretization/staggered/freeflow/velocityoutput.hh:56
Velocity output for implicit (porous media) models.
Definition: io/velocityoutput.hh:41
std::vector< Dune::FieldVector< Scalar, dimWorld > > VelocityVector
Definition: io/velocityoutput.hh:50
Default velocity output policy for porous media models.