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>
47template<
class Gr
idGeometry>
50 using GridView =
typename GridGeometry::GridView;
51 static constexpr int dim = GridView::dimension;
58 const std::string&
name,
60 Dune::VTK::DataMode dm = Dune::VTK::conforming,
68 const auto precisionString = getParamFromGroup<std::string>(
paramGroup,
"Vtk.Precision",
"Float32");
71 writer_ = std::make_shared<Dune::VTKWriter<GridView>>(
gridGeometry.gridView(), dm, coordPrecision);
72 sequenceWriter_ = std::make_unique<Dune::VTKSequenceWriter<GridView>>(writer_,
name);
73 addProcessRank_ = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank",
true);
80 {
return paramGroup_; }
91 template<
typename Vector>
93 const std::string&
name,
107 template<
typename Vector>
109 const std::string&
name,
114 const auto nComp = getNumberOfComponents_(v);
116 const auto numElemDofs =
gridGeometry().elementMapper().size();
117 const auto numVertexDofs =
gridGeometry().vertexMapper().size();
122 if(numElemDofs == numVertexDofs)
123 DUNE_THROW(Dune::InvalidStateException,
"Automatic deduction of FieldType failed. Please explicitly specify FieldType::element or FieldType::vertex.");
125 if(v.size() == numElemDofs)
127 else if(v.size() == numVertexDofs)
130 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
136 if(v.size() != numElemDofs)
137 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
140 if(v.size() != numVertexDofs)
141 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
159 for (
auto i = 0UL; i < fields_.size(); ++i)
161 if (fields_[i].
name() == field.name() && fields_[i].codim() == field.codim())
163 fields_[i] = std::move(field);
164 std::cout << Fmt::format(
165 "VtkOutputModule: Replaced field \"{}\" (codim {}). "
166 "A field by the same name & codim had already been registered previously.\n",
167 field.name(), field.codim()
174 fields_.push_back(std::move(field));
182 void write(
double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
187 if (dm_ == Dune::VTK::conforming)
189 else if (dm_ == Dune::VTK::nonconforming)
192 DUNE_THROW(Dune::NotImplemented,
"Output for provided VTK data mode");
197 std::cout << Fmt::format(
"Writing output for problem \"{}\". Took {:.2g} seconds.\n", name_, timer.elapsed());
204 const std::string&
name()
const {
return name_; }
205 Dune::VTK::DataMode
dataMode()
const {
return dm_; }
206 Dumux::Vtk::Precision
precision()
const {
return precision_; }
208 Dune::VTKWriter<GridView>&
writer() {
return *writer_; }
211 const std::vector<Field>&
fields()
const {
return fields_; }
245 std::vector<int> rank;
248 if (!fields_.empty() || addProcessRank_)
250 const auto numCells = gridGeometry_.gridView().size(0);
255 rank.resize(numCells);
257 for (
const auto& element : elements(gridGeometry_.gridView(), Dune::Partitions::interior))
259 const auto eIdxGlobal = gridGeometry_.elementMapper().index(element);
260 rank[eIdxGlobal] = gridGeometry_.gridView().comm().rank();
270 this->
addCellData(
Field(gridGeometry_.gridView(), gridGeometry_.elementMapper(), rank,
"process rank", 1, 0));
273 for (
auto&& field : fields_)
275 if (field.codim() == 0)
277 else if (field.codim() == dim)
280 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
294 this->addedCellData_.clear();
295 this->addedVertexData_.clear();
301 DUNE_THROW(Dune::NotImplemented,
"Non-conforming VTK output");
305 template<
class Vector>
306 std::size_t getNumberOfComponents_(
const Vector& v)
308 if constexpr (Dune::IsIndexable<decltype(std::declval<Vector>()[0])>::value)
314 const GridGeometry& gridGeometry_;
316 const std::string paramGroup_;
317 Dune::VTK::DataMode dm_;
319 Dumux::Vtk::Precision precision_;
321 std::shared_ptr<Dune::VTKWriter<GridView>> writer_;
322 std::unique_ptr<Dune::VTKSequenceWriter<GridView>> sequenceWriter_;
324 std::vector<Field> fields_;
326 bool addProcessRank_ =
true;
341template<
class Gr
idVariables,
class SolutionVector>
345 using GridGeometry =
typename GridVariables::GridGeometry;
347 using VV =
typename GridVariables::VolumeVariables;
348 using Scalar =
typename GridVariables::Scalar;
350 using GridView =
typename GridGeometry::GridView;
353 dim = GridView::dimension,
354 dimWorld = GridView::dimensionworld
357 using Element =
typename GridView::template Codim<0>::Entity;
358 using VolVarsVector = Dune::FieldVector<Scalar, dimWorld>;
365 struct VolVarScalarDataInfo { std::function<Scalar(
const VV&)> get; std::string
name; Dumux::Vtk::Precision precision_; };
366 struct VolVarVectorDataInfo { std::function<VolVarsVector(
const VV&)> get; std::string
name; Dumux::Vtk::Precision precision_; };
377 const SolutionVector&
sol,
378 const std::string&
name,
380 Dune::VTK::DataMode dm = Dune::VTK::conforming,
387 enableVelocityOutput_ = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddVelocity",
false);
388 addProcessRank_ = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank",
true);
409 const std::string&
name)
413 for (
auto i = 0UL; i < volVarScalarDataInfo_.size(); ++i)
415 if (volVarScalarDataInfo_[i].
name ==
name)
417 volVarScalarDataInfo_[i] = VolVarScalarDataInfo{f,
name, this->
precision()};
418 std::cout << Fmt::format(
419 "VtkOutputModule: Replaced volume variable output \"{}\". "
420 "A field by the same name had already been registered previously.\n",
428 volVarScalarDataInfo_.push_back(VolVarScalarDataInfo{f,
name, this->
precision()});
435 template<
class VVV = VolVarsVector,
typename std::enable_if_t<(VVV::dimension > 1),
int> = 0>
437 const std::string&
name)
441 for (
auto i = 0UL; i < volVarVectorDataInfo_.size(); ++i)
443 if (volVarVectorDataInfo_[i].
name ==
name)
445 volVarVectorDataInfo_[i] = VolVarVectorDataInfo{f,
name, this->
precision()};
446 std::cout << Fmt::format(
447 "VtkOutputModule: Replaced volume variable output \"{}\". "
448 "A field by the same name had already been registered previously.\n",
456 volVarVectorDataInfo_.push_back(VolVarVectorDataInfo{f,
name, this->
precision()});
461 const auto&
problem()
const {
return gridVariables_.curGridVolVars().problem(); }
463 const GridGeometry&
gridGeometry()
const {
return gridVariables_.gridGeometry(); }
464 const SolutionVector&
sol()
const {
return sol_; }
475 void writeConforming_(
double time, Dune::VTK::OutputType type)
override
477 const Dune::VTK::DataMode dm = Dune::VTK::conforming;
484 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
487 std::vector<double> rank;
490 std::vector<std::vector<Scalar>> volVarScalarData;
491 std::vector<std::vector<VolVarsVector>> volVarVectorData;
494 if (!volVarScalarDataInfo_.empty()
495 || !volVarVectorDataInfo_.empty()
496 || !this->fields().empty()
497 || velocityOutput_->enableOutput()
500 const auto numCells =
gridGeometry().gridView().size(0);
501 const auto numDofs = numDofs_();
504 if (!volVarScalarDataInfo_.empty())
505 volVarScalarData.resize(volVarScalarDataInfo_.size(), std::vector<Scalar>(numDofs));
506 if (!volVarVectorDataInfo_.empty())
507 volVarVectorData.resize(volVarVectorDataInfo_.size(), std::vector<VolVarsVector>(numDofs));
509 if (velocityOutput_->enableOutput())
511 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
514 velocity[phaseIdx].resize(numCells);
516 velocity[phaseIdx].resize(numDofs);
519 if(isBox && dim == 1)
520 velocity[phaseIdx].resize(numCells);
522 velocity[phaseIdx].resize(numDofs);
528 if (addProcessRank_) rank.resize(numCells);
531 auto elemVolVars =
localView(gridVariables_.curGridVolVars());
532 for (
const auto& element : elements(
gridGeometry().gridView()))
534 if (!velocityOutput_->enableOutput() &&
535 element.partitionType() != Dune::PartitionType::InteriorEntity)
539 const auto eIdxGlobal =
gridGeometry().elementMapper().index(element);
542 if (velocityOutput_->enableOutput())
544 fvGeometry.bind(element);
545 elemVolVars.bind(element, fvGeometry, sol_);
549 fvGeometry.bindElement(element);
550 elemVolVars.bindElement(element, fvGeometry, sol_);
554 if (velocityOutput_->enableOutput())
556 const auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
558 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
559 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
561 else if (
element.partitionType() != Dune::PartitionType::InteriorEntity)
566 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
568 for (
const auto& scv :
scvs(fvGeometry))
570 const auto dofIdxGlobal = scv.dofIndex();
571 const auto& volVars = elemVolVars[scv];
574 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
575 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo_[i].get(volVars);
578 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
579 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo_[i].get(volVars);
585 rank[eIdxGlobal] =
static_cast<double>(
gridGeometry().gridView().comm().rank());
593 if constexpr (isBox || isPQ1Bubble || isPQ2)
595 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
597 volVarScalarDataInfo_[i].name, 1, dim, dm, this->
precision()) );
598 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
600 volVarVectorDataInfo_[i].name, dimWorld, dim, dm, this->
precision()) );
602 if constexpr (isPQ1Bubble)
604 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
606 volVarScalarDataInfo_[i].name, 1, 0,dm, this->
precision()) );
607 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
609 volVarVectorDataInfo_[i].name, dimWorld, 0,dm, this->
precision()) );
615 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
617 volVarScalarDataInfo_[i].name, 1, 0,dm, this->
precision()) );
618 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
620 volVarVectorDataInfo_[i].name, dimWorld, 0,dm, this->
precision()) );
624 if (velocityOutput_->enableOutput())
629 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
631 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
637 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
639 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
649 for (
auto&& field : this->
fields())
651 if (field.codim() == 0)
653 else if (field.codim() == dim)
656 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
670 this->addedCellData_.clear();
671 this->addedVertexData_.clear();
675 void writeNonConforming_(
double time, Dune::VTK::OutputType type)
override
677 const Dune::VTK::DataMode dm = Dune::VTK::nonconforming;
680 if(!isBox && !isDiamond)
681 DUNE_THROW(Dune::NotImplemented,
682 "Non-conforming output for discretization scheme " << GridGeometry::discMethod
690 if (enableVelocityOutput_ && !velocityOutput_->enableOutput())
691 std::cerr <<
"Warning! Velocity output was enabled in the input file"
692 <<
" but no velocity output policy was set for the VTK output module:"
693 <<
" There will be no velocity output."
694 <<
" Use the addVelocityOutput member function of the VTK output module." << std::endl;
696 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
699 std::vector<double> rank;
702 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
703 using VectorDataContainer = std::vector< std::vector<VolVarsVector> >;
704 std::vector< ScalarDataContainer > volVarScalarData;
705 std::vector< VectorDataContainer > volVarVectorData;
708 if (!volVarScalarDataInfo_.empty()
709 || !volVarVectorDataInfo_.empty()
710 || !this->fields().empty()
711 || velocityOutput_->enableOutput()
714 const auto numCells =
gridGeometry().gridView().size(0);
715 const auto outputSize = numDofs_();
718 if (!volVarScalarDataInfo_.empty())
719 volVarScalarData.resize(volVarScalarDataInfo_.size(), ScalarDataContainer(numCells));
720 if (!volVarVectorDataInfo_.empty())
721 volVarVectorData.resize(volVarVectorDataInfo_.size(), VectorDataContainer(numCells));
723 if (velocityOutput_->enableOutput())
725 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
727 if((isBox && dim == 1) || isDiamond)
728 velocity[phaseIdx].resize(numCells);
730 velocity[phaseIdx].resize(outputSize);
735 if (addProcessRank_) rank.resize(numCells);
739 auto elemVolVars =
localView(gridVariables_.curGridVolVars());
740 for (
const auto& element : elements(
gridGeometry().gridView()))
742 if (!velocityOutput_->enableOutput() &&
743 element.partitionType() != Dune::PartitionType::InteriorEntity)
747 const auto eIdxGlobal =
gridGeometry().elementMapper().index(element);
750 if (velocityOutput_->enableOutput())
752 fvGeometry.bind(element);
753 elemVolVars.bind(element, fvGeometry, sol_);
757 fvGeometry.bindElement(element);
758 elemVolVars.bindElement(element, fvGeometry, sol_);
761 const auto numLocalDofs = fvGeometry.numScv();
763 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
764 volVarScalarData[i][eIdxGlobal].resize(numLocalDofs);
765 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
766 volVarVectorData[i][eIdxGlobal].resize(numLocalDofs);
769 if (velocityOutput_->enableOutput())
771 const auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
772 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
773 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
775 else if (
element.partitionType() != Dune::PartitionType::InteriorEntity)
780 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
782 for (
const auto& scv :
scvs(fvGeometry))
784 const auto& volVars = elemVolVars[scv];
787 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
788 volVarScalarData[i][eIdxGlobal][scv.localDofIndex()] = volVarScalarDataInfo_[i].get(volVars);
791 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
792 volVarVectorData[i][eIdxGlobal][scv.localDofIndex()] = volVarVectorDataInfo_[i].get(volVars);
798 rank[eIdxGlobal] =
static_cast<double>(
gridGeometry().gridView().comm().rank());
806 static constexpr int dofLocCodim = isDiamond ? 1 : dim;
807 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
810 volVarScalarData[i], volVarScalarDataInfo_[i].
name,
814 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
817 volVarVectorData[i], volVarVectorDataInfo_[i].
name,
818 dimWorld, dofLocCodim, dm, this->
precision()
822 if (velocityOutput_->enableOutput())
825 if (dim > 1 && !isDiamond)
826 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
829 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
830 dimWorld, dofLocCodim, dm, this->precision()
835 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
838 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
839 dimWorld, 0, dm, this->precision()
848 rank,
"process rank", 1, 0
852 for (
const auto& field : this->
fields())
854 if (field.codim() == 0)
856 else if (field.codim() == dim || field.codim() == 1)
859 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
872 this->addedCellData_.clear();
873 this->addedVertexData_.clear();
877 std::size_t numDofs_()
const
881 if constexpr (isBox || isDiamond || isPQ1Bubble || isPQ2)
887 const GridVariables& gridVariables_;
888 const SolutionVector& sol_;
890 std::vector<VolVarScalarDataInfo> volVarScalarDataInfo_;
891 std::vector<VolVarVectorDataInfo> volVarVectorDataInfo_;
893 std::shared_ptr<VelocityOutput> velocityOutput_;
894 bool enableVelocityOutput_ =
false;
895 bool addProcessRank_ =
true;
Velocity output for implicit (porous media) models.
Definition: io/velocityoutput.hh:29
std::vector< Dune::FieldVector< Scalar, dimWorld > > VelocityVector
Definition: io/velocityoutput.hh:38
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:49
Dumux::Vtk::Precision precision() const
Definition: io/vtkoutputmodule.hh:206
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:79
void addCellData(const Field &field)
Definition: io/vtkoutputmodule.hh:213
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Definition: io/vtkoutputmodule.hh:182
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:57
Dune::VTK::DataMode dataMode() const
Definition: io/vtkoutputmodule.hh:205
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:208
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:209
const std::string & name() const
Definition: io/vtkoutputmodule.hh:204
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:211
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:299
std::vector< std::string > addedCellData_
Definition: io/vtkoutputmodule.hh:233
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:108
bool verbose() const
Definition: io/vtkoutputmodule.hh:203
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:55
std::vector< std::string > addedVertexData_
Definition: io/vtkoutputmodule.hh:234
virtual void writeConforming_(double time, Dune::VTK::OutputType type)
Assembles the fields and adds them to the writer (conforming output)
Definition: io/vtkoutputmodule.hh:238
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:92
void addVertexData(const Field &field)
Definition: io/vtkoutputmodule.hh:222
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:201
void addField(Field &&field)
Add a scalar or vector valued vtk field.
Definition: io/vtkoutputmodule.hh:155
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:343
const auto & problem() const
Definition: io/vtkoutputmodule.hh:461
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:462
void addVolumeVariable(std::function< VolVarsVector(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:436
VV VolumeVariables
export type of the volume variables for the outputfields
Definition: io/vtkoutputmodule.hh:374
void addVolumeVariable(std::function< Scalar(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:408
void addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:402
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:470
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:372
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:464
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:463
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:376
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:466
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:467
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.