version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_IO_VTK_OUTPUT_MODULE_HH
13#define DUMUX_IO_VTK_OUTPUT_MODULE_HH
14
15#include <functional>
16#include <memory>
17#include <string>
18
19#include <dune/common/timer.hh>
20#include <dune/common/fvector.hh>
21#include <dune/common/typetraits.hh>
22
23#include <dune/geometry/type.hh>
24#include <dune/geometry/multilineargeometry.hh>
25
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>
31
33#include <dumux/io/format.hh>
35
38#include "velocityoutput.hh"
39
40namespace Dumux {
41
47template<class GridGeometry>
49{
50 using GridView = typename GridGeometry::GridView;
51 static constexpr int dim = GridView::dimension;
52
53public:
55 using Field = Vtk::template Field<GridView>;
56
57 VtkOutputModuleBase(const GridGeometry& gridGeometry,
58 const std::string& name,
59 const std::string& paramGroup = "",
60 Dune::VTK::DataMode dm = Dune::VTK::conforming,
61 bool verbose = true)
62 : gridGeometry_(gridGeometry)
63 , name_(name)
64 , paramGroup_(paramGroup)
65 , dm_(dm)
66 , verbose_(gridGeometry.gridView().comm().rank() == 0 && verbose)
67 {
68 const auto precisionString = getParamFromGroup<std::string>(paramGroup, "Vtk.Precision", "Float32");
69 precision_ = Dumux::Vtk::stringToPrecision(precisionString);
70 const auto coordPrecision = Dumux::Vtk::stringToPrecision(getParamFromGroup<std::string>(paramGroup, "Vtk.CoordPrecision", precisionString));
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);
74 }
75
76 virtual ~VtkOutputModuleBase() = default;
77
79 const std::string& paramGroup() const
80 { return paramGroup_; }
81
91 template<typename Vector>
92 void addField(const Vector& v,
93 const std::string& name,
95 { addField(v, name, this->precision(), fieldType); }
96
107 template<typename Vector>
108 void addField(const Vector& v,
109 const std::string& name,
110 Dumux::Vtk::Precision precision,
112 {
113 // Deduce the number of components from the given vector type
114 const auto nComp = getNumberOfComponents_(v);
115
116 const auto numElemDofs = gridGeometry().elementMapper().size();
117 const auto numVertexDofs = gridGeometry().vertexMapper().size();
118
119 // Automatically deduce the field type ...
120 if(fieldType == Vtk::FieldType::automatic)
121 {
122 if(numElemDofs == numVertexDofs)
123 DUNE_THROW(Dune::InvalidStateException, "Automatic deduction of FieldType failed. Please explicitly specify FieldType::element or FieldType::vertex.");
124
125 if(v.size() == numElemDofs)
126 fieldType = Vtk::FieldType::element;
127 else if(v.size() == numVertexDofs)
128 fieldType = Vtk::FieldType::vertex;
129 else
130 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
131 }
132 // ... or check if the user-specified type matches the size of v
133 else
134 {
135 if(fieldType == Vtk::FieldType::element)
136 if(v.size() != numElemDofs)
137 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
138
139 if(fieldType == Vtk::FieldType::vertex)
140 if(v.size() != numVertexDofs)
141 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
142 }
143
144 // add the appropriate field
145 if (fieldType == Vtk::FieldType::element)
146 fields_.emplace_back(gridGeometry_.gridView(), gridGeometry_.elementMapper(), v, name, nComp, 0, dm_, precision);
147 else
148 fields_.emplace_back(gridGeometry_.gridView(), gridGeometry_.vertexMapper(), v, name, nComp, dim, dm_, precision);
149 }
150
155 void addField(Field&& field)
156 {
157 fields_.push_back(std::move(field));
158 }
159
165 void write(double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
166 {
167 Dune::Timer timer;
168
169 // write to file depending on data mode
170 if (dm_ == Dune::VTK::conforming)
171 writeConforming_(time, type);
172 else if (dm_ == Dune::VTK::nonconforming)
173 writeNonConforming_(time, type);
174 else
175 DUNE_THROW(Dune::NotImplemented, "Output for provided VTK data mode");
176
178 timer.stop();
179 if (verbose_)
180 std::cout << Fmt::format("Writing output for problem \"{}\". Took {:.2g} seconds.\n", name_, timer.elapsed());
181 }
182
183protected:
184 const GridGeometry& gridGeometry() const { return gridGeometry_; }
185
186 bool verbose() const { return verbose_; }
187 const std::string& name() const { return name_; }
188 Dune::VTK::DataMode dataMode() const { return dm_; }
189 Dumux::Vtk::Precision precision() const { return precision_; }
190
191 Dune::VTKWriter<GridView>& writer() { return *writer_; }
192 Dune::VTKSequenceWriter<GridView>& sequenceWriter() { return *sequenceWriter_; }
193
194 const std::vector<Field>& fields() const { return fields_; }
195
196private:
198 virtual void writeConforming_(double time, Dune::VTK::OutputType type)
199 {
203
204 // process rank
205 std::vector<int> rank;
206
208 if (!fields_.empty() || addProcessRank_)
209 {
210 const auto numCells = gridGeometry_.gridView().size(0);
211
212 // maybe allocate space for the process rank
213 if (addProcessRank_)
214 {
215 rank.resize(numCells);
216
217 for (const auto& element : elements(gridGeometry_.gridView(), Dune::Partitions::interior))
218 {
219 const auto eIdxGlobal = gridGeometry_.elementMapper().index(element);
220 rank[eIdxGlobal] = gridGeometry_.gridView().comm().rank();
221 }
222 }
223
227
228 // the process rank
229 if (addProcessRank_)
230 this->sequenceWriter().addCellData(Field(gridGeometry_.gridView(), gridGeometry_.elementMapper(), rank, "process rank", 1, 0).get());
231
232 // also register additional (non-standardized) user fields if any
233 for (auto&& field : fields_)
234 {
235 if (field.codim() == 0)
236 this->sequenceWriter().addCellData(field.get());
237 else if (field.codim() == dim)
238 this->sequenceWriter().addVertexData(field.get());
239 else
240 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
241 }
242 }
243
247 this->sequenceWriter().write(time, type);
248
252 this->writer().clear();
253 }
254
256 virtual void writeNonConforming_(double time, Dune::VTK::OutputType type)
257 {
258 DUNE_THROW(Dune::NotImplemented, "Non-conforming VTK output");
259 }
260
262 template<class Vector>
263 std::size_t getNumberOfComponents_(const Vector& v)
264 {
265 if constexpr (Dune::IsIndexable<decltype(std::declval<Vector>()[0])>::value)
266 return v[0].size();
267 else
268 return 1;
269 }
270
271 const GridGeometry& gridGeometry_;
272 std::string name_;
273 const std::string paramGroup_;
274 Dune::VTK::DataMode dm_;
275 bool verbose_;
276 Dumux::Vtk::Precision precision_;
277
278 std::shared_ptr<Dune::VTKWriter<GridView>> writer_;
279 std::unique_ptr<Dune::VTKSequenceWriter<GridView>> sequenceWriter_;
280
281 std::vector<Field> fields_;
282
283 bool addProcessRank_ = true;
284};
285
298template<class GridVariables, class SolutionVector>
299class VtkOutputModule : public VtkOutputModuleBase<typename GridVariables::GridGeometry>
300{
302 using GridGeometry = typename GridVariables::GridGeometry;
303
304 using VV = typename GridVariables::VolumeVariables;
305 using Scalar = typename GridVariables::Scalar;
306
307 using GridView = typename GridGeometry::GridView;
308
309 enum {
310 dim = GridView::dimension,
311 dimWorld = GridView::dimensionworld
312 };
313
314 using Element = typename GridView::template Codim<0>::Entity;
315 using VolVarsVector = Dune::FieldVector<Scalar, dimWorld>;
316
317 static constexpr bool isBox = GridGeometry::discMethod == DiscretizationMethods::box;
318 static constexpr bool isDiamond = GridGeometry::discMethod == DiscretizationMethods::fcdiamond;
319 static constexpr bool isPQ1Bubble = GridGeometry::discMethod == DiscretizationMethods::pq1bubble;
320
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_; };
323
325
326public:
328 using Field = Vtk::template Field<GridView>;
330 using VolumeVariables = VV;
331
332 VtkOutputModule(const GridVariables& gridVariables,
333 const SolutionVector& sol,
334 const std::string& name,
335 const std::string& paramGroup = "",
336 Dune::VTK::DataMode dm = Dune::VTK::conforming,
337 bool verbose = true)
339 , gridVariables_(gridVariables)
340 , sol_(sol)
341 , velocityOutput_(std::make_shared<VelocityOutputType>())
342 {
343 enableVelocityOutput_ = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddVelocity", false);
344 addProcessRank_ = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddProcessRank", true);
345 }
346
351
358 void addVelocityOutput(std::shared_ptr<VelocityOutputType> velocityOutput)
359 { velocityOutput_ = velocityOutput; }
360
364 void addVolumeVariable(std::function<Scalar(const VolumeVariables&)>&& f,
365 const std::string& name)
366 {
367 volVarScalarDataInfo_.push_back(VolVarScalarDataInfo{f, name, this->precision()});
368 }
369
374 template<class VVV = VolVarsVector, typename std::enable_if_t<(VVV::dimension > 1), int> = 0>
375 void addVolumeVariable(std::function<VolVarsVector(const VolumeVariables&)>&& f,
376 const std::string& name)
377 {
378 volVarVectorDataInfo_.push_back(VolVarVectorDataInfo{f, name, this->precision()});
379 }
380
381protected:
382 // some return functions for differing implementations to use
383 const auto& problem() const { return gridVariables_.curGridVolVars().problem(); }
384 const GridVariables& gridVariables() const { return gridVariables_; }
385 const GridGeometry& gridGeometry() const { return gridVariables_.gridGeometry(); }
386 const SolutionVector& sol() const { return sol_; }
387
388 const std::vector<VolVarScalarDataInfo>& volVarScalarDataInfo() const { return volVarScalarDataInfo_; }
389 const std::vector<VolVarVectorDataInfo>& volVarVectorDataInfo() const { return volVarVectorDataInfo_; }
390
392 const VelocityOutput& velocityOutput() const { return *velocityOutput_; }
393
394private:
395
397 void writeConforming_(double time, Dune::VTK::OutputType type) override
398 {
399 const Dune::VTK::DataMode dm = Dune::VTK::conforming;
403
404 // instantiate the velocity output
405 using VelocityVector = typename VelocityOutput::VelocityVector;
406 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
407
408 // process rank
409 std::vector<double> rank;
410
411 // volume variable data
412 std::vector<std::vector<Scalar>> volVarScalarData;
413 std::vector<std::vector<VolVarsVector>> volVarVectorData;
414
416 if (!volVarScalarDataInfo_.empty()
417 || !volVarVectorDataInfo_.empty()
418 || !this->fields().empty()
419 || velocityOutput_->enableOutput()
420 || addProcessRank_)
421 {
422 const auto numCells = gridGeometry().gridView().size(0);
423 const auto numDofs = numDofs_();
424
425 // get fields for all volume variables
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));
430
431 if (velocityOutput_->enableOutput())
432 {
433 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
434 {
435 if (velocityOutput_->fieldType() == VelocityOutput::FieldType::element)
436 velocity[phaseIdx].resize(numCells);
437 else if (velocityOutput_->fieldType() == VelocityOutput::FieldType::vertex)
438 velocity[phaseIdx].resize(numDofs);
439 else
440 {
441 if(isBox && dim == 1)
442 velocity[phaseIdx].resize(numCells);
443 else
444 velocity[phaseIdx].resize(numDofs);
445 }
446 }
447 }
448
449 // maybe allocate space for the process rank
450 if (addProcessRank_) rank.resize(numCells);
451
452 auto fvGeometry = localView(gridGeometry());
453 auto elemVolVars = localView(gridVariables_.curGridVolVars());
454 for (const auto& element : elements(gridGeometry().gridView(), Dune::Partitions::interior))
455 {
456 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
457 // If velocity output is enabled we need to bind to the whole stencil
458 // otherwise element-local data is sufficient
459 if (velocityOutput_->enableOutput())
460 {
461 fvGeometry.bind(element);
462 elemVolVars.bind(element, fvGeometry, sol_);
463 }
464 else
465 {
466 fvGeometry.bindElement(element);
467 elemVolVars.bindElement(element, fvGeometry, sol_);
468 }
469
470 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
471 {
472 for (const auto& scv : scvs(fvGeometry))
473 {
474 const auto dofIdxGlobal = scv.dofIndex();
475 const auto& volVars = elemVolVars[scv];
476
477 // get the scalar-valued data
478 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
479 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo_[i].get(volVars);
480
481 // get the vector-valued data
482 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
483 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo_[i].get(volVars);
484 }
485 }
486
487 // velocity output
488 if (velocityOutput_->enableOutput())
489 {
490 const auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
491
492 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
493 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
494 }
495
497 if (addProcessRank_)
498 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
499 }
500
504
505 // volume variables if any
506 if constexpr (isBox || isPQ1Bubble)
507 {
508 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
509 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarScalarData[i],
510 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/dim, dm, this->precision()).get() );
511 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
512 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarVectorData[i],
513 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()).get() );
514
515 if constexpr (isPQ1Bubble)
516 {
517 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
518 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarScalarData[i],
519 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0,dm, this->precision()).get() );
520 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
521 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarVectorData[i],
522 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0,dm, this->precision()).get() );
523 }
524
525 }
526 else
527 {
528 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
529 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarScalarData[i],
530 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0,dm, this->precision()).get() );
531 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
532 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarVectorData[i],
533 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0,dm, this->precision()).get() );
534 }
535
536 // the velocity field
537 if (velocityOutput_->enableOutput())
538 {
539 if (velocityOutput_->fieldType() == VelocityOutput::FieldType::vertex
540 || ( (velocityOutput_->fieldType() == VelocityOutput::FieldType::automatic) && dim > 1 && isBox ))
541 {
542 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
543 this->sequenceWriter().addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
544 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
545 /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()).get() );
546 }
547 // cell-centered models
548 else
549 {
550 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
551 this->sequenceWriter().addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
552 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
553 /*numComp*/dimWorld, /*codim*/0, dm, this->precision()).get() );
554 }
555 }
556
557 // the process rank
558 if (addProcessRank_)
559 this->sequenceWriter().addCellData(Field(gridGeometry().gridView(), gridGeometry().elementMapper(), rank, "process rank", 1, 0).get());
560
561 // also register additional (non-standardized) user fields if any
562 for (auto&& field : this->fields())
563 {
564 if (field.codim() == 0)
565 this->sequenceWriter().addCellData(field.get());
566 else if (field.codim() == dim)
567 this->sequenceWriter().addVertexData(field.get());
568 else
569 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
570 }
571 }
572
576 this->sequenceWriter().write(time, type);
577
581 this->writer().clear();
582 }
583
585 void writeNonConforming_(double time, Dune::VTK::OutputType type) override
586 {
587 const Dune::VTK::DataMode dm = Dune::VTK::nonconforming;
588
589 // only supports finite-element-like discretization schemes
590 if(!isBox && !isDiamond)
591 DUNE_THROW(Dune::NotImplemented,
592 "Non-conforming output for discretization scheme " << GridGeometry::discMethod
593 );
594
598
599 // check the velocity output
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;
605 using VelocityVector = typename VelocityOutput::VelocityVector;
606 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
607
608 // process rank
609 std::vector<double> rank;
610
611 // volume variable data (indexing: volvardata/element/localindex)
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;
616
618 if (!volVarScalarDataInfo_.empty()
619 || !volVarVectorDataInfo_.empty()
620 || !this->fields().empty()
621 || velocityOutput_->enableOutput()
622 || addProcessRank_)
623 {
624 const auto numCells = gridGeometry().gridView().size(0);
625 const auto outputSize = numDofs_();
626
627 // get fields for all volume variables
628 if (!volVarScalarDataInfo_.empty())
629 volVarScalarData.resize(volVarScalarDataInfo_.size(), ScalarDataContainer(numCells));
630 if (!volVarVectorDataInfo_.empty())
631 volVarVectorData.resize(volVarVectorDataInfo_.size(), VectorDataContainer(numCells));
632
633 if (velocityOutput_->enableOutput())
634 {
635 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
636 {
637 if((isBox && dim == 1) || isDiamond)
638 velocity[phaseIdx].resize(numCells);
639 else
640 velocity[phaseIdx].resize(outputSize);
641 }
642 }
643
644 // maybe allocate space for the process rank
645 if (addProcessRank_) rank.resize(numCells);
646
647 // now we go element-local to extract values at local dof locations
648 auto fvGeometry = localView(gridGeometry());
649 auto elemVolVars = localView(gridVariables_.curGridVolVars());
650 for (const auto& element : elements(gridGeometry().gridView(), Dune::Partitions::interior))
651 {
652 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
653 // If velocity output is enabled we need to bind to the whole stencil
654 // otherwise element-local data is sufficient
655 if (velocityOutput_->enableOutput())
656 {
657 fvGeometry.bind(element);
658 elemVolVars.bind(element, fvGeometry, sol_);
659 }
660 else
661 {
662 fvGeometry.bindElement(element);
663 elemVolVars.bindElement(element, fvGeometry, sol_);
664 }
665
666 const auto numLocalDofs = fvGeometry.numScv();
667 // resize element-local data containers
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);
672
673 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
674 {
675 for (const auto& scv : scvs(fvGeometry))
676 {
677 const auto& volVars = elemVolVars[scv];
678
679 // get the scalar-valued data
680 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
681 volVarScalarData[i][eIdxGlobal][scv.localDofIndex()] = volVarScalarDataInfo_[i].get(volVars);
682
683 // get the vector-valued data
684 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
685 volVarVectorData[i][eIdxGlobal][scv.localDofIndex()] = volVarVectorDataInfo_[i].get(volVars);
686 }
687 }
688
689 // velocity output
690 if (velocityOutput_->enableOutput())
691 {
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);
695 }
696
698 if (addProcessRank_)
699 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
700 }
701
705
706 // volume variables if any
707 static constexpr int dofLocCodim = isDiamond ? 1 : dim;
708 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
709 this->sequenceWriter().addVertexData(Field(
710 gridGeometry().gridView(), gridGeometry().elementMapper(),
711 volVarScalarData[i], volVarScalarDataInfo_[i].name,
712 /*numComp*/1, /*codim*/dofLocCodim, /*nonconforming*/dm, this->precision()
713 ).get());
714
715 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
716 this->sequenceWriter().addVertexData(Field(
717 gridGeometry().gridView(), gridGeometry().elementMapper(),
718 volVarVectorData[i], volVarVectorDataInfo_[i].name,
719 /*numComp*/dimWorld, /*codim*/dofLocCodim, /*nonconforming*/dm, this->precision()
720 ).get());
721
722 // the velocity field
723 if (velocityOutput_->enableOutput())
724 {
725 // node-wise velocities
726 if (dim > 1 && !isDiamond)
727 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
728 this->sequenceWriter().addVertexData(Field(
729 gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
730 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
731 /*numComp*/dimWorld, /*codim*/dofLocCodim, dm, this->precision()
732 ).get());
733
734 // cell-wise velocities
735 else
736 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
737 this->sequenceWriter().addCellData(Field(
738 gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
739 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
740 /*numComp*/dimWorld, /*codim*/0, dm, this->precision()
741 ).get());
742 }
743 }
744
745 // the process rank
746 if (addProcessRank_)
747 this->sequenceWriter().addCellData(Field(
748 gridGeometry().gridView(), gridGeometry().elementMapper(),
749 rank, "process rank", /*numComp*/1, /*codim*/0
750 ).get());
751
752 // also register additional (non-standardized) user fields if any
753 for (const auto& field : this->fields())
754 {
755 if (field.codim() == 0)
756 this->sequenceWriter().addCellData(field.get());
757 else if (field.codim() == dim || field.codim() == 1)
758 this->sequenceWriter().addVertexData(field.get());
759 else
760 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
761 }
762
766 this->sequenceWriter().write(time, type);
767
771 this->writer().clear();
772 }
773
775 std::size_t numDofs_() const
776 {
777 // TODO this should actually always be dofMapper.size()
778 // maybe some discretizations needs special treatment (?)
779 if constexpr (isBox || isDiamond || isPQ1Bubble)
780 return gridGeometry().dofMapper().size();
781 else
782 return gridGeometry().elementMapper().size();
783 }
784
785 const GridVariables& gridVariables_;
786 const SolutionVector& sol_;
787
788 std::vector<VolVarScalarDataInfo> volVarScalarDataInfo_;
789 std::vector<VolVarVectorDataInfo> volVarVectorDataInfo_;
790
791 std::shared_ptr<VelocityOutput> velocityOutput_;
792 bool enableVelocityOutput_ = false;
793 bool addProcessRank_ = true;
794};
795
796} // end namespace Dumux
797
798#endif
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 &paramGroup="", 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 &paramGroup="", 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.
Formatting based on the fmt-library which implements std::format of C++20.
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
Precision stringToPrecision(std::string_view precisionName)
Maps a string (e.g. from input) to a Dune precision type.
Definition: precision.hh:27
FieldType
Identifier for vtk field types.
Definition: fieldtype.hh:22
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
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
TODO: docme!