3.2-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 DUMUX_IO_VTK_OUTPUT_MODULE_HH
25#define DUMUX_IO_VTK_OUTPUT_MODULE_HH
26
27#include <functional>
28#include <memory>
29#include <string>
30
31#include <dune/common/timer.hh>
32#include <dune/common/fvector.hh>
33#include <dune/common/typetraits.hh>
34
35#include <dune/geometry/type.hh>
36#include <dune/geometry/multilineargeometry.hh>
37
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>
43
47
48#include "vtkfunction.hh"
49#include "velocityoutput.hh"
50
51namespace Dumux {
52
58template<class GridGeometry>
60{
61 using GridView = typename GridGeometry::GridView;
62 using Field = Vtk::template Field<GridView>;
63 static constexpr int dim = GridView::dimension;
64
65public:
67 enum class FieldType : unsigned int
68 {
69 element, vertex, automatic
70 };
71
72 VtkOutputModuleBase(const GridGeometry& gridGeometry,
73 const std::string& name,
74 const std::string& paramGroup = "",
75 Dune::VTK::DataMode dm = Dune::VTK::conforming,
76 bool verbose = true)
77 : gridGeometry_(gridGeometry)
78 , name_(name)
79 , paramGroup_(paramGroup)
80 , dm_(dm)
81 , verbose_(gridGeometry.gridView().comm().rank() == 0 && verbose)
82 {
83 const auto precisionString = getParamFromGroup<std::string>(paramGroup, "Vtk.Precision", "Float32");
84 precision_ = Dumux::Vtk::stringToPrecision(precisionString);
85 const auto coordPrecision = Dumux::Vtk::stringToPrecision(getParamFromGroup<std::string>(paramGroup, "Vtk.CoordPrecision", precisionString));
86#if DUNE_VERSION_LT(DUNE_GRID, 2, 7)
87 if (precision_ != Dumux::Vtk::Precision::float32 || coordPrecision != Dumux::Vtk::Precision::float32)
88 std::cerr << "Warning: Specifying VTK output precision other than Float32 is only supported in Dune 2.7 and newer. "
89 << "Ignoring parameter and defaulting to Float32." << std::endl;
90 writer_ = std::make_shared<Dune::VTKWriter<GridView>>(gridGeometry.gridView(), dm);
91#else
92 writer_ = std::make_shared<Dune::VTKWriter<GridView>>(gridGeometry.gridView(), dm, coordPrecision);
93#endif
94 sequenceWriter_ = std::make_unique<Dune::VTKSequenceWriter<GridView>>(writer_, name);
95 }
96
97 virtual ~VtkOutputModuleBase() = default;
98
100 const std::string& paramGroup() const
101 { return paramGroup_; }
102
112 template<typename Vector>
113 void addField(const Vector& v,
114 const std::string& name,
116 { addField(v, name, this->precision(), fieldType); }
117
128 template<typename Vector>
129 void addField(const Vector& v,
130 const std::string& name,
131 Dumux::Vtk::Precision precision,
133 {
134 // Deduce the number of components from the given vector type
135 const auto nComp = getNumberOfComponents_(v);
136
137 const auto numElemDofs = gridGeometry().elementMapper().size();
138 const auto numVertexDofs = gridGeometry().vertexMapper().size();
139
140 // Automatically deduce the field type ...
141 if(fieldType == FieldType::automatic)
142 {
143 if(numElemDofs == numVertexDofs)
144 DUNE_THROW(Dune::InvalidStateException, "Automatic deduction of FieldType failed. Please explicitly specify FieldType::element or FieldType::vertex.");
145
146 if(v.size() == numElemDofs)
147 fieldType = FieldType::element;
148 else if(v.size() == numVertexDofs)
149 fieldType = FieldType::vertex;
150 else
151 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
152 }
153 // ... or check if the user-specified type matches the size of v
154 else
155 {
156 if(fieldType == FieldType::element)
157 if(v.size() != numElemDofs)
158 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
159
160 if(fieldType == FieldType::vertex)
161 if(v.size() != numVertexDofs)
162 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
163 }
164
165 // add the appropriate field
166 if (fieldType == FieldType::element)
167 fields_.emplace_back(gridGeometry_.gridView(), gridGeometry_.elementMapper(), v, name, nComp, 0, dm_, precision);
168 else
169 fields_.emplace_back(gridGeometry_.gridView(), gridGeometry_.vertexMapper(), v, name, nComp, dim, dm_, precision);
170 }
171
177 void write(double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
178 {
179 Dune::Timer timer;
180
181 // write to file depending on data mode
182 if (dm_ == Dune::VTK::conforming)
183 writeConforming_(time, type);
184 else if (dm_ == Dune::VTK::nonconforming)
185 writeNonConforming_(time, type);
186 else
187 DUNE_THROW(Dune::NotImplemented, "Output for provided VTK data mode");
188
190 timer.stop();
191 if (verbose_)
192 {
193 std::cout << "Writing output for problem \"" << name_ << "\". Took " << timer.elapsed() << " seconds." << std::endl;
194 }
195 }
196
197protected:
198 const GridGeometry& gridGeometry() const { return gridGeometry_; }
199
200 bool verbose() const { return verbose_; }
201 const std::string& name() const { return name_; }
202 Dune::VTK::DataMode dataMode() const { return dm_; }
203 Dumux::Vtk::Precision precision() const { return precision_; }
204
205 Dune::VTKWriter<GridView>& writer() { return *writer_; }
206 Dune::VTKSequenceWriter<GridView>& sequenceWriter() { return *sequenceWriter_; }
207
208 const std::vector<Field>& fields() const { return fields_; }
209
210private:
212 virtual void writeConforming_(double time, Dune::VTK::OutputType type)
213 {
217
218 // process rank
219 static bool addProcessRank = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddProcessRank");
220 std::vector<int> rank;
221
223 if (!fields_.empty() || addProcessRank)
224 {
225 const auto numCells = gridGeometry_.gridView().size(0);
226
227 // maybe allocate space for the process rank
228 if (addProcessRank)
229 {
230 rank.resize(numCells);
231
232 for (const auto& element : elements(gridGeometry_.gridView(), Dune::Partitions::interior))
233 {
234 const auto eIdxGlobal = gridGeometry_.elementMapper().index(element);
235 rank[eIdxGlobal] = gridGeometry_.gridView().comm().rank();
236 }
237 }
238
242
243 // the process rank
244 if (addProcessRank)
245 this->sequenceWriter().addCellData(Field(gridGeometry_.gridView(), gridGeometry_.elementMapper(), rank, "process rank", 1, 0).get());
246
247 // also register additional (non-standardized) user fields if any
248 for (auto&& field : fields_)
249 {
250 if (field.codim() == 0)
251 this->sequenceWriter().addCellData(field.get());
252 else if (field.codim() == dim)
253 this->sequenceWriter().addVertexData(field.get());
254 else
255 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
256 }
257 }
258
262 this->sequenceWriter().write(time, type);
263
267 this->writer().clear();
268 }
269
271 virtual void writeNonConforming_(double time, Dune::VTK::OutputType type)
272 {
273 DUNE_THROW(Dune::NotImplemented, "Non-conforming VTK output");
274 }
275
277 template<class Vector>
278 std::size_t getNumberOfComponents_(const Vector& v)
279 {
280 if constexpr (IsIndexable<decltype(std::declval<Vector>()[0])>::value)
281 return v[0].size();
282 else
283 return 1;
284 }
285
286 const GridGeometry& gridGeometry_;
287 std::string name_;
288 const std::string paramGroup_;
289 Dune::VTK::DataMode dm_;
290 bool verbose_;
291 Dumux::Vtk::Precision precision_;
292
293 std::shared_ptr<Dune::VTKWriter<GridView>> writer_;
294 std::unique_ptr<Dune::VTKSequenceWriter<GridView>> sequenceWriter_;
295
296 std::vector<Field> fields_;
297};
298
311template<class GridVariables, class SolutionVector>
312class VtkOutputModule : public VtkOutputModuleBase<typename GridVariables::GridGeometry>
313{
315 using GridGeometry = typename GridVariables::GridGeometry;
316
317 using VV = typename GridVariables::VolumeVariables;
318 using Scalar = typename GridVariables::Scalar;
319
320 using GridView = typename GridGeometry::GridView;
321
322 enum {
323 dim = GridView::dimension,
324 dimWorld = GridView::dimensionworld
325 };
326
327 using Element = typename GridView::template Codim<0>::Entity;
328 using VolVarsVector = Dune::FieldVector<Scalar, dimWorld>;
329
330 static constexpr bool isBox = GridGeometry::discMethod == DiscretizationMethod::box;
331 static constexpr int dofCodim = isBox ? dim : 0;
332
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_; };
335 using Field = Vtk::template Field<GridView>;
336
338
339public:
341 using VolumeVariables = VV;
342
343 VtkOutputModule(const GridVariables& gridVariables,
344 const SolutionVector& sol,
345 const std::string& name,
346 const std::string& paramGroup = "",
347 Dune::VTK::DataMode dm = Dune::VTK::conforming,
348 bool verbose = true)
350 , gridVariables_(gridVariables)
351 , sol_(sol)
352 , velocityOutput_(std::make_shared<VelocityOutputType>())
353 {}
354
359
366 void addVelocityOutput(std::shared_ptr<VelocityOutputType> velocityOutput)
367 { velocityOutput_ = velocityOutput; }
368
372 void addVolumeVariable(std::function<Scalar(const VolumeVariables&)>&& f,
373 const std::string& name)
374 {
375 volVarScalarDataInfo_.push_back(VolVarScalarDataInfo{f, name, this->precision()});
376 }
377
382 template<class VVV = VolVarsVector, typename std::enable_if_t<(VVV::dimension > 1), int> = 0>
383 void addVolumeVariable(std::function<VolVarsVector(const VolumeVariables&)>&& f,
384 const std::string& name)
385 {
386 volVarVectorDataInfo_.push_back(VolVarVectorDataInfo{f, name, this->precision()});
387 }
388
389protected:
390 // some return functions for differing implementations to use
391 const auto& problem() const { return gridVariables_.curGridVolVars().problem(); }
392 const GridVariables& gridVariables() const { return gridVariables_; }
393 const GridGeometry& gridGeometry() const { return gridVariables_.gridGeometry(); }
394 const SolutionVector& sol() const { return sol_; }
395
396 const std::vector<VolVarScalarDataInfo>& volVarScalarDataInfo() const { return volVarScalarDataInfo_; }
397 const std::vector<VolVarVectorDataInfo>& volVarVectorDataInfo() const { return volVarVectorDataInfo_; }
398
400 const VelocityOutput& velocityOutput() const { return *velocityOutput_; }
401
402private:
403
405 void writeConforming_(double time, Dune::VTK::OutputType type) override
406 {
407 const Dune::VTK::DataMode dm = Dune::VTK::conforming;
411
412 // instantiate the velocity output
413 using VelocityVector = typename VelocityOutput::VelocityVector;
414 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
415
416 // process rank
417 static bool addProcessRank = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddProcessRank");
418 std::vector<double> rank;
419
420 // volume variable data
421 std::vector<std::vector<Scalar>> volVarScalarData;
422 std::vector<std::vector<VolVarsVector>> volVarVectorData;
423
425 if (!volVarScalarDataInfo_.empty()
426 || !volVarVectorDataInfo_.empty()
427 || !this->fields().empty()
428 || velocityOutput_->enableOutput()
429 || addProcessRank)
430 {
431 const auto numCells = gridGeometry().gridView().size(0);
432 const auto numDofs = numDofs_();
433
434 // get fields for all volume variables
435 if (!volVarScalarDataInfo_.empty())
436 volVarScalarData.resize(volVarScalarDataInfo_.size(), std::vector<Scalar>(numDofs));
437 if (!volVarVectorDataInfo_.empty())
438 volVarVectorData.resize(volVarVectorDataInfo_.size(), std::vector<VolVarsVector>(numDofs));
439
440 if (velocityOutput_->enableOutput())
441 {
442 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
443 {
444 if(isBox && dim == 1)
445 velocity[phaseIdx].resize(numCells);
446 else
447 velocity[phaseIdx].resize(numDofs);
448 }
449 }
450
451 // maybe allocate space for the process rank
452 if (addProcessRank) rank.resize(numCells);
453
454 for (const auto& element : elements(gridGeometry().gridView(), Dune::Partitions::interior))
455 {
456 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
457
458 auto fvGeometry = localView(gridGeometry());
459 auto elemVolVars = localView(gridVariables_.curGridVolVars());
460
461 // If velocity output is enabled we need to bind to the whole stencil
462 // otherwise element-local data is sufficient
463 if (velocityOutput_->enableOutput())
464 {
465 fvGeometry.bind(element);
466 elemVolVars.bind(element, fvGeometry, sol_);
467 }
468 else
469 {
470 fvGeometry.bindElement(element);
471 elemVolVars.bindElement(element, fvGeometry, sol_);
472 }
473
474 if (!volVarScalarDataInfo_.empty()
475 || !volVarVectorDataInfo_.empty())
476 {
477 for (auto&& scv : scvs(fvGeometry))
478 {
479 const auto dofIdxGlobal = scv.dofIndex();
480 const auto& volVars = elemVolVars[scv];
481
482 // get the scalar-valued data
483 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
484 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo_[i].get(volVars);
485
486 // get the vector-valued data
487 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
488 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo_[i].get(volVars);
489 }
490 }
491
492 // velocity output
493 if (velocityOutput_->enableOutput())
494 {
495 auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache());
496 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
497
498 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
499 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
500 }
501
503 if (addProcessRank)
504 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
505 }
506
510
511 // volume variables if any
512 if (isBox)
513 {
514 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
515 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), volVarScalarData[i],
516 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/dim, dm, this->precision()).get() );
517 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
518 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), volVarVectorData[i],
519 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()).get() );
520 }
521 else
522 {
523 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
524 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarScalarData[i],
525 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0,dm, this->precision()).get() );
526 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
527 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarVectorData[i],
528 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0,dm, this->precision()).get() );
529 }
530
531 // the velocity field
532 if (velocityOutput_->enableOutput())
533 {
534 if (isBox && dim > 1)
535 {
536 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
537 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
538 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
539 /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()).get() );
540 }
541 // cell-centered models
542 else
543 {
544 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
545 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
546 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
547 /*numComp*/dimWorld, /*codim*/0, dm, this->precision()).get() );
548 }
549 }
550
551 // the process rank
552 if (addProcessRank)
553 this->sequenceWriter().addCellData(Field(gridGeometry().gridView(), gridGeometry().elementMapper(), rank, "process rank", 1, 0).get());
554
555 // also register additional (non-standardized) user fields if any
556 for (auto&& field : this->fields())
557 {
558 if (field.codim() == 0)
559 this->sequenceWriter().addCellData(field.get());
560 else if (field.codim() == dim)
561 this->sequenceWriter().addVertexData(field.get());
562 else
563 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
564 }
565 }
566
570 this->sequenceWriter().write(time, type);
571
575 this->writer().clear();
576 }
577
579 void writeNonConforming_(double time, Dune::VTK::OutputType type) override
580 {
581 const Dune::VTK::DataMode dm = Dune::VTK::nonconforming;
582
583 if(!isBox)
584 DUNE_THROW(Dune::InvalidStateException, "Non-conforming output makes no sense for cell-centered schemes!");
585
589
590 // check the velocity output
591 bool enableVelocityOutput = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddVelocity");
592 if (enableVelocityOutput == true && !velocityOutput_->enableOutput())
593 std::cerr << "Warning! Velocity output was enabled in the input file"
594 << " but no velocity output policy was set for the VTK output module:"
595 << " There will be no velocity output."
596 << " Use the addVelocityOutput member function of the VTK output module." << std::endl;
597 using VelocityVector = typename VelocityOutput::VelocityVector;
598 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
599
600 // process rank
601 static bool addProcessRank = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddProcessRank");
602 std::vector<double> rank;
603
604 // volume variable data (indexing: volvardata/element/localcorner)
605 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
606 using VectorDataContainer = std::vector< std::vector<VolVarsVector> >;
607 std::vector< ScalarDataContainer > volVarScalarData;
608 std::vector< VectorDataContainer > volVarVectorData;
609
611 if (!volVarScalarDataInfo_.empty()
612 || !volVarVectorDataInfo_.empty()
613 || !this->fields().empty()
614 || velocityOutput_->enableOutput()
615 || addProcessRank)
616 {
617 const auto numCells = gridGeometry().gridView().size(0);
618 const auto numDofs = numDofs_();
619
620 // get fields for all volume variables
621 if (!volVarScalarDataInfo_.empty())
622 volVarScalarData.resize(volVarScalarDataInfo_.size(), ScalarDataContainer(numCells));
623 if (!volVarVectorDataInfo_.empty())
624 volVarVectorData.resize(volVarVectorDataInfo_.size(), VectorDataContainer(numCells));
625
626 if (velocityOutput_->enableOutput())
627 {
628 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
629 {
630 if(isBox && dim == 1)
631 velocity[phaseIdx].resize(numCells);
632 else
633 velocity[phaseIdx].resize(numDofs);
634 }
635 }
636
637 // maybe allocate space for the process rank
638 if (addProcessRank) rank.resize(numCells);
639
640 for (const auto& element : elements(gridGeometry().gridView(), Dune::Partitions::interior))
641 {
642 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
643 const auto numCorners = element.subEntities(dim);
644
645 auto fvGeometry = localView(gridGeometry());
646 auto elemVolVars = localView(gridVariables_.curGridVolVars());
647
648 // resize element-local data containers
649 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
650 volVarScalarData[i][eIdxGlobal].resize(numCorners);
651 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
652 volVarVectorData[i][eIdxGlobal].resize(numCorners);
653
654 // If velocity output is enabled we need to bind to the whole stencil
655 // otherwise element-local data is sufficient
656 if (velocityOutput_->enableOutput())
657 {
658 fvGeometry.bind(element);
659 elemVolVars.bind(element, fvGeometry, sol_);
660 }
661 else
662 {
663 fvGeometry.bindElement(element);
664 elemVolVars.bindElement(element, fvGeometry, sol_);
665 }
666
667 if (!volVarScalarDataInfo_.empty()
668 || !volVarVectorDataInfo_.empty())
669 {
670 for (auto&& scv : scvs(fvGeometry))
671 {
672 const auto& volVars = elemVolVars[scv];
673
674 // get the scalar-valued data
675 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
676 volVarScalarData[i][scv.elementIndex()][scv.localDofIndex()] = volVarScalarDataInfo_[i].get(volVars);
677
678 // get the vector-valued data
679 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
680 volVarVectorData[i][scv.elementIndex()][scv.localDofIndex()] = volVarVectorDataInfo_[i].get(volVars);
681 }
682 }
683
684 // velocity output
685 if (velocityOutput_->enableOutput())
686 {
687 auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache());
688 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
689
690 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
691 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
692 }
693
695 if (addProcessRank)
696 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
697 }
698
702
703 // volume variables if any
704 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
705 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarScalarData[i],
706 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/dim, /*nonconforming*/dm, this->precision()).get() );
707
708 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
709 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarVectorData[i],
710 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/dim, /*nonconforming*/dm, this->precision()).get() );
711
712 // the velocity field
713 if (velocityOutput_->enableOutput())
714 {
715 // node-wise velocities
716 if (dim > 1)
717 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
718 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
719 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
720 /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()).get() );
721
722 // cell-wise velocities
723 else
724 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
725 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
726 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
727 /*numComp*/dimWorld, /*codim*/0,dm, this->precision()).get());
728 }
729
730 // the process rank
731 if (addProcessRank)
732 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), rank, "process rank", 1, 0).get() );
733
734 // also register additional (non-standardized) user fields if any
735 for (auto&& field : this->fields())
736 {
737 if (field.codim() == 0)
738 this->sequenceWriter().addCellData(field.get());
739 else if (field.codim() == dim)
740 this->sequenceWriter().addVertexData(field.get());
741 else
742 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
743 }
744 }
745
749 this->sequenceWriter().write(time, type);
750
754 this->writer().clear();
755 }
756
758 std::size_t numDofs_() const { return dofCodim == dim ? gridGeometry().vertexMapper().size() : gridGeometry().elementMapper().size(); }
759
760 const GridVariables& gridVariables_;
761 const SolutionVector& sol_;
762
763 std::vector<VolVarScalarDataInfo> volVarScalarDataInfo_;
764 std::vector<VolVarVectorDataInfo> volVarVectorDataInfo_;
765
766 std::shared_ptr<VelocityOutput> velocityOutput_;
767};
768
769} // end namespace Dumux
770
771#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
Precision stringToPrecision(std::string_view precisionName)
Maps a string (e.g. from input) to a Dune precision type.
Definition: vtkprecision.hh:53
Definition: adapt.hh:29
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:60
Dumux::Vtk::Precision precision() const
Definition: io/vtkoutputmodule.hh:203
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:100
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 &paramGroup="", Dune::VTK::DataMode dm=Dune::VTK::conforming, bool verbose=true)
Definition: io/vtkoutputmodule.hh:72
FieldType
export field type
Definition: io/vtkoutputmodule.hh:68
Dune::VTK::DataMode dataMode() const
Definition: io/vtkoutputmodule.hh:202
void addField(const Vector &v, const std::string &name, Dumux::Vtk::Precision precision, FieldType fieldType=FieldType::automatic)
Add a scalar or vector valued vtk field.
Definition: io/vtkoutputmodule.hh:129
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:205
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:206
const std::string & name() const
Definition: io/vtkoutputmodule.hh:201
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:208
virtual ~VtkOutputModuleBase()=default
bool verbose() const
Definition: io/vtkoutputmodule.hh:200
void addField(const Vector &v, const std::string &name, FieldType fieldType=FieldType::automatic)
Add a scalar or vector valued vtk field.
Definition: io/vtkoutputmodule.hh:113
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:198
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:313
const auto & problem() const
Definition: io/vtkoutputmodule.hh:391
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:392
void addVolumeVariable(std::function< VolVarsVector(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:383
VV VolumeVariables
export type of the volume variables for the outputfields
Definition: io/vtkoutputmodule.hh:341
void addVolumeVariable(std::function< Scalar(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:372
void addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:366
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:400
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:394
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:393
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:343
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:396
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:397
Velocity output for porous media models.