12#ifndef DUMUX_IO_VTK_OUTPUT_MODULE_HH
13#define DUMUX_IO_VTK_OUTPUT_MODULE_HH
19#include <dune/common/timer.hh>
20#include <dune/common/fvector.hh>
21#include <dune/common/typetraits.hh>
23#include <dune/geometry/type.hh>
24#include <dune/geometry/multilineargeometry.hh>
26#include <dune/grid/common/mcmgmapper.hh>
27#include <dune/grid/common/partitionset.hh>
28#include <dune/grid/io/file/vtk/vtkwriter.hh>
29#include <dune/grid/io/file/vtk/vtksequencewriter.hh>
30#include <dune/grid/common/partitionset.hh>
32#include <dumux/common/concepts/variables_.hh>
48template<
class Gr
idGeometry>
51 using GridView =
typename GridGeometry::GridView;
52 static constexpr int dim = GridView::dimension;
59 const std::string&
name,
61 Dune::VTK::DataMode dm = Dune::VTK::conforming,
69 const auto precisionString = getParamFromGroup<std::string>(
paramGroup,
"Vtk.Precision",
"Float32");
72 writer_ = std::make_shared<Dune::VTKWriter<GridView>>(
gridGeometry.gridView(), dm, coordPrecision);
73 sequenceWriter_ = std::make_unique<Dune::VTKSequenceWriter<GridView>>(writer_,
name);
74 addProcessRank_ = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank",
true);
81 {
return paramGroup_; }
92 template<
typename Vector>
94 const std::string&
name,
108 template<
typename Vector>
110 const std::string&
name,
115 const auto nComp = getNumberOfComponents_(v);
117 const auto numElemDofs =
gridGeometry().elementMapper().size();
118 const auto numVertexDofs =
gridGeometry().vertexMapper().size();
123 if(numElemDofs == numVertexDofs)
124 DUNE_THROW(Dune::InvalidStateException,
"Automatic deduction of FieldType failed. Please explicitly specify FieldType::element or FieldType::vertex.");
126 if(v.size() == numElemDofs)
128 else if(v.size() == numVertexDofs)
131 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
137 if(v.size() != numElemDofs)
138 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
141 if(v.size() != numVertexDofs)
142 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
160 for (
auto i = 0UL; i < fields_.size(); ++i)
162 if (fields_[i].
name() == field.name() && fields_[i].codim() == field.codim())
164 fields_[i] = std::move(field);
165 std::cout << Fmt::format(
166 "VtkOutputModule: Replaced field \"{}\" (codim {}). "
167 "A field by the same name & codim had already been registered previously.\n",
168 field.name(), field.codim()
175 fields_.push_back(std::move(field));
183 void write(
double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
188 if (dm_ == Dune::VTK::conforming)
190 else if (dm_ == Dune::VTK::nonconforming)
193 DUNE_THROW(Dune::NotImplemented,
"Output for provided VTK data mode");
198 std::cout << Fmt::format(
"Writing output for problem \"{}\". Took {:.2g} seconds.\n", name_, timer.elapsed());
205 const std::string&
name()
const {
return name_; }
206 Dune::VTK::DataMode
dataMode()
const {
return dm_; }
207 Dumux::Vtk::Precision
precision()
const {
return precision_; }
209 Dune::VTKWriter<GridView>&
writer() {
return *writer_; }
212 const std::vector<Field>&
fields()
const {
return fields_; }
246 std::vector<int> rank;
249 if (!fields_.empty() || addProcessRank_)
251 const auto numCells = gridGeometry_.gridView().size(0);
256 rank.resize(numCells);
258 for (
const auto& element : elements(gridGeometry_.gridView(), Dune::Partitions::interior))
260 const auto eIdxGlobal = gridGeometry_.elementMapper().index(element);
261 rank[eIdxGlobal] = gridGeometry_.gridView().comm().rank();
271 this->
addCellData(
Field(gridGeometry_.gridView(), gridGeometry_.elementMapper(), rank,
"process rank", 1, 0));
274 for (
auto&& field : fields_)
276 if (field.codim() == 0)
278 else if (field.codim() == dim)
281 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
295 this->addedCellData_.clear();
296 this->addedVertexData_.clear();
302 DUNE_THROW(Dune::NotImplemented,
"Non-conforming VTK output");
306 template<
class Vector>
307 std::size_t getNumberOfComponents_(
const Vector& v)
309 if constexpr (Dune::IsIndexable<decltype(std::declval<Vector>()[0])>::value)
315 const GridGeometry& gridGeometry_;
317 const std::string paramGroup_;
318 Dune::VTK::DataMode dm_;
320 Dumux::Vtk::Precision precision_;
322 std::shared_ptr<Dune::VTKWriter<GridView>> writer_;
323 std::unique_ptr<Dune::VTKSequenceWriter<GridView>> sequenceWriter_;
325 std::vector<Field> fields_;
327 bool addProcessRank_ =
true;
342template<
class Gr
idVariables,
class SolutionVector>
346 using GridGeometry =
typename GridVariables::GridGeometry;
348 using VV = Concept::Variables_t<GridVariables>;
349 using Scalar =
typename GridVariables::Scalar;
351 using GridView =
typename GridGeometry::GridView;
354 dim = GridView::dimension,
355 dimWorld = GridView::dimensionworld
358 using Element =
typename GridView::template Codim<0>::Entity;
359 using VolVarsVector = Dune::FieldVector<Scalar, dimWorld>;
366 struct VolVarScalarDataInfo { std::function<Scalar(
const VV&)> get; std::string
name; Dumux::Vtk::Precision precision_; };
367 struct VolVarVectorDataInfo { std::function<VolVarsVector(
const VV&)> get; std::string
name; Dumux::Vtk::Precision precision_; };
378 const SolutionVector&
sol,
379 const std::string&
name,
381 Dune::VTK::DataMode dm = Dune::VTK::conforming,
388 enableVelocityOutput_ = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddVelocity",
false);
389 addProcessRank_ = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank",
true);
410 const std::string&
name)
414 for (
auto i = 0UL; i < volVarScalarDataInfo_.size(); ++i)
416 if (volVarScalarDataInfo_[i].
name ==
name)
418 volVarScalarDataInfo_[i] = VolVarScalarDataInfo{f,
name, this->
precision()};
419 std::cout << Fmt::format(
420 "VtkOutputModule: Replaced volume variable output \"{}\". "
421 "A field by the same name had already been registered previously.\n",
429 volVarScalarDataInfo_.push_back(VolVarScalarDataInfo{f,
name, this->
precision()});
436 template<
class VVV = VolVarsVector,
typename std::enable_if_t<(VVV::dimension > 1),
int> = 0>
438 const std::string&
name)
442 for (
auto i = 0UL; i < volVarVectorDataInfo_.size(); ++i)
444 if (volVarVectorDataInfo_[i].
name ==
name)
446 volVarVectorDataInfo_[i] = VolVarVectorDataInfo{f,
name, this->
precision()};
447 std::cout << Fmt::format(
448 "VtkOutputModule: Replaced volume variable output \"{}\". "
449 "A field by the same name had already been registered previously.\n",
457 volVarVectorDataInfo_.push_back(VolVarVectorDataInfo{f,
name, this->
precision()});
464 const GridGeometry&
gridGeometry()
const {
return gridVariables_.gridGeometry(); }
465 const SolutionVector&
sol()
const {
return sol_; }
475 if constexpr (Concept::FVGridVariables<GridVariables>)
476 return gridVariables_.curGridVolVars();
478 return gridVariables_.curGridVars();
484 void writeConforming_(
double time, Dune::VTK::OutputType type)
override
486 const Dune::VTK::DataMode dm = Dune::VTK::conforming;
493 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
496 std::vector<double> rank;
499 std::vector<std::vector<Scalar>> volVarScalarData;
500 std::vector<std::vector<VolVarsVector>> volVarVectorData;
503 if (!volVarScalarDataInfo_.empty()
504 || !volVarVectorDataInfo_.empty()
505 || !this->fields().empty()
506 || velocityOutput_->enableOutput()
509 const auto numCells =
gridGeometry().gridView().size(0);
510 const auto numDofs = numDofs_();
513 if (!volVarScalarDataInfo_.empty())
514 volVarScalarData.resize(volVarScalarDataInfo_.size(), std::vector<Scalar>(numDofs));
515 if (!volVarVectorDataInfo_.empty())
516 volVarVectorData.resize(volVarVectorDataInfo_.size(), std::vector<VolVarsVector>(numDofs));
518 if (velocityOutput_->enableOutput())
520 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
523 velocity[phaseIdx].resize(numCells);
525 velocity[phaseIdx].resize(numDofs);
528 if(isBox && dim == 1)
529 velocity[phaseIdx].resize(numCells);
531 velocity[phaseIdx].resize(numDofs);
537 if (addProcessRank_) rank.resize(numCells);
541 for (
const auto& element : elements(
gridGeometry().gridView()))
543 if (!velocityOutput_->enableOutput() &&
544 element.partitionType() != Dune::PartitionType::InteriorEntity)
548 const auto eIdxGlobal =
gridGeometry().elementMapper().index(element);
551 if (velocityOutput_->enableOutput())
553 fvGeometry.bind(element);
554 elemVolVars.bind(element, fvGeometry, sol_);
558 fvGeometry.bindElement(element);
559 elemVolVars.bindElement(element, fvGeometry, sol_);
563 if (velocityOutput_->enableOutput())
565 if constexpr (Concept::FVGridVariables<GridVariables>)
567 const auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
569 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
570 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
574 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
575 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, phaseIdx);
578 else if (
element.partitionType() != Dune::PartitionType::InteriorEntity)
583 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
585 for (
const auto& scv :
scvs(fvGeometry))
587 const auto dofIdxGlobal = scv.dofIndex();
588 const auto& volVars = elemVolVars[scv];
591 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
592 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo_[i].get(volVars);
595 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
596 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo_[i].get(volVars);
602 rank[eIdxGlobal] =
static_cast<double>(
gridGeometry().gridView().comm().rank());
610 if constexpr (isBox || isPQ1Bubble || isPQ2)
612 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
614 volVarScalarDataInfo_[i].name, 1, dim, dm, this->
precision()) );
615 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
617 volVarVectorDataInfo_[i].name, dimWorld, dim, dm, this->
precision()) );
619 if constexpr (isPQ1Bubble)
621 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
623 volVarScalarDataInfo_[i].name, 1, 0,dm, this->
precision()) );
624 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
626 volVarVectorDataInfo_[i].name, dimWorld, 0,dm, this->
precision()) );
632 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
634 volVarScalarDataInfo_[i].name, 1, 0,dm, this->
precision()) );
635 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
637 volVarVectorDataInfo_[i].name, dimWorld, 0,dm, this->
precision()) );
641 if (velocityOutput_->enableOutput())
646 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
648 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
654 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
656 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
666 for (
auto&& field : this->
fields())
668 if (field.codim() == 0)
670 else if (field.codim() == dim)
673 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
687 this->addedCellData_.clear();
688 this->addedVertexData_.clear();
692 void writeNonConforming_(
double time, Dune::VTK::OutputType type)
override
694 const Dune::VTK::DataMode dm = Dune::VTK::nonconforming;
697 if(!isBox && !isDiamond)
698 DUNE_THROW(Dune::NotImplemented,
699 "Non-conforming output for discretization scheme " << GridGeometry::discMethod
707 if (enableVelocityOutput_ && !velocityOutput_->enableOutput())
708 std::cerr <<
"Warning! Velocity output was enabled in the input file"
709 <<
" but no velocity output policy was set for the VTK output module:"
710 <<
" There will be no velocity output."
711 <<
" Use the addVelocityOutput member function of the VTK output module." << std::endl;
713 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
716 std::vector<double> rank;
719 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
720 using VectorDataContainer = std::vector< std::vector<VolVarsVector> >;
721 std::vector< ScalarDataContainer > volVarScalarData;
722 std::vector< VectorDataContainer > volVarVectorData;
725 if (!volVarScalarDataInfo_.empty()
726 || !volVarVectorDataInfo_.empty()
727 || !this->fields().empty()
728 || velocityOutput_->enableOutput()
731 const auto numCells =
gridGeometry().gridView().size(0);
732 const auto outputSize = numDofs_();
735 if (!volVarScalarDataInfo_.empty())
736 volVarScalarData.resize(volVarScalarDataInfo_.size(), ScalarDataContainer(numCells));
737 if (!volVarVectorDataInfo_.empty())
738 volVarVectorData.resize(volVarVectorDataInfo_.size(), VectorDataContainer(numCells));
740 if (velocityOutput_->enableOutput())
742 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
744 if((isBox && dim == 1) || isDiamond)
745 velocity[phaseIdx].resize(numCells);
747 velocity[phaseIdx].resize(outputSize);
752 if (addProcessRank_) rank.resize(numCells);
757 for (
const auto& element : elements(
gridGeometry().gridView()))
759 if (!velocityOutput_->enableOutput() &&
760 element.partitionType() != Dune::PartitionType::InteriorEntity)
764 const auto eIdxGlobal =
gridGeometry().elementMapper().index(element);
767 if (velocityOutput_->enableOutput())
769 fvGeometry.bind(element);
770 elemVolVars.bind(element, fvGeometry, sol_);
774 fvGeometry.bindElement(element);
775 elemVolVars.bindElement(element, fvGeometry, sol_);
778 const auto numLocalDofs = fvGeometry.numScv();
780 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
781 volVarScalarData[i][eIdxGlobal].resize(numLocalDofs);
782 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
783 volVarVectorData[i][eIdxGlobal].resize(numLocalDofs);
786 if (velocityOutput_->enableOutput())
788 if constexpr (Concept::FVGridVariables<GridVariables>)
790 const auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
791 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
792 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
796 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
797 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, phaseIdx);
800 else if (
element.partitionType() != Dune::PartitionType::InteriorEntity)
805 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
807 for (
const auto& scv :
scvs(fvGeometry))
809 const auto& volVars = elemVolVars[scv];
812 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
813 volVarScalarData[i][eIdxGlobal][scv.localDofIndex()] = volVarScalarDataInfo_[i].get(volVars);
816 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
817 volVarVectorData[i][eIdxGlobal][scv.localDofIndex()] = volVarVectorDataInfo_[i].get(volVars);
823 rank[eIdxGlobal] =
static_cast<double>(
gridGeometry().gridView().comm().rank());
831 static constexpr int dofLocCodim = isDiamond ? 1 : dim;
832 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
835 volVarScalarData[i], volVarScalarDataInfo_[i].
name,
839 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
842 volVarVectorData[i], volVarVectorDataInfo_[i].
name,
843 dimWorld, dofLocCodim, dm, this->
precision()
847 if (velocityOutput_->enableOutput())
850 if (dim > 1 && !isDiamond)
851 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
854 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
855 dimWorld, dofLocCodim, dm, this->precision()
860 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
863 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
864 dimWorld, 0, dm, this->precision()
873 rank,
"process rank", 1, 0
877 for (
const auto& field : this->
fields())
879 if (field.codim() == 0)
881 else if (field.codim() == dim || field.codim() == 1)
884 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
897 this->addedCellData_.clear();
898 this->addedVertexData_.clear();
902 std::size_t numDofs_()
const
906 if constexpr (isBox || isDiamond || isPQ1Bubble || isPQ2)
912 const GridVariables& gridVariables_;
913 const SolutionVector& sol_;
915 std::vector<VolVarScalarDataInfo> volVarScalarDataInfo_;
916 std::vector<VolVarVectorDataInfo> volVarVectorDataInfo_;
918 std::shared_ptr<VelocityOutput> velocityOutput_;
919 bool enableVelocityOutput_ =
false;
920 bool addProcessRank_ =
true;
Definition: io/velocityoutput.hh:31
std::vector< Dune::FieldVector< Scalar, dimWorld > > VelocityVector
Definition: io/velocityoutput.hh:40
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:50
Dumux::Vtk::Precision precision() const
Definition: io/vtkoutputmodule.hh:207
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:80
void addCellData(const Field &field)
Definition: io/vtkoutputmodule.hh:214
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Definition: io/vtkoutputmodule.hh:183
VtkOutputModuleBase(const GridGeometry &gridGeometry, const std::string &name, const std::string ¶mGroup="", Dune::VTK::DataMode dm=Dune::VTK::conforming, bool verbose=true)
Definition: io/vtkoutputmodule.hh:58
Dune::VTK::DataMode dataMode() const
Definition: io/vtkoutputmodule.hh:206
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:209
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:210
const std::string & name() const
Definition: io/vtkoutputmodule.hh:205
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:212
virtual ~VtkOutputModuleBase()=default
virtual void writeNonConforming_(double time, Dune::VTK::OutputType type)
Assembles the fields and adds them to the writer (nonconforming output)
Definition: io/vtkoutputmodule.hh:300
std::vector< std::string > addedCellData_
Definition: io/vtkoutputmodule.hh:234
void addField(const Vector &v, const std::string &name, Dumux::Vtk::Precision precision, Vtk::FieldType fieldType=Vtk::FieldType::automatic)
Add a scalar or vector valued vtk field.
Definition: io/vtkoutputmodule.hh:109
bool verbose() const
Definition: io/vtkoutputmodule.hh:204
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:56
std::vector< std::string > addedVertexData_
Definition: io/vtkoutputmodule.hh:235
virtual void writeConforming_(double time, Dune::VTK::OutputType type)
Assembles the fields and adds them to the writer (conforming output)
Definition: io/vtkoutputmodule.hh:239
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:93
void addVertexData(const Field &field)
Definition: io/vtkoutputmodule.hh:223
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:202
void addField(Field &&field)
Add a scalar or vector valued vtk field.
Definition: io/vtkoutputmodule.hh:156
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:344
const auto & problem() const
Definition: io/vtkoutputmodule.hh:462
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:463
void addVolumeVariable(std::function< VolVarsVector(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:437
VV VolumeVariables
export type of the volume variables for the outputfields
Definition: io/vtkoutputmodule.hh:375
void addVolumeVariable(std::function< Scalar(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:409
void addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:403
const auto & curGridVariables_() const
Definition: io/vtkoutputmodule.hh:473
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:471
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:373
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:465
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:464
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)
Definition: io/vtkoutputmodule.hh:377
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:467
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:468
Vtk field types available in Dumux.
Dune style VTK functions.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
The available discretization methods in Dumux.
constexpr FCDiamond fcdiamond
Definition: method.hh:162
constexpr PQ2 pq2
Definition: method.hh:157
constexpr Box box
Definition: method.hh:156
constexpr PQ1Bubble pq1bubble
Definition: method.hh:158
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.