25#ifndef DUMUX_PNM_VTK_OUTPUT_MODULE_HH
26#define DUMUX_PNM_VTK_OUTPUT_MODULE_HH
37template<
class Gr
idVariables,
class FluxVariables,
class SolutionVector>
41 using Scalar =
typename GridVariables::Scalar;
42 using FluxVarsCache =
typename GridVariables::GridFluxVariablesCache::FluxVariablesCache;
44 struct ThroatFluxDataInfo { std::function<Scalar(
const FluxVariables&,
const FluxVarsCache&)> get; std::string
name; };
50 const SolutionVector&
sol,
51 const std::string&
name,
53 Dune::VTK::DataMode dm = Dune::VTK::conforming,
57 if constexpr (GridVariables::VolumeVariables::numFluidPhases() >= 1)
68 void addFluxVariable(std::function<Scalar(
const FluxVariables&,
const FluxVarsCache&)>&& f,
const std::string&
name)
70 throatFluxDataInfo_.push_back(ThroatFluxDataInfo{f,
name});
71 const auto numElems =
problem().gridGeometry().gridView().size(0);
72 throatFluxData_.push_back(std::vector<Scalar>(numElems));
77 void write(
double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
79 const auto gridView =
problem().gridGeometry().gridView();
80 const auto numElems = gridView.size(0);
83 for (
auto& data : throatFluxData_)
84 data.resize(numElems);
87 for (
const auto& element : elements(gridView, Dune::Partitions::interior))
91 fvElementGeometry.bindElement(element);
93 const auto eIdx =
problem().gridGeometry().elementMapper().index(element);
98 elemVolVars.bind(element, fvElementGeometry, this->
sol());
99 elemFluxVarsCache.bind(element, fvElementGeometry, elemVolVars);
102 std::size_t dataIdx = 0;
103 for (
auto&& scvf : scvfs(fvElementGeometry))
105 if (!scvf.boundary())
107 FluxVariables fluxVars;
108 fluxVars.init(
problem(), element, fvElementGeometry, elemVolVars, scvf, elemFluxVarsCache);
111 for(
auto& data : throatFluxData_)
112 data[eIdx] = throatFluxDataInfo_[dataIdx++].get(fluxVars, elemFluxVarsCache[scvf]);
121 auto clearAndShrink = [] (
auto& data)
124 data.shrink_to_fit();
127 for (
auto& data : throatFluxData_)
128 clearAndShrink(data);
135 std::vector<ThroatFluxDataInfo> throatFluxDataInfo_;
136 std::list<std::vector<Scalar>> throatFluxData_;
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Definition: discretization/porenetwork/fvelementgeometry.hh:33
Velocity output for implicit (porous media) models.
Definition: io/velocityoutput.hh:41
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:61
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:96
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Definition: io/vtkoutputmodule.hh:173
const std::string & name() const
Definition: io/vtkoutputmodule.hh:195
bool verbose() const
Definition: io/vtkoutputmodule.hh:194
void addField(const Vector &v, const std::string &name, Vtk::FieldType fieldType=Vtk::FieldType::automatic)
Add a scalar or vector valued vtk field.
Definition: io/vtkoutputmodule.hh:109
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:307
const auto & problem() const
Definition: io/vtkoutputmodule.hh:385
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:386
void addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:360
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:388
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:387
Adds vtk output fields specific to pore-network models.
Definition: pnmvtkoutputmodule.hh:39
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Gather and process all required data and write them to a vtk file.
Definition: pnmvtkoutputmodule.hh:77
VtkOutputModule(const GridVariables &gridVariables, const SolutionVector &sol, const std::string &name, const std::string ¶mGroup="", Dune::VTK::DataMode dm=Dune::VTK::conforming, bool verbose=true)
The constructor.
Definition: pnmvtkoutputmodule.hh:49
void addFluxVariable(std::function< Scalar(const FluxVariables &, const FluxVarsCache &)> &&f, const std::string &name)
Definition: pnmvtkoutputmodule.hh:68
const auto & problem() const
Return a reference to the problem.
Definition: pnmvtkoutputmodule.hh:132
Velocity output for porous media models.
A VTK output module to simplify writing dumux simulation data to VTK format.