version 3.10-dev
pnmvtkoutputmodule.hh
Go to the documentation of this file.
1// // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_PNM_VTK_OUTPUT_MODULE_HH
14#define DUMUX_PNM_VTK_OUTPUT_MODULE_HH
15
17#include "velocityoutput.hh"
18
19namespace Dumux::PoreNetwork {
20
25template<class GridVariables, class FluxVariables, class SolutionVector>
26class VtkOutputModule : public Dumux::VtkOutputModule<GridVariables, SolutionVector>
27{
29 using Scalar = typename GridVariables::Scalar;
30 using FluxVarsCache = typename GridVariables::GridFluxVariablesCache::FluxVariablesCache;
31
32 struct ThroatFluxDataInfo { std::function<Scalar(const FluxVariables&, const FluxVarsCache&)> get; std::string name; };
33
34public:
35
37 VtkOutputModule(const GridVariables& gridVariables,
38 const SolutionVector& sol,
39 const std::string& name,
40 const std::string& paramGroup = "",
41 Dune::VTK::DataMode dm = Dune::VTK::conforming,
42 bool verbose = true)
44 {
45 if constexpr (GridVariables::VolumeVariables::numFluidPhases() >= 1)
46 {
47 // enable velocity output per default
49 this->addVelocityOutput(std::make_shared<VelocityOutput>(gridVariables));
50 }
51 }
52
56 void addFluxVariable(std::function<Scalar(const FluxVariables&, const FluxVarsCache&)>&& f, const std::string& name)
57 {
58 throatFluxDataInfo_.push_back(ThroatFluxDataInfo{f, name});
59 const auto numElems = problem().gridGeometry().gridView().size(0);
60 throatFluxData_.push_back(std::vector<Scalar>(numElems));
61 ParentType::addField(throatFluxData_.back(), throatFluxDataInfo_.back().name, Vtk::FieldType::element);
62 }
63
65 void write(double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
66 {
67 const auto gridView = problem().gridGeometry().gridView();
68 const auto numElems = gridView.size(0);
69
70 // resize all fields
71 for (auto& data : throatFluxData_)
72 data.resize(numElems);
73
74 auto fvElementGeometry = localView(problem().gridGeometry());
75 auto elemVolVars = localView(this->gridVariables().curGridVolVars());
76 auto elemFluxVarsCache = localView(this->gridVariables().gridFluxVarsCache());
77 // iterate over all elements
78 for (const auto& element : elements(gridView, Dune::Partitions::interior))
79 {
80 const auto eIdx = problem().gridGeometry().elementMapper().index(element);
81
82 // make sure FVElementGeometry & vol vars are bound to the element
83 fvElementGeometry.bindElement(element);
84 elemVolVars.bind(element, fvElementGeometry, this->sol());
85 elemFluxVarsCache.bind(element, fvElementGeometry, elemVolVars);
86
87 // treat the throat flux related data
88 std::size_t dataIdx = 0;
89 for (auto&& scvf : scvfs(fvElementGeometry))
90 {
91 if (!scvf.boundary())
92 {
93 FluxVariables fluxVars;
94 fluxVars.init(problem(), element, fvElementGeometry, elemVolVars, scvf, elemFluxVarsCache);
95
96 dataIdx = 0;
97 for(auto& data : throatFluxData_)
98 data[eIdx] = throatFluxDataInfo_[dataIdx++].get(fluxVars, elemFluxVarsCache[scvf]);
99 }
100 }
101 }
102
103 // call the ParentType's write method to write out all data
104 ParentType::write(time, type);
105
106 // empty the data containers in order to save some memory
107 auto clearAndShrink = [] (auto& data)
108 {
109 data.clear();
110 data.shrink_to_fit();
111 };
112
113 for (auto& data : throatFluxData_)
114 clearAndShrink(data);
115 }
116
118 const auto& problem() const { return ParentType::problem(); }
119
120private:
121 std::vector<ThroatFluxDataInfo> throatFluxDataInfo_;
122 std::list<std::vector<Scalar>> throatFluxData_;
123};
124
125} //namespace Dumux::PoreNetwork
126
127#endif
Adds vtk output fields specific to pore-network models.
Definition: pnmvtkoutputmodule.hh:27
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:65
VtkOutputModule(const GridVariables &gridVariables, const SolutionVector &sol, const std::string &name, const std::string &paramGroup="", Dune::VTK::DataMode dm=Dune::VTK::conforming, bool verbose=true)
The constructor.
Definition: pnmvtkoutputmodule.hh:37
void addFluxVariable(std::function< Scalar(const FluxVariables &, const FluxVarsCache &)> &&f, const std::string &name)
Definition: pnmvtkoutputmodule.hh:56
const auto & problem() const
Return a reference to the problem.
Definition: pnmvtkoutputmodule.hh:118
Velocity output for implicit (porous media) models.
Definition: io/velocityoutput.hh:29
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:49
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
const std::string & name() const
Definition: io/vtkoutputmodule.hh:187
bool verbose() const
Definition: io/vtkoutputmodule.hh:186
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
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 addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:358
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:386
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:385
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: discretization/porenetwork/fvelementgeometry.hh:24
TODO: docme!