3.6-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
45#include <dumux/io/format.hh>
47
50#include "velocityoutput.hh"
51
52namespace Dumux {
53
59template<class GridGeometry>
61{
62 using GridView = typename GridGeometry::GridView;
63 static constexpr int dim = GridView::dimension;
64
65public:
67 using Field = Vtk::template Field<GridView>;
68
69 VtkOutputModuleBase(const GridGeometry& gridGeometry,
70 const std::string& name,
71 const std::string& paramGroup = "",
72 Dune::VTK::DataMode dm = Dune::VTK::conforming,
73 bool verbose = true)
74 : gridGeometry_(gridGeometry)
75 , name_(name)
76 , paramGroup_(paramGroup)
77 , dm_(dm)
78 , verbose_(gridGeometry.gridView().comm().rank() == 0 && verbose)
79 {
80 const auto precisionString = getParamFromGroup<std::string>(paramGroup, "Vtk.Precision", "Float32");
81 precision_ = Dumux::Vtk::stringToPrecision(precisionString);
82 const auto coordPrecision = Dumux::Vtk::stringToPrecision(getParamFromGroup<std::string>(paramGroup, "Vtk.CoordPrecision", precisionString));
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);
86 }
87
88 virtual ~VtkOutputModuleBase() = default;
89
91 const std::string& paramGroup() const
92 { return paramGroup_; }
93
103 template<typename Vector>
104 void addField(const Vector& v,
105 const std::string& name,
107 { addField(v, name, this->precision(), fieldType); }
108
119 template<typename Vector>
120 void addField(const Vector& v,
121 const std::string& name,
122 Dumux::Vtk::Precision precision,
124 {
125 // Deduce the number of components from the given vector type
126 const auto nComp = getNumberOfComponents_(v);
127
128 const auto numElemDofs = gridGeometry().elementMapper().size();
129 const auto numVertexDofs = gridGeometry().vertexMapper().size();
130
131 // Automatically deduce the field type ...
132 if(fieldType == Vtk::FieldType::automatic)
133 {
134 if(numElemDofs == numVertexDofs)
135 DUNE_THROW(Dune::InvalidStateException, "Automatic deduction of FieldType failed. Please explicitly specify FieldType::element or FieldType::vertex.");
136
137 if(v.size() == numElemDofs)
138 fieldType = Vtk::FieldType::element;
139 else if(v.size() == numVertexDofs)
140 fieldType = Vtk::FieldType::vertex;
141 else
142 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
143 }
144 // ... or check if the user-specified type matches the size of v
145 else
146 {
147 if(fieldType == Vtk::FieldType::element)
148 if(v.size() != numElemDofs)
149 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
150
151 if(fieldType == Vtk::FieldType::vertex)
152 if(v.size() != numVertexDofs)
153 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
154 }
155
156 // add the appropriate field
157 if (fieldType == Vtk::FieldType::element)
158 fields_.emplace_back(gridGeometry_.gridView(), gridGeometry_.elementMapper(), v, name, nComp, 0, dm_, precision);
159 else
160 fields_.emplace_back(gridGeometry_.gridView(), gridGeometry_.vertexMapper(), v, name, nComp, dim, dm_, precision);
161 }
162
167 void addField(Field&& field)
168 {
169 fields_.push_back(std::move(field));
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 std::cout << Fmt::format("Writing output for problem \"{}\". Took {:.2g} seconds.\n", name_, timer.elapsed());
193 }
194
195protected:
196 const GridGeometry& gridGeometry() const { return gridGeometry_; }
197
198 bool verbose() const { return verbose_; }
199 const std::string& name() const { return name_; }
200 Dune::VTK::DataMode dataMode() const { return dm_; }
201 Dumux::Vtk::Precision precision() const { return precision_; }
202
203 Dune::VTKWriter<GridView>& writer() { return *writer_; }
204 Dune::VTKSequenceWriter<GridView>& sequenceWriter() { return *sequenceWriter_; }
205
206 const std::vector<Field>& fields() const { return fields_; }
207
208private:
210 virtual void writeConforming_(double time, Dune::VTK::OutputType type)
211 {
215
216 // process rank
217 std::vector<int> rank;
218
220 if (!fields_.empty() || addProcessRank_)
221 {
222 const auto numCells = gridGeometry_.gridView().size(0);
223
224 // maybe allocate space for the process rank
225 if (addProcessRank_)
226 {
227 rank.resize(numCells);
228
229 for (const auto& element : elements(gridGeometry_.gridView(), Dune::Partitions::interior))
230 {
231 const auto eIdxGlobal = gridGeometry_.elementMapper().index(element);
232 rank[eIdxGlobal] = gridGeometry_.gridView().comm().rank();
233 }
234 }
235
239
240 // the process rank
241 if (addProcessRank_)
242 this->sequenceWriter().addCellData(Field(gridGeometry_.gridView(), gridGeometry_.elementMapper(), rank, "process rank", 1, 0).get());
243
244 // also register additional (non-standardized) user fields if any
245 for (auto&& field : fields_)
246 {
247 if (field.codim() == 0)
248 this->sequenceWriter().addCellData(field.get());
249 else if (field.codim() == dim)
250 this->sequenceWriter().addVertexData(field.get());
251 else
252 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
253 }
254 }
255
259 this->sequenceWriter().write(time, type);
260
264 this->writer().clear();
265 }
266
268 virtual void writeNonConforming_(double time, Dune::VTK::OutputType type)
269 {
270 DUNE_THROW(Dune::NotImplemented, "Non-conforming VTK output");
271 }
272
274 template<class Vector>
275 std::size_t getNumberOfComponents_(const Vector& v)
276 {
277 if constexpr (Dune::IsIndexable<decltype(std::declval<Vector>()[0])>::value)
278 return v[0].size();
279 else
280 return 1;
281 }
282
283 const GridGeometry& gridGeometry_;
284 std::string name_;
285 const std::string paramGroup_;
286 Dune::VTK::DataMode dm_;
287 bool verbose_;
288 Dumux::Vtk::Precision precision_;
289
290 std::shared_ptr<Dune::VTKWriter<GridView>> writer_;
291 std::unique_ptr<Dune::VTKSequenceWriter<GridView>> sequenceWriter_;
292
293 std::vector<Field> fields_;
294
295 bool addProcessRank_ = true;
296};
297
310template<class GridVariables, class SolutionVector>
311class VtkOutputModule : public VtkOutputModuleBase<typename GridVariables::GridGeometry>
312{
314 using GridGeometry = typename GridVariables::GridGeometry;
315
316 using VV = typename GridVariables::VolumeVariables;
317 using Scalar = typename GridVariables::Scalar;
318
319 using GridView = typename GridGeometry::GridView;
320
321 enum {
322 dim = GridView::dimension,
323 dimWorld = GridView::dimensionworld
324 };
325
326 using Element = typename GridView::template Codim<0>::Entity;
327 using VolVarsVector = Dune::FieldVector<Scalar, dimWorld>;
328
329 static constexpr bool isBox = GridGeometry::discMethod == DiscretizationMethods::box;
330 static constexpr bool isDiamond = GridGeometry::discMethod == DiscretizationMethods::fcdiamond;
331 static constexpr bool isPQ1Bubble = GridGeometry::discMethod == DiscretizationMethods::pq1bubble;
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
337
338public:
340 using Field = Vtk::template Field<GridView>;
342 using VolumeVariables = VV;
343
344 VtkOutputModule(const GridVariables& gridVariables,
345 const SolutionVector& sol,
346 const std::string& name,
347 const std::string& paramGroup = "",
348 Dune::VTK::DataMode dm = Dune::VTK::conforming,
349 bool verbose = true)
351 , gridVariables_(gridVariables)
352 , sol_(sol)
353 , velocityOutput_(std::make_shared<VelocityOutputType>())
354 {
355 enableVelocityOutput_ = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddVelocity", false);
356 addProcessRank_ = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddProcessRank", true);
357 }
358
363
370 void addVelocityOutput(std::shared_ptr<VelocityOutputType> velocityOutput)
371 { velocityOutput_ = velocityOutput; }
372
376 void addVolumeVariable(std::function<Scalar(const VolumeVariables&)>&& f,
377 const std::string& name)
378 {
379 volVarScalarDataInfo_.push_back(VolVarScalarDataInfo{f, name, this->precision()});
380 }
381
386 template<class VVV = VolVarsVector, typename std::enable_if_t<(VVV::dimension > 1), int> = 0>
387 void addVolumeVariable(std::function<VolVarsVector(const VolumeVariables&)>&& f,
388 const std::string& name)
389 {
390 volVarVectorDataInfo_.push_back(VolVarVectorDataInfo{f, name, this->precision()});
391 }
392
393protected:
394 // some return functions for differing implementations to use
395 const auto& problem() const { return gridVariables_.curGridVolVars().problem(); }
396 const GridVariables& gridVariables() const { return gridVariables_; }
397 const GridGeometry& gridGeometry() const { return gridVariables_.gridGeometry(); }
398 const SolutionVector& sol() const { return sol_; }
399
400 const std::vector<VolVarScalarDataInfo>& volVarScalarDataInfo() const { return volVarScalarDataInfo_; }
401 const std::vector<VolVarVectorDataInfo>& volVarVectorDataInfo() const { return volVarVectorDataInfo_; }
402
404 const VelocityOutput& velocityOutput() const { return *velocityOutput_; }
405
406private:
407
409 void writeConforming_(double time, Dune::VTK::OutputType type) override
410 {
411 const Dune::VTK::DataMode dm = Dune::VTK::conforming;
415
416 // instantiate the velocity output
417 using VelocityVector = typename VelocityOutput::VelocityVector;
418 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
419
420 // process rank
421 std::vector<double> rank;
422
423 // volume variable data
424 std::vector<std::vector<Scalar>> volVarScalarData;
425 std::vector<std::vector<VolVarsVector>> volVarVectorData;
426
428 if (!volVarScalarDataInfo_.empty()
429 || !volVarVectorDataInfo_.empty()
430 || !this->fields().empty()
431 || velocityOutput_->enableOutput()
432 || addProcessRank_)
433 {
434 const auto numCells = gridGeometry().gridView().size(0);
435 const auto numDofs = numDofs_();
436
437 // get fields for all volume variables
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));
442
443 if (velocityOutput_->enableOutput())
444 {
445 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
446 {
447 if (velocityOutput_->fieldType() == VelocityOutput::FieldType::element)
448 velocity[phaseIdx].resize(numCells);
449 else if (velocityOutput_->fieldType() == VelocityOutput::FieldType::vertex)
450 velocity[phaseIdx].resize(numDofs);
451 else
452 {
453 if(isBox && dim == 1)
454 velocity[phaseIdx].resize(numCells);
455 else
456 velocity[phaseIdx].resize(numDofs);
457 }
458 }
459 }
460
461 // maybe allocate space for the process rank
462 if (addProcessRank_) rank.resize(numCells);
463
464 auto fvGeometry = localView(gridGeometry());
465 auto elemVolVars = localView(gridVariables_.curGridVolVars());
466 for (const auto& element : elements(gridGeometry().gridView(), Dune::Partitions::interior))
467 {
468 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
469 // If velocity output is enabled we need to bind to the whole stencil
470 // otherwise element-local data is sufficient
471 if (velocityOutput_->enableOutput())
472 {
473 fvGeometry.bind(element);
474 elemVolVars.bind(element, fvGeometry, sol_);
475 }
476 else
477 {
478 fvGeometry.bindElement(element);
479 elemVolVars.bindElement(element, fvGeometry, sol_);
480 }
481
482 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
483 {
484 for (const auto& scv : scvs(fvGeometry))
485 {
486 const auto dofIdxGlobal = scv.dofIndex();
487 const auto& volVars = elemVolVars[scv];
488
489 // get the scalar-valued data
490 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
491 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo_[i].get(volVars);
492
493 // get the vector-valued data
494 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
495 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo_[i].get(volVars);
496 }
497 }
498
499 // velocity output
500 if (velocityOutput_->enableOutput())
501 {
502 const auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
503
504 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
505 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
506 }
507
509 if (addProcessRank_)
510 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
511 }
512
516
517 // volume variables if any
518 if constexpr (isBox || isPQ1Bubble)
519 {
520 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
521 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarScalarData[i],
522 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/dim, dm, this->precision()).get() );
523 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
524 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarVectorData[i],
525 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()).get() );
526
527 if constexpr (isPQ1Bubble)
528 {
529 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
530 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarScalarData[i],
531 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0,dm, this->precision()).get() );
532 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
533 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarVectorData[i],
534 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0,dm, this->precision()).get() );
535 }
536
537 }
538 else
539 {
540 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
541 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarScalarData[i],
542 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0,dm, this->precision()).get() );
543 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
544 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarVectorData[i],
545 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0,dm, this->precision()).get() );
546 }
547
548 // the velocity field
549 if (velocityOutput_->enableOutput())
550 {
551 if (velocityOutput_->fieldType() == VelocityOutput::FieldType::vertex
552 || ( (velocityOutput_->fieldType() == VelocityOutput::FieldType::automatic) && dim > 1 && isBox ))
553 {
554 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
555 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
556 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
557 /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()).get() );
558 }
559 // cell-centered models
560 else
561 {
562 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
563 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
564 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
565 /*numComp*/dimWorld, /*codim*/0, dm, this->precision()).get() );
566 }
567 }
568
569 // the process rank
570 if (addProcessRank_)
571 this->sequenceWriter().addCellData(Field(gridGeometry().gridView(), gridGeometry().elementMapper(), rank, "process rank", 1, 0).get());
572
573 // also register additional (non-standardized) user fields if any
574 for (auto&& field : this->fields())
575 {
576 if (field.codim() == 0)
577 this->sequenceWriter().addCellData(field.get());
578 else if (field.codim() == dim)
579 this->sequenceWriter().addVertexData(field.get());
580 else
581 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
582 }
583 }
584
588 this->sequenceWriter().write(time, type);
589
593 this->writer().clear();
594 }
595
597 void writeNonConforming_(double time, Dune::VTK::OutputType type) override
598 {
599 const Dune::VTK::DataMode dm = Dune::VTK::nonconforming;
600
601 // only supports finite-element-like discretization schemes
602 if(!isBox && !isDiamond)
603 DUNE_THROW(Dune::NotImplemented,
604 "Non-conforming output for discretization scheme " << GridGeometry::discMethod
605 );
606
610
611 // check the velocity output
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;
617 using VelocityVector = typename VelocityOutput::VelocityVector;
618 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
619
620 // process rank
621 std::vector<double> rank;
622
623 // volume variable data (indexing: volvardata/element/localindex)
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;
628
630 if (!volVarScalarDataInfo_.empty()
631 || !volVarVectorDataInfo_.empty()
632 || !this->fields().empty()
633 || velocityOutput_->enableOutput()
634 || addProcessRank_)
635 {
636 const auto numCells = gridGeometry().gridView().size(0);
637 const auto outputSize = numDofs_();
638
639 // get fields for all volume variables
640 if (!volVarScalarDataInfo_.empty())
641 volVarScalarData.resize(volVarScalarDataInfo_.size(), ScalarDataContainer(numCells));
642 if (!volVarVectorDataInfo_.empty())
643 volVarVectorData.resize(volVarVectorDataInfo_.size(), VectorDataContainer(numCells));
644
645 if (velocityOutput_->enableOutput())
646 {
647 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
648 {
649 if((isBox && dim == 1) || isDiamond)
650 velocity[phaseIdx].resize(numCells);
651 else
652 velocity[phaseIdx].resize(outputSize);
653 }
654 }
655
656 // maybe allocate space for the process rank
657 if (addProcessRank_) rank.resize(numCells);
658
659 // now we go element-local to extract values at local dof locations
660 auto fvGeometry = localView(gridGeometry());
661 auto elemVolVars = localView(gridVariables_.curGridVolVars());
662 for (const auto& element : elements(gridGeometry().gridView(), Dune::Partitions::interior))
663 {
664 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
665 // If velocity output is enabled we need to bind to the whole stencil
666 // otherwise element-local data is sufficient
667 if (velocityOutput_->enableOutput())
668 {
669 fvGeometry.bind(element);
670 elemVolVars.bind(element, fvGeometry, sol_);
671 }
672 else
673 {
674 fvGeometry.bindElement(element);
675 elemVolVars.bindElement(element, fvGeometry, sol_);
676 }
677
678 const auto numLocalDofs = fvGeometry.numScv();
679 // resize element-local data containers
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);
684
685 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
686 {
687 for (const auto& scv : scvs(fvGeometry))
688 {
689 const auto& volVars = elemVolVars[scv];
690
691 // get the scalar-valued data
692 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
693 volVarScalarData[i][eIdxGlobal][scv.localDofIndex()] = volVarScalarDataInfo_[i].get(volVars);
694
695 // get the vector-valued data
696 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
697 volVarVectorData[i][eIdxGlobal][scv.localDofIndex()] = volVarVectorDataInfo_[i].get(volVars);
698 }
699 }
700
701 // velocity output
702 if (velocityOutput_->enableOutput())
703 {
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);
707 }
708
710 if (addProcessRank_)
711 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
712 }
713
717
718 // volume variables if any
719 static constexpr int dofLocCodim = isDiamond ? 1 : dim;
720 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
721 this->sequenceWriter().addVertexData(Field(
722 gridGeometry().gridView(), gridGeometry().elementMapper(),
723 volVarScalarData[i], volVarScalarDataInfo_[i].name,
724 /*numComp*/1, /*codim*/dofLocCodim, /*nonconforming*/dm, this->precision()
725 ).get());
726
727 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
728 this->sequenceWriter().addVertexData(Field(
729 gridGeometry().gridView(), gridGeometry().elementMapper(),
730 volVarVectorData[i], volVarVectorDataInfo_[i].name,
731 /*numComp*/dimWorld, /*codim*/dofLocCodim, /*nonconforming*/dm, this->precision()
732 ).get());
733
734 // the velocity field
735 if (velocityOutput_->enableOutput())
736 {
737 // node-wise velocities
738 if (dim > 1 && !isDiamond)
739 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
740 this->sequenceWriter().addVertexData(Field(
741 gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
742 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
743 /*numComp*/dimWorld, /*codim*/dofLocCodim, dm, this->precision()
744 ).get());
745
746 // cell-wise velocities
747 else
748 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
749 this->sequenceWriter().addCellData(Field(
750 gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
751 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
752 /*numComp*/dimWorld, /*codim*/0, dm, this->precision()
753 ).get());
754 }
755 }
756
757 // the process rank
758 if (addProcessRank_)
759 this->sequenceWriter().addCellData(Field(
760 gridGeometry().gridView(), gridGeometry().elementMapper(),
761 rank, "process rank", /*numComp*/1, /*codim*/0
762 ).get());
763
764 // also register additional (non-standardized) user fields if any
765 for (const auto& field : this->fields())
766 {
767 if (field.codim() == 0)
768 this->sequenceWriter().addCellData(field.get());
769 else if (field.codim() == dim || field.codim() == 1)
770 this->sequenceWriter().addVertexData(field.get());
771 else
772 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
773 }
774
778 this->sequenceWriter().write(time, type);
779
783 this->writer().clear();
784 }
785
787 std::size_t numDofs_() const
788 {
789 // TODO this should actually always be dofMapper.size()
790 // maybe some discretizations needs special treatment (?)
791 if constexpr (isBox || isDiamond || isPQ1Bubble)
792 return gridGeometry().dofMapper().size();
793 else
794 return gridGeometry().elementMapper().size();
795 }
796
797 const GridVariables& gridVariables_;
798 const SolutionVector& sol_;
799
800 std::vector<VolVarScalarDataInfo> volVarScalarDataInfo_;
801 std::vector<VolVarVectorDataInfo> volVarVectorDataInfo_;
802
803 std::shared_ptr<VelocityOutput> velocityOutput_;
804 bool enableVelocityOutput_ = false;
805 bool addProcessRank_ = true;
806};
807
808} // end namespace Dumux
809
810#endif
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 &paramGroup="", 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 &paramGroup="", 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
TODO: docme!