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!");
146 fields_.emplace_back(gridGeometry_.gridView(), gridGeometry_.elementMapper(), v,
name, nComp, 0, dm_,
precision);
148 fields_.emplace_back(gridGeometry_.gridView(), gridGeometry_.vertexMapper(), v,
name, nComp, dim, dm_,
precision);
157 fields_.push_back(std::move(field));
165 void write(
double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
170 if (dm_ == Dune::VTK::conforming)
171 writeConforming_(time, type);
172 else if (dm_ == Dune::VTK::nonconforming)
173 writeNonConforming_(time, type);
175 DUNE_THROW(Dune::NotImplemented,
"Output for provided VTK data mode");
180 std::cout << Fmt::format(
"Writing output for problem \"{}\". Took {:.2g} seconds.\n", name_, timer.elapsed());
187 const std::string&
name()
const {
return name_; }
188 Dune::VTK::DataMode
dataMode()
const {
return dm_; }
189 Dumux::Vtk::Precision
precision()
const {
return precision_; }
191 Dune::VTKWriter<GridView>&
writer() {
return *writer_; }
194 const std::vector<Field>&
fields()
const {
return fields_; }
198 virtual void writeConforming_(
double time, Dune::VTK::OutputType type)
205 std::vector<int> rank;
208 if (!fields_.empty() || addProcessRank_)
210 const auto numCells = gridGeometry_.gridView().size(0);
215 rank.resize(numCells);
217 for (
const auto& element : elements(gridGeometry_.gridView(), Dune::Partitions::interior))
219 const auto eIdxGlobal = gridGeometry_.elementMapper().index(element);
220 rank[eIdxGlobal] = gridGeometry_.gridView().comm().rank();
230 this->
sequenceWriter().addCellData(
Field(gridGeometry_.gridView(), gridGeometry_.elementMapper(), rank,
"process rank", 1, 0).get());
233 for (
auto&& field : fields_)
235 if (field.codim() == 0)
237 else if (field.codim() == dim)
240 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
256 virtual void writeNonConforming_(
double time, Dune::VTK::OutputType type)
258 DUNE_THROW(Dune::NotImplemented,
"Non-conforming VTK output");
262 template<
class Vector>
263 std::size_t getNumberOfComponents_(
const Vector& v)
265 if constexpr (Dune::IsIndexable<decltype(std::declval<Vector>()[0])>::value)
271 const GridGeometry& gridGeometry_;
273 const std::string paramGroup_;
274 Dune::VTK::DataMode dm_;
276 Dumux::Vtk::Precision precision_;
278 std::shared_ptr<Dune::VTKWriter<GridView>> writer_;
279 std::unique_ptr<Dune::VTKSequenceWriter<GridView>> sequenceWriter_;
281 std::vector<Field> fields_;
283 bool addProcessRank_ =
true;
298template<
class Gr
idVariables,
class SolutionVector>
302 using GridGeometry =
typename GridVariables::GridGeometry;
304 using VV =
typename GridVariables::VolumeVariables;
305 using Scalar =
typename GridVariables::Scalar;
307 using GridView =
typename GridGeometry::GridView;
310 dim = GridView::dimension,
311 dimWorld = GridView::dimensionworld
314 using Element =
typename GridView::template Codim<0>::Entity;
315 using VolVarsVector = Dune::FieldVector<Scalar, dimWorld>;
321 struct VolVarScalarDataInfo { std::function<Scalar(
const VV&)> get; std::string
name; Dumux::Vtk::Precision precision_; };
322 struct VolVarVectorDataInfo { std::function<VolVarsVector(
const VV&)> get; std::string
name; Dumux::Vtk::Precision precision_; };
333 const SolutionVector&
sol,
334 const std::string&
name,
336 Dune::VTK::DataMode dm = Dune::VTK::conforming,
343 enableVelocityOutput_ = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddVelocity",
false);
344 addProcessRank_ = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank",
true);
365 const std::string&
name)
367 volVarScalarDataInfo_.push_back(VolVarScalarDataInfo{f,
name, this->
precision()});
374 template<
class VVV = VolVarsVector,
typename std::enable_if_t<(VVV::dimension > 1),
int> = 0>
376 const std::string&
name)
378 volVarVectorDataInfo_.push_back(VolVarVectorDataInfo{f,
name, this->
precision()});
383 const auto&
problem()
const {
return gridVariables_.curGridVolVars().problem(); }
385 const GridGeometry&
gridGeometry()
const {
return gridVariables_.gridGeometry(); }
386 const SolutionVector&
sol()
const {
return sol_; }
397 void writeConforming_(
double time, Dune::VTK::OutputType type)
override
399 const Dune::VTK::DataMode dm = Dune::VTK::conforming;
406 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
409 std::vector<double> rank;
412 std::vector<std::vector<Scalar>> volVarScalarData;
413 std::vector<std::vector<VolVarsVector>> volVarVectorData;
416 if (!volVarScalarDataInfo_.empty()
417 || !volVarVectorDataInfo_.empty()
418 || !this->fields().empty()
419 || velocityOutput_->enableOutput()
422 const auto numCells =
gridGeometry().gridView().size(0);
423 const auto numDofs = numDofs_();
426 if (!volVarScalarDataInfo_.empty())
427 volVarScalarData.resize(volVarScalarDataInfo_.size(), std::vector<Scalar>(numDofs));
428 if (!volVarVectorDataInfo_.empty())
429 volVarVectorData.resize(volVarVectorDataInfo_.size(), std::vector<VolVarsVector>(numDofs));
431 if (velocityOutput_->enableOutput())
433 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
436 velocity[phaseIdx].resize(numCells);
438 velocity[phaseIdx].resize(numDofs);
441 if(isBox && dim == 1)
442 velocity[phaseIdx].resize(numCells);
444 velocity[phaseIdx].resize(numDofs);
450 if (addProcessRank_) rank.resize(numCells);
453 auto elemVolVars =
localView(gridVariables_.curGridVolVars());
454 for (
const auto& element : elements(
gridGeometry().gridView(), Dune::Partitions::interior))
456 const auto eIdxGlobal =
gridGeometry().elementMapper().index(element);
459 if (velocityOutput_->enableOutput())
461 fvGeometry.bind(element);
462 elemVolVars.bind(element, fvGeometry, sol_);
466 fvGeometry.bindElement(element);
467 elemVolVars.bindElement(element, fvGeometry, sol_);
470 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
472 for (
const auto& scv : scvs(fvGeometry))
474 const auto dofIdxGlobal = scv.dofIndex();
475 const auto& volVars = elemVolVars[scv];
478 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
479 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo_[i].get(volVars);
482 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
483 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo_[i].get(volVars);
488 if (velocityOutput_->enableOutput())
490 const auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
492 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
493 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
498 rank[eIdxGlobal] =
static_cast<double>(
gridGeometry().gridView().comm().rank());
506 if constexpr (isBox || isPQ1Bubble)
508 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
510 volVarScalarDataInfo_[i].name, 1, dim, dm, this->
precision()).get() );
511 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
513 volVarVectorDataInfo_[i].name, dimWorld, dim, dm, this->
precision()).get() );
515 if constexpr (isPQ1Bubble)
517 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
519 volVarScalarDataInfo_[i].name, 1, 0,dm, this->
precision()).get() );
520 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
522 volVarVectorDataInfo_[i].name, dimWorld, 0,dm, this->
precision()).get() );
528 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
530 volVarScalarDataInfo_[i].name, 1, 0,dm, this->
precision()).get() );
531 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
533 volVarVectorDataInfo_[i].name, dimWorld, 0,dm, this->
precision()).get() );
537 if (velocityOutput_->enableOutput())
542 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
544 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
545 dimWorld, dim, dm, this->
precision()).get() );
550 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
552 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
553 dimWorld, 0, dm, this->
precision()).get() );
562 for (
auto&& field : this->
fields())
564 if (field.codim() == 0)
566 else if (field.codim() == dim)
569 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
585 void writeNonConforming_(
double time, Dune::VTK::OutputType type)
override
587 const Dune::VTK::DataMode dm = Dune::VTK::nonconforming;
590 if(!isBox && !isDiamond)
591 DUNE_THROW(Dune::NotImplemented,
592 "Non-conforming output for discretization scheme " << GridGeometry::discMethod
600 if (enableVelocityOutput_ && !velocityOutput_->enableOutput())
601 std::cerr <<
"Warning! Velocity output was enabled in the input file"
602 <<
" but no velocity output policy was set for the VTK output module:"
603 <<
" There will be no velocity output."
604 <<
" Use the addVelocityOutput member function of the VTK output module." << std::endl;
606 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
609 std::vector<double> rank;
612 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
613 using VectorDataContainer = std::vector< std::vector<VolVarsVector> >;
614 std::vector< ScalarDataContainer > volVarScalarData;
615 std::vector< VectorDataContainer > volVarVectorData;
618 if (!volVarScalarDataInfo_.empty()
619 || !volVarVectorDataInfo_.empty()
620 || !this->fields().empty()
621 || velocityOutput_->enableOutput()
624 const auto numCells =
gridGeometry().gridView().size(0);
625 const auto outputSize = numDofs_();
628 if (!volVarScalarDataInfo_.empty())
629 volVarScalarData.resize(volVarScalarDataInfo_.size(), ScalarDataContainer(numCells));
630 if (!volVarVectorDataInfo_.empty())
631 volVarVectorData.resize(volVarVectorDataInfo_.size(), VectorDataContainer(numCells));
633 if (velocityOutput_->enableOutput())
635 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
637 if((isBox && dim == 1) || isDiamond)
638 velocity[phaseIdx].resize(numCells);
640 velocity[phaseIdx].resize(outputSize);
645 if (addProcessRank_) rank.resize(numCells);
649 auto elemVolVars =
localView(gridVariables_.curGridVolVars());
650 for (
const auto& element : elements(
gridGeometry().gridView(), Dune::Partitions::interior))
652 const auto eIdxGlobal =
gridGeometry().elementMapper().index(element);
655 if (velocityOutput_->enableOutput())
657 fvGeometry.bind(element);
658 elemVolVars.bind(element, fvGeometry, sol_);
662 fvGeometry.bindElement(element);
663 elemVolVars.bindElement(element, fvGeometry, sol_);
666 const auto numLocalDofs = fvGeometry.numScv();
668 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
669 volVarScalarData[i][eIdxGlobal].resize(numLocalDofs);
670 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
671 volVarVectorData[i][eIdxGlobal].resize(numLocalDofs);
673 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
675 for (
const auto& scv : scvs(fvGeometry))
677 const auto& volVars = elemVolVars[scv];
680 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
681 volVarScalarData[i][eIdxGlobal][scv.localDofIndex()] = volVarScalarDataInfo_[i].get(volVars);
684 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
685 volVarVectorData[i][eIdxGlobal][scv.localDofIndex()] = volVarVectorDataInfo_[i].get(volVars);
690 if (velocityOutput_->enableOutput())
692 const auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
693 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
694 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
699 rank[eIdxGlobal] =
static_cast<double>(
gridGeometry().gridView().comm().rank());
707 static constexpr int dofLocCodim = isDiamond ? 1 : dim;
708 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
711 volVarScalarData[i], volVarScalarDataInfo_[i].
name,
715 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
718 volVarVectorData[i], volVarVectorDataInfo_[i].
name,
719 dimWorld, dofLocCodim, dm, this->
precision()
723 if (velocityOutput_->enableOutput())
726 if (dim > 1 && !isDiamond)
727 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
730 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
731 dimWorld, dofLocCodim, dm, this->precision()
736 for (
int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
739 "velocity_" + velocityOutput_->phaseName(phaseIdx) +
" (m/s)",
740 dimWorld, 0, dm, this->precision()
749 rank,
"process rank", 1, 0
753 for (
const auto& field : this->
fields())
755 if (field.codim() == 0)
757 else if (field.codim() == dim || field.codim() == 1)
760 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
775 std::size_t numDofs_()
const
779 if constexpr (isBox || isDiamond || isPQ1Bubble)
785 const GridVariables& gridVariables_;
786 const SolutionVector& sol_;
788 std::vector<VolVarScalarDataInfo> volVarScalarDataInfo_;
789 std::vector<VolVarVectorDataInfo> volVarVectorDataInfo_;
791 std::shared_ptr<VelocityOutput> velocityOutput_;
792 bool enableVelocityOutput_ =
false;
793 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:189
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:79
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Definition: io/vtkoutputmodule.hh:165
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:188
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:191
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:192
const std::string & name() const
Definition: io/vtkoutputmodule.hh:187
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:194
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:108
bool verbose() const
Definition: io/vtkoutputmodule.hh:186
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:55
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
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:184
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:300
const auto & problem() const
Definition: io/vtkoutputmodule.hh:383
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:384
void addVolumeVariable(std::function< VolVarsVector(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:375
VV VolumeVariables
export type of the volume variables for the outputfields
Definition: io/vtkoutputmodule.hh:330
void addVolumeVariable(std::function< Scalar(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:364
void addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:358
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:392
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:328
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:386
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:385
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:332
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:388
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:389
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:152
constexpr Box box
Definition: method.hh:147
constexpr PQ1Bubble pq1bubble
Definition: method.hh:148
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.