24#ifndef DUMUX_IO_VTK_OUTPUT_MODULE_HH
25#define DUMUX_IO_VTK_OUTPUT_MODULE_HH
31#include <dune/common/timer.hh>
32#include <dune/common/fvector.hh>
33#include <dune/common/typetraits.hh>
35#include <dune/geometry/type.hh>
36#include <dune/geometry/multilineargeometry.hh>
38#include <dune/grid/common/mcmgmapper.hh>
39#include <dune/grid/common/partitionset.hh>
40#include <dune/grid/io/file/vtk/vtkwriter.hh>
41#include <dune/grid/io/file/vtk/vtksequencewriter.hh>
42#include <dune/grid/common/partitionset.hh>
59template<
class Gr
idGeometry>
62 using GridView =
typename GridGeometry::GridView;
63 static constexpr int dim = GridView::dimension;
70 const std::string&
name,
72 Dune::VTK::DataMode dm = Dune::VTK::conforming,
80 const auto precisionString = getParamFromGroup<std::string>(
paramGroup,
"Vtk.Precision",
"Float32");
83 writer_ = std::make_shared<Dune::VTKWriter<GridView>>(
gridGeometry.gridView(), dm, coordPrecision);
84 sequenceWriter_ = std::make_unique<Dune::VTKSequenceWriter<GridView>>(writer_,
name);
85 addProcessRank_ = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank",
true);
92 {
return paramGroup_; }
103 template<
typename Vector>
105 const std::string&
name,
119 template<
typename Vector>
121 const std::string&
name,
126 const auto nComp = getNumberOfComponents_(v);
128 const auto numElemDofs =
gridGeometry().elementMapper().size();
129 const auto numVertexDofs =
gridGeometry().vertexMapper().size();
134 if(numElemDofs == numVertexDofs)
135 DUNE_THROW(Dune::InvalidStateException,
"Automatic deduction of FieldType failed. Please explicitly specify FieldType::element or FieldType::vertex.");
137 if(v.size() == numElemDofs)
139 else if(v.size() == numVertexDofs)
142 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
148 if(v.size() != numElemDofs)
149 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
152 if(v.size() != numVertexDofs)
153 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
158 fields_.emplace_back(gridGeometry_.gridView(), gridGeometry_.elementMapper(), v,
name, nComp, 0, dm_,
precision);
160 fields_.emplace_back(gridGeometry_.gridView(), gridGeometry_.vertexMapper(), v,
name, nComp, dim, dm_,
precision);
169 fields_.push_back(std::move(field));
177 void write(
double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
182 if (dm_ == Dune::VTK::conforming)
183 writeConforming_(time, type);
184 else if (dm_ == Dune::VTK::nonconforming)
185 writeNonConforming_(time, type);
187 DUNE_THROW(Dune::NotImplemented,
"Output for provided VTK data mode");
192 std::cout << Fmt::format(
"Writing output for problem \"{}\". Took {:.2g} seconds.\n", name_, timer.elapsed());
199 const std::string&
name()
const {
return name_; }
200 Dune::VTK::DataMode
dataMode()
const {
return dm_; }
201 Dumux::Vtk::Precision
precision()
const {
return precision_; }
203 Dune::VTKWriter<GridView>&
writer() {
return *writer_; }
206 const std::vector<Field>&
fields()
const {
return fields_; }
210 virtual void writeConforming_(
double time, Dune::VTK::OutputType type)
217 std::vector<int> rank;
220 if (!fields_.empty() || addProcessRank_)
222 const auto numCells = gridGeometry_.gridView().size(0);
227 rank.resize(numCells);
229 for (
const auto& element : elements(gridGeometry_.gridView(), Dune::Partitions::interior))
231 const auto eIdxGlobal = gridGeometry_.elementMapper().index(element);
232 rank[eIdxGlobal] = gridGeometry_.gridView().comm().rank();
242 this->
sequenceWriter().addCellData(
Field(gridGeometry_.gridView(), gridGeometry_.elementMapper(), rank,
"process rank", 1, 0).get());
245 for (
auto&& field : fields_)
247 if (field.codim() == 0)
249 else if (field.codim() == dim)
252 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
268 virtual void writeNonConforming_(
double time, Dune::VTK::OutputType type)
270 DUNE_THROW(Dune::NotImplemented,
"Non-conforming VTK output");
274 template<
class Vector>
275 std::size_t getNumberOfComponents_(
const Vector& v)
277 if constexpr (Dune::IsIndexable<decltype(std::declval<Vector>()[0])>::value)
283 const GridGeometry& gridGeometry_;
285 const std::string paramGroup_;
286 Dune::VTK::DataMode dm_;
288 Dumux::Vtk::Precision precision_;
290 std::shared_ptr<Dune::VTKWriter<GridView>> writer_;
291 std::unique_ptr<Dune::VTKSequenceWriter<GridView>> sequenceWriter_;
293 std::vector<Field> fields_;
295 bool addProcessRank_ =
true;
310template<
class Gr
idVariables,
class SolutionVector>
314 using GridGeometry =
typename GridVariables::GridGeometry;
316 using VV =
typename GridVariables::VolumeVariables;
317 using Scalar =
typename GridVariables::Scalar;
319 using GridView =
typename GridGeometry::GridView;
322 dim = GridView::dimension,
323 dimWorld = GridView::dimensionworld
326 using Element =
typename GridView::template Codim<0>::Entity;
327 using VolVarsVector = Dune::FieldVector<Scalar, dimWorld>;
333 struct VolVarScalarDataInfo { std::function<Scalar(
const VV&)> get; std::string
name; Dumux::Vtk::Precision precision_; };
334 struct VolVarVectorDataInfo { std::function<VolVarsVector(
const VV&)> get; std::string
name; Dumux::Vtk::Precision precision_; };
345 const SolutionVector&
sol,
346 const std::string&
name,
348 Dune::VTK::DataMode dm = Dune::VTK::conforming,
355 enableVelocityOutput_ = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddVelocity",
false);
356 addProcessRank_ = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank",
true);
377 const std::string&
name)
379 volVarScalarDataInfo_.push_back(VolVarScalarDataInfo{f,
name, this->
precision()});
386 template<
class VVV = VolVarsVector,
typename std::enable_if_t<(VVV::dimension > 1),
int> = 0>
388 const std::string&
name)
390 volVarVectorDataInfo_.push_back(VolVarVectorDataInfo{f,
name, this->
precision()});
395 const auto&
problem()
const {
return gridVariables_.curGridVolVars().problem(); }
397 const GridGeometry&
gridGeometry()
const {
return gridVariables_.gridGeometry(); }
398 const SolutionVector&
sol()
const {
return sol_; }
409 void writeConforming_(
double time, Dune::VTK::OutputType type)
override
411 const Dune::VTK::DataMode dm = Dune::VTK::conforming;
418 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
421 std::vector<double> rank;
424 std::vector<std::vector<Scalar>> volVarScalarData;
425 std::vector<std::vector<VolVarsVector>> volVarVectorData;
428 if (!volVarScalarDataInfo_.empty()
429 || !volVarVectorDataInfo_.empty()
430 || !this->fields().empty()
431 || velocityOutput_->enableOutput()
434 const auto numCells =
gridGeometry().gridView().size(0);
435 const auto numDofs = numDofs_();
438 if (!volVarScalarDataInfo_.empty())
439 volVarScalarData.resize(volVarScalarDataInfo_.size(), std::vector<Scalar>(numDofs));
440 if (!volVarVectorDataInfo_.empty())
441 volVarVectorData.resize(volVarVectorDataInfo_.size(), std::vector<VolVarsVector>(numDofs));
443 if (velocityOutput_->enableOutput())
445 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
448 velocity[phaseIdx].resize(numCells);
450 velocity[phaseIdx].resize(numDofs);
453 if(isBox && dim == 1)
454 velocity[phaseIdx].resize(numCells);
456 velocity[phaseIdx].resize(numDofs);
462 if (addProcessRank_) rank.resize(numCells);
465 auto elemVolVars =
localView(gridVariables_.curGridVolVars());
466 for (
const auto& element : elements(
gridGeometry().gridView(), Dune::Partitions::interior))
468 const auto eIdxGlobal =
gridGeometry().elementMapper().index(element);
471 if (velocityOutput_->enableOutput())
473 fvGeometry.bind(element);
474 elemVolVars.bind(element, fvGeometry, sol_);
478 fvGeometry.bindElement(element);
479 elemVolVars.bindElement(element, fvGeometry, sol_);
482 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
484 for (
const auto& scv : scvs(fvGeometry))
486 const auto dofIdxGlobal = scv.dofIndex();
487 const auto& volVars = elemVolVars[scv];
490 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
491 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo_[i].get(volVars);
494 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
495 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo_[i].get(volVars);
500 if (velocityOutput_->enableOutput())
502 const auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
504 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
505 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
510 rank[eIdxGlobal] =
static_cast<double>(
gridGeometry().gridView().comm().rank());
518 if constexpr (isBox || isPQ1Bubble)
520 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
522 volVarScalarDataInfo_[i].name, 1, dim, dm, this->
precision()).get() );
523 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
525 volVarVectorDataInfo_[i].name, dimWorld, dim, dm, this->
precision()).get() );
527 if constexpr (isPQ1Bubble)
529 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
531 volVarScalarDataInfo_[i].name, 1, 0,dm, this->
precision()).get() );
532 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
534 volVarVectorDataInfo_[i].name, dimWorld, 0,dm, this->
precision()).get() );
540 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
542 volVarScalarDataInfo_[i].name, 1, 0,dm, this->
precision()).get() );
543 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
545 volVarVectorDataInfo_[i].name, dimWorld, 0,dm, this->
precision()).get() );
549 if (velocityOutput_->enableOutput())
554 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
556 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
557 dimWorld, dim, dm, this->
precision()).get() );
562 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
564 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
565 dimWorld, 0, dm, this->
precision()).get() );
574 for (
auto&& field : this->
fields())
576 if (field.codim() == 0)
578 else if (field.codim() == dim)
581 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
597 void writeNonConforming_(
double time, Dune::VTK::OutputType type)
override
599 const Dune::VTK::DataMode dm = Dune::VTK::nonconforming;
602 if(!isBox && !isDiamond)
603 DUNE_THROW(Dune::NotImplemented,
604 "Non-conforming output for discretization scheme " << GridGeometry::discMethod
612 if (enableVelocityOutput_ && !velocityOutput_->enableOutput())
613 std::cerr <<
"Warning! Velocity output was enabled in the input file"
614 <<
" but no velocity output policy was set for the VTK output module:"
615 <<
" There will be no velocity output."
616 <<
" Use the addVelocityOutput member function of the VTK output module." << std::endl;
618 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
621 std::vector<double> rank;
624 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
625 using VectorDataContainer = std::vector< std::vector<VolVarsVector> >;
626 std::vector< ScalarDataContainer > volVarScalarData;
627 std::vector< VectorDataContainer > volVarVectorData;
630 if (!volVarScalarDataInfo_.empty()
631 || !volVarVectorDataInfo_.empty()
632 || !this->fields().empty()
633 || velocityOutput_->enableOutput()
636 const auto numCells =
gridGeometry().gridView().size(0);
637 const auto outputSize = numDofs_();
640 if (!volVarScalarDataInfo_.empty())
641 volVarScalarData.resize(volVarScalarDataInfo_.size(), ScalarDataContainer(numCells));
642 if (!volVarVectorDataInfo_.empty())
643 volVarVectorData.resize(volVarVectorDataInfo_.size(), VectorDataContainer(numCells));
645 if (velocityOutput_->enableOutput())
647 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
649 if((isBox && dim == 1) || isDiamond)
650 velocity[phaseIdx].resize(numCells);
652 velocity[phaseIdx].resize(outputSize);
657 if (addProcessRank_) rank.resize(numCells);
661 auto elemVolVars =
localView(gridVariables_.curGridVolVars());
662 for (
const auto& element : elements(
gridGeometry().gridView(), Dune::Partitions::interior))
664 const auto eIdxGlobal =
gridGeometry().elementMapper().index(element);
667 if (velocityOutput_->enableOutput())
669 fvGeometry.bind(element);
670 elemVolVars.bind(element, fvGeometry, sol_);
674 fvGeometry.bindElement(element);
675 elemVolVars.bindElement(element, fvGeometry, sol_);
678 const auto numLocalDofs = fvGeometry.numScv();
680 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
681 volVarScalarData[i][eIdxGlobal].resize(numLocalDofs);
682 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
683 volVarVectorData[i][eIdxGlobal].resize(numLocalDofs);
685 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
687 for (
const auto& scv : scvs(fvGeometry))
689 const auto& volVars = elemVolVars[scv];
692 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
693 volVarScalarData[i][eIdxGlobal][scv.localDofIndex()] = volVarScalarDataInfo_[i].get(volVars);
696 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
697 volVarVectorData[i][eIdxGlobal][scv.localDofIndex()] = volVarVectorDataInfo_[i].get(volVars);
702 if (velocityOutput_->enableOutput())
704 const auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
705 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
706 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
711 rank[eIdxGlobal] =
static_cast<double>(
gridGeometry().gridView().comm().rank());
719 static constexpr int dofLocCodim = isDiamond ? 1 : dim;
720 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
723 volVarScalarData[i], volVarScalarDataInfo_[i].
name,
727 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
730 volVarVectorData[i], volVarVectorDataInfo_[i].
name,
731 dimWorld, dofLocCodim, dm, this->
precision()
735 if (velocityOutput_->enableOutput())
738 if (dim > 1 && !isDiamond)
739 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
742 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
743 dimWorld, dofLocCodim, dm, this->precision()
748 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
751 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
752 dimWorld, 0, dm, this->precision()
761 rank,
"process rank", 1, 0
765 for (
const auto& field : this->
fields())
767 if (field.codim() == 0)
769 else if (field.codim() == dim || field.codim() == 1)
772 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
787 std::size_t numDofs_()
const
791 if constexpr (isBox || isDiamond || isPQ1Bubble)
797 const GridVariables& gridVariables_;
798 const SolutionVector& sol_;
800 std::vector<VolVarScalarDataInfo> volVarScalarDataInfo_;
801 std::vector<VolVarVectorDataInfo> volVarVectorDataInfo_;
803 std::shared_ptr<VelocityOutput> velocityOutput_;
804 bool enableVelocityOutput_ =
false;
805 bool addProcessRank_ =
true;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Vtk field types available in Dumux.
Dune style VTK functions.
Formatting based on the fmt-library which implements std::format of C++20.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Precision stringToPrecision(std::string_view precisionName)
Maps a string (e.g. from input) to a Dune precision type.
Definition: precision.hh:39
FieldType
Identifier for vtk field types.
Definition: fieldtype.hh:34
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
constexpr FCDiamond fcdiamond
Definition: method.hh:141
constexpr Box box
Definition: method.hh:136
constexpr PQ1Bubble pq1bubble
Definition: method.hh:137
Velocity output for implicit (porous media) models.
Definition: io/velocityoutput.hh:41
std::vector< Dune::FieldVector< Scalar, dimWorld > > VelocityVector
Definition: io/velocityoutput.hh:50
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:61
Dumux::Vtk::Precision precision() const
Definition: io/vtkoutputmodule.hh:201
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:91
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Definition: io/vtkoutputmodule.hh:177
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:69
Dune::VTK::DataMode dataMode() const
Definition: io/vtkoutputmodule.hh:200
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:203
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:204
const std::string & name() const
Definition: io/vtkoutputmodule.hh:199
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:206
virtual ~VtkOutputModuleBase()=default
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:120
bool verbose() const
Definition: io/vtkoutputmodule.hh:198
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:67
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:104
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:196
void addField(Field &&field)
Add a scalar or vector valued vtk field.
Definition: io/vtkoutputmodule.hh:167
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:312
const auto & problem() const
Definition: io/vtkoutputmodule.hh:395
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:396
void addVolumeVariable(std::function< VolVarsVector(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:387
VV VolumeVariables
export type of the volume variables for the outputfields
Definition: io/vtkoutputmodule.hh:342
void addVolumeVariable(std::function< Scalar(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:376
void addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:370
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:404
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:340
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:398
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:397
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:344
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:400
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:401