3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
io/vtkoutputmodule.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef VTK_OUTPUT_MODULE_HH
25#define VTK_OUTPUT_MODULE_HH
26
27#include <functional>
28
29#include <dune/common/timer.hh>
30#include <dune/common/fvector.hh>
31#include <dune/common/typetraits.hh>
32
33#include <dune/geometry/type.hh>
34#include <dune/geometry/multilineargeometry.hh>
35
36#include <dune/grid/common/mcmgmapper.hh>
37#include <dune/grid/common/partitionset.hh>
38#include <dune/grid/io/file/vtk/vtkwriter.hh>
39#include <dune/grid/io/file/vtk/vtksequencewriter.hh>
40#include <dune/grid/common/partitionset.hh>
41
45
46#include "vtkfunction.hh"
47#include "velocityoutput.hh"
48
49namespace Dumux {
50
64template<class GridVariables, class SolutionVector>
66{
67 using GridGeometry = typename GridVariables::GridGeometry;
68
69 using VV = typename GridVariables::VolumeVariables;
70 using Scalar = typename GridVariables::Scalar;
71
72 using GridView = typename GridGeometry::GridView;
73
74 enum {
75 dim = GridView::dimension,
76 dimWorld = GridView::dimensionworld
77 };
78
79 using Element = typename GridView::template Codim<0>::Entity;
80 using VolVarsVector = Dune::FieldVector<Scalar, dimWorld>;
81
82 static constexpr bool isBox = GridGeometry::discMethod == DiscretizationMethod::box;
83 static constexpr int dofCodim = isBox ? dim : 0;
84
85 struct VolVarScalarDataInfo { std::function<Scalar(const VV&)> get; std::string name; };
86 struct VolVarVectorDataInfo { std::function<VolVarsVector(const VV&)> get; std::string name; };
87 using Field = Vtk::template Field<GridView>;
88
90
91public:
93 using VolumeVariables = VV;
94
96 enum class FieldType : unsigned int
97 {
98 element, vertex, automatic
99 };
100
101 VtkOutputModule(const GridVariables& gridVariables,
102 const SolutionVector& sol,
103 const std::string& name,
104 const std::string& paramGroup = "",
105 Dune::VTK::DataMode dm = Dune::VTK::conforming,
106 bool verbose = true)
107 : gridVariables_(gridVariables)
108 , sol_(sol)
109 , name_(name)
110 , paramGroup_(paramGroup)
111 , verbose_(gridVariables.gridGeometry().gridView().comm().rank() == 0 && verbose)
112 , dm_(dm)
113 , writer_(std::make_shared<Dune::VTKWriter<GridView>>(gridVariables.gridGeometry().gridView(), dm))
114 , sequenceWriter_(writer_, name)
115 , velocityOutput_(std::make_shared<VelocityOutputType>())
116 {}
117
119 const std::string& paramGroup() const
120 { return paramGroup_; }
121
126
133 void addVelocityOutput(std::shared_ptr<VelocityOutputType> velocityOutput)
134 { velocityOutput_ = velocityOutput; }
135
139 void addVolumeVariable(std::function<Scalar(const VolumeVariables&)>&& f, const std::string& name)
140 {
141 volVarScalarDataInfo_.push_back(VolVarScalarDataInfo{f, name});
142 }
143
148 template<class VVV = VolVarsVector, typename std::enable_if_t<(VVV::dimension > 1), int> = 0>
149 void addVolumeVariable(std::function<VolVarsVector(const VolumeVariables&)>&& f, const std::string& name)
150 {
151 volVarVectorDataInfo_.push_back(VolVarVectorDataInfo{f, name});
152 }
153
160 template<typename Vector>
161 void addField(const Vector& v, const std::string& name, FieldType fieldType = FieldType::automatic)
162 {
163 // Deduce the number of components from the given vector type
164 const auto nComp = getNumberOfComponents_(v);
165
166 const auto numElemDofs = gridGeometry().elementMapper().size();
167 const auto numVertexDofs = gridGeometry().vertexMapper().size();
168
169 // Automatically deduce the field type ...
170 if(fieldType == FieldType::automatic)
171 {
172 if(numElemDofs == numVertexDofs)
173 DUNE_THROW(Dune::InvalidStateException, "Automatic deduction of FieldType failed. Please explicitly specify FieldType::element or FieldType::vertex.");
174
175 if(v.size() == numElemDofs)
176 fieldType = FieldType::element;
177 else if(v.size() == numVertexDofs)
178 fieldType = FieldType::vertex;
179 else
180 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
181 }
182 // ... or check if the user-specified type matches the size of v
183 else
184 {
185 if(fieldType == FieldType::element)
186 if(v.size() != numElemDofs)
187 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
188
189 if(fieldType == FieldType::vertex)
190 if(v.size() != numVertexDofs)
191 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
192 }
193
194 // add the appropriate field
195 if (fieldType == FieldType::element)
196 fields_.emplace_back(gridGeometry().gridView(), gridGeometry().elementMapper(), v, name, nComp, 0);
197 else
198 fields_.emplace_back(gridGeometry().gridView(), gridGeometry().vertexMapper(), v, name, nComp, dim);
199 }
200
204
210 void write(double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
211 {
212 Dune::Timer timer;
213
214 // write to file depending on data mode
215 if (dm_ == Dune::VTK::conforming)
216 writeConforming_(time, type);
217 else if (dm_ == Dune::VTK::nonconforming)
218 writeNonConforming_(time, type);
219 else
220 DUNE_THROW(Dune::NotImplemented, "Output for provided vtk data mode");
221
223 timer.stop();
224 if (verbose_)
225 {
226 std::cout << "Writing output for problem \"" << name_ << "\". Took " << timer.elapsed() << " seconds." << std::endl;
227 }
228 }
229
230protected:
231 // some return functions for differing implementations to use
232 const auto& problem() const { return gridVariables_.curGridVolVars().problem(); }
233 const GridGeometry& gridGeometry() const { return gridVariables_.gridGeometry(); }
234 const GridVariables& gridVariables() const { return gridVariables_; }
235 const SolutionVector& sol() const { return sol_; }
236
237 bool verbose() const { return verbose_; }
238 const std::string& name() const { return name_; }
239 Dune::VTK::DataMode dataMode() const { return dm_; }
240
241 Dune::VTKWriter<GridView>& writer() { return *writer_; }
242 Dune::VTKSequenceWriter<GridView>& sequenceWriter() { return sequenceWriter_; }
243
244 const std::vector<VolVarScalarDataInfo>& volVarScalarDataInfo() const { return volVarScalarDataInfo_; }
245 const std::vector<VolVarVectorDataInfo>& volVarVectorDataInfo() const { return volVarVectorDataInfo_; }
246 const std::vector<Field>& fields() const { return fields_; }
247
249 const VelocityOutput& velocityOutput() const { return *velocityOutput_; }
250
251private:
252
254 void writeConforming_(double time, Dune::VTK::OutputType type)
255 {
259
260 // instantiate the velocity output
261 using VelocityVector = typename VelocityOutput::VelocityVector;
262 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
263
264 // process rank
265 static bool addProcessRank = getParamFromGroup<bool>(paramGroup_, "Vtk.AddProcessRank");
266 std::vector<double> rank;
267
268 // volume variable data
269 std::vector<std::vector<Scalar>> volVarScalarData;
270 std::vector<std::vector<VolVarsVector>> volVarVectorData;
271
273 if (!volVarScalarDataInfo_.empty()
274 || !volVarVectorDataInfo_.empty()
275 || !fields_.empty()
276 || velocityOutput_->enableOutput()
277 || addProcessRank)
278 {
279 const auto numCells = gridGeometry().gridView().size(0);
280 const auto numDofs = numDofs_();
281
282 // get fields for all volume variables
283 if (!volVarScalarDataInfo_.empty())
284 volVarScalarData.resize(volVarScalarDataInfo_.size(), std::vector<Scalar>(numDofs));
285 if (!volVarVectorDataInfo_.empty())
286 volVarVectorData.resize(volVarVectorDataInfo_.size(), std::vector<VolVarsVector>(numDofs));
287
288 if (velocityOutput_->enableOutput())
289 {
290 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
291 {
292 if(isBox && dim == 1)
293 velocity[phaseIdx].resize(numCells);
294 else
295 velocity[phaseIdx].resize(numDofs);
296 }
297 }
298
299 // maybe allocate space for the process rank
300 if (addProcessRank) rank.resize(numCells);
301
302 for (const auto& element : elements(gridGeometry().gridView(), Dune::Partitions::interior))
303 {
304 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
305
306 auto fvGeometry = localView(gridGeometry());
307 auto elemVolVars = localView(gridVariables_.curGridVolVars());
308
309 // If velocity output is enabled we need to bind to the whole stencil
310 // otherwise element-local data is sufficient
311 if (velocityOutput_->enableOutput())
312 {
313 fvGeometry.bind(element);
314 elemVolVars.bind(element, fvGeometry, sol_);
315 }
316 else
317 {
318 fvGeometry.bindElement(element);
319 elemVolVars.bindElement(element, fvGeometry, sol_);
320 }
321
322 if (!volVarScalarDataInfo_.empty()
323 || !volVarVectorDataInfo_.empty())
324 {
325 for (auto&& scv : scvs(fvGeometry))
326 {
327 const auto dofIdxGlobal = scv.dofIndex();
328 const auto& volVars = elemVolVars[scv];
329
330 // get the scalar-valued data
331 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
332 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo_[i].get(volVars);
333
334 // get the vector-valued data
335 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
336 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo_[i].get(volVars);
337 }
338 }
339
340 // velocity output
341 if (velocityOutput_->enableOutput())
342 {
343 auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache());
344 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
345
346 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
347 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
348 }
349
351 if (addProcessRank)
352 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
353 }
354
358
359 // volume variables if any
360 if (isBox)
361 {
362 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
363 sequenceWriter_.addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), volVarScalarData[i],
364 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/dim).get() );
365 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
366 sequenceWriter_.addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), volVarVectorData[i],
367 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/dim).get() );
368 }
369 else
370 {
371 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
372 sequenceWriter_.addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarScalarData[i],
373 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0).get() );
374 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
375 sequenceWriter_.addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarVectorData[i],
376 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0).get() );
377 }
378
379 // the velocity field
380 if (velocityOutput_->enableOutput())
381 {
382 if (isBox && dim > 1)
383 {
384 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
385 sequenceWriter_.addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
386 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
387 /*numComp*/dimWorld, /*codim*/dim).get() );
388 }
389 // cell-centered models
390 else
391 {
392 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
393 sequenceWriter_.addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
394 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
395 /*numComp*/dimWorld, /*codim*/0).get() );
396 }
397 }
398
399 // the process rank
400 if (addProcessRank)
401 sequenceWriter_.addCellData(Field(gridGeometry().gridView(), gridGeometry().elementMapper(), rank, "process rank", 1, 0).get());
402
403 // also register additional (non-standardized) user fields if any
404 for (auto&& field : fields_)
405 {
406 if (field.codim() == 0)
407 sequenceWriter_.addCellData(field.get());
408 else if (field.codim() == dim)
409 sequenceWriter_.addVertexData(field.get());
410 else
411 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
412 }
413 }
414
418 sequenceWriter_.write(time, type);
419
423 writer_->clear();
424 }
425
427 void writeNonConforming_(double time, Dune::VTK::OutputType type)
428 {
429 if(!isBox)
430 DUNE_THROW(Dune::InvalidStateException, "Non-conforming output makes no sense for cell-centered schemes!");
431
435
436 // check the velocity output
437 bool enableVelocityOutput = getParamFromGroup<bool>(paramGroup_, "Vtk.AddVelocity");
438 if (enableVelocityOutput == true && !velocityOutput_->enableOutput())
439 std::cerr << "Warning! Velocity output was enabled in the input file"
440 << " but no velocity output policy was set for the VTK output module:"
441 << " There will be no velocity output."
442 << " Use the addVelocityOutput member function of the VTK output module." << std::endl;
443 using VelocityVector = typename VelocityOutput::VelocityVector;
444 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
445
446 // process rank
447 static bool addProcessRank = getParamFromGroup<bool>(paramGroup_, "Vtk.AddProcessRank");
448 std::vector<double> rank;
449
450 // volume variable data (indexing: volvardata/element/localcorner)
451 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
452 using VectorDataContainer = std::vector< std::vector<VolVarsVector> >;
453 std::vector< ScalarDataContainer > volVarScalarData;
454 std::vector< VectorDataContainer > volVarVectorData;
455
457 if (!volVarScalarDataInfo_.empty()
458 || !volVarVectorDataInfo_.empty()
459 || !fields_.empty()
460 || velocityOutput_->enableOutput()
461 || addProcessRank)
462 {
463 const auto numCells = gridGeometry().gridView().size(0);
464 const auto numDofs = numDofs_();
465
466 // get fields for all volume variables
467 if (!volVarScalarDataInfo_.empty())
468 volVarScalarData.resize(volVarScalarDataInfo_.size(), ScalarDataContainer(numCells));
469 if (!volVarVectorDataInfo_.empty())
470 volVarVectorData.resize(volVarVectorDataInfo_.size(), VectorDataContainer(numCells));
471
472 if (velocityOutput_->enableOutput())
473 {
474 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
475 {
476 if(isBox && dim == 1)
477 velocity[phaseIdx].resize(numCells);
478 else
479 velocity[phaseIdx].resize(numDofs);
480 }
481 }
482
483 // maybe allocate space for the process rank
484 if (addProcessRank) rank.resize(numCells);
485
486 for (const auto& element : elements(gridGeometry().gridView(), Dune::Partitions::interior))
487 {
488 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
489 const auto numCorners = element.subEntities(dim);
490
491 auto fvGeometry = localView(gridGeometry());
492 auto elemVolVars = localView(gridVariables_.curGridVolVars());
493
494 // resize element-local data containers
495 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
496 volVarScalarData[i][eIdxGlobal].resize(numCorners);
497 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
498 volVarVectorData[i][eIdxGlobal].resize(numCorners);
499
500 // If velocity output is enabled we need to bind to the whole stencil
501 // otherwise element-local data is sufficient
502 if (velocityOutput_->enableOutput())
503 {
504 fvGeometry.bind(element);
505 elemVolVars.bind(element, fvGeometry, sol_);
506 }
507 else
508 {
509 fvGeometry.bindElement(element);
510 elemVolVars.bindElement(element, fvGeometry, sol_);
511 }
512
513 if (!volVarScalarDataInfo_.empty()
514 || !volVarVectorDataInfo_.empty())
515 {
516 for (auto&& scv : scvs(fvGeometry))
517 {
518 const auto& volVars = elemVolVars[scv];
519
520 // get the scalar-valued data
521 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
522 volVarScalarData[i][scv.elementIndex()][scv.localDofIndex()] = volVarScalarDataInfo_[i].get(volVars);
523
524 // get the vector-valued data
525 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
526 volVarVectorData[i][scv.elementIndex()][scv.localDofIndex()] = volVarVectorDataInfo_[i].get(volVars);
527 }
528 }
529
530 // velocity output
531 if (velocityOutput_->enableOutput())
532 {
533 auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache());
534 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
535
536 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
537 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
538 }
539
541 if (addProcessRank)
542 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
543 }
544
548
549 // volume variables if any
550 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
551 sequenceWriter_.addVertexData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarScalarData[i],
552 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/dim, /*nonconforming*/dm_).get() );
553
554 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
555 sequenceWriter_.addVertexData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarVectorData[i],
556 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/dim, /*nonconforming*/dm_).get() );
557
558 // the velocity field
559 if (velocityOutput_->enableOutput())
560 {
561 // node-wise velocities
562 if (dim > 1)
563 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
564 sequenceWriter_.addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
565 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
566 /*numComp*/dimWorld, /*codim*/dim).get() );
567
568 // cell-wise velocities
569 else
570 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
571 sequenceWriter_.addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
572 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
573 /*numComp*/dimWorld, /*codim*/0).get() );
574 }
575
576 // the process rank
577 if (addProcessRank)
578 sequenceWriter_.addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), rank, "process rank", 1, 0).get() );
579
580 // also register additional (non-standardized) user fields if any
581 for (auto&& field : fields_)
582 {
583 if (field.codim() == 0)
584 sequenceWriter_.addCellData(field.get());
585 else if (field.codim() == dim)
586 sequenceWriter_.addVertexData(field.get());
587 else
588 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
589 }
590 }
591
595 sequenceWriter_.write(time, type);
596
600 writer_->clear();
601 }
602
604 template<class Vector, typename std::enable_if_t<IsIndexable<decltype(std::declval<Vector>()[0])>::value, int> = 0>
605 std::size_t getNumberOfComponents_(const Vector& v) { return v[0].size(); }
606
608 template<class Vector, typename std::enable_if_t<!IsIndexable<decltype(std::declval<Vector>()[0])>::value, int> = 0>
609 std::size_t getNumberOfComponents_(const Vector& v) { return 1; }
610
612 std::size_t numDofs_() const { return dofCodim == dim ? gridGeometry().vertexMapper().size() : gridGeometry().elementMapper().size(); }
613
614 const GridVariables& gridVariables_;
615 const SolutionVector& sol_;
616
617 std::string name_;
618 std::string paramGroup_;
619 bool verbose_;
620 Dune::VTK::DataMode dm_;
621
622 std::shared_ptr<Dune::VTKWriter<GridView>> writer_;
623 Dune::VTKSequenceWriter<GridView> sequenceWriter_;
624
625 std::vector<VolVarScalarDataInfo> volVarScalarDataInfo_;
626 std::vector<VolVarVectorDataInfo> volVarVectorDataInfo_;
627
628 std::vector<Field> fields_;
629 std::shared_ptr<VelocityOutput> velocityOutput_;
630};
631
632} // end namespace Dumux
633
634#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods 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:38
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Definition: common/properties/model.hh:34
void resize(const Vector &v, std::size_t size)
Definition: test_isvalid.cc:26
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:66
const auto & problem() const
Definition: io/vtkoutputmodule.hh:232
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:234
void addVolumeVariable(std::function< VolVarsVector(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:149
VV VolumeVariables
export type of the volume variables for the outputfields
Definition: io/vtkoutputmodule.hh:93
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Writing data.
Definition: io/vtkoutputmodule.hh:210
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:119
void addVolumeVariable(std::function< Scalar(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:139
bool verbose() const
Definition: io/vtkoutputmodule.hh:237
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:242
void addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:133
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:249
void addField(const Vector &v, const std::string &name, FieldType fieldType=FieldType::automatic)
Definition: io/vtkoutputmodule.hh:161
const std::string & name() const
Definition: io/vtkoutputmodule.hh:238
Dune::VTK::DataMode dataMode() const
Definition: io/vtkoutputmodule.hh:239
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:241
FieldType
export field type
Definition: io/vtkoutputmodule.hh:97
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:235
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:233
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:246
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)
Definition: io/vtkoutputmodule.hh:101
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:244
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:245
Velocity output for porous media models.