version 3.11-dev
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-FileCopyrightText: 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 addField(Field(gridGeometry_.gridView(), gridGeometry_.elementMapper(), v, name, nComp, 0, dm_, precision));
147 else
148 addField(Field(gridGeometry_.gridView(), gridGeometry_.vertexMapper(), v, name, nComp, dim, dm_, precision));
149 }
150
155 void addField(Field&& field)
156 {
157 // data arrays in the vtk output have to be unique within cell and point data
158 // look if we have a field by the same name with the same codim, if so, replace it
159 for (auto i = 0UL; i < fields_.size(); ++i)
160 {
161 if (fields_[i].name() == field.name() && fields_[i].codim() == field.codim())
162 {
163 fields_[i] = std::move(field);
164 std::cout << Fmt::format(
165 "VtkOutputModule: Replaced field \"{}\" (codim {}). "
166 "A field by the same name & codim had already been registered previously.\n",
167 field.name(), field.codim()
168 );
169 return;
170 }
171 }
172
173 // otherwise add it to the end of the fields
174 fields_.push_back(std::move(field));
175 }
176
182 void write(double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
183 {
184 Dune::Timer timer;
185
186 // write to file depending on data mode
187 if (dm_ == Dune::VTK::conforming)
188 writeConforming_(time, type);
189 else if (dm_ == Dune::VTK::nonconforming)
190 writeNonConforming_(time, type);
191 else
192 DUNE_THROW(Dune::NotImplemented, "Output for provided VTK data mode");
193
195 timer.stop();
196 if (verbose_)
197 std::cout << Fmt::format("Writing output for problem \"{}\". Took {:.2g} seconds.\n", name_, timer.elapsed());
198 }
199
200protected:
201 const GridGeometry& gridGeometry() const { return gridGeometry_; }
202
203 bool verbose() const { return verbose_; }
204 const std::string& name() const { return name_; }
205 Dune::VTK::DataMode dataMode() const { return dm_; }
206 Dumux::Vtk::Precision precision() const { return precision_; }
207
208 Dune::VTKWriter<GridView>& writer() { return *writer_; }
209 Dune::VTKSequenceWriter<GridView>& sequenceWriter() { return *sequenceWriter_; }
210
211 const std::vector<Field>& fields() const { return fields_; }
212
213 void addCellData(const Field& field)
214 {
215 if (std::ranges::find(addedCellData_, field.name()) == addedCellData_.end())
216 {
217 this->sequenceWriter().addCellData(field.get());
218 addedCellData_.push_back(field.name());
219 }
220 }
221
222 void addVertexData(const Field& field)
223 {
224 if (std::ranges::find(addedVertexData_, field.name()) == addedVertexData_.end())
225 {
226 this->sequenceWriter().addVertexData(field.get());
227 addedVertexData_.push_back(field.name());
228 }
229 }
230
231 // keep track of what has been already added because Dune::VTK::Writer
232 // potentially adds the same field twice which is not allowed in VTK/Paraview
233 std::vector<std::string> addedCellData_;
234 std::vector<std::string> addedVertexData_;
235
236private:
238 virtual void writeConforming_(double time, Dune::VTK::OutputType type)
239 {
243
244 // process rank
245 std::vector<int> rank;
246
248 if (!fields_.empty() || addProcessRank_)
249 {
250 const auto numCells = gridGeometry_.gridView().size(0);
251
252 // maybe allocate space for the process rank
253 if (addProcessRank_)
254 {
255 rank.resize(numCells);
256
257 for (const auto& element : elements(gridGeometry_.gridView(), Dune::Partitions::interior))
258 {
259 const auto eIdxGlobal = gridGeometry_.elementMapper().index(element);
260 rank[eIdxGlobal] = gridGeometry_.gridView().comm().rank();
261 }
262 }
263
267
268 // the process rank
269 if (addProcessRank_)
270 this->addCellData(Field(gridGeometry_.gridView(), gridGeometry_.elementMapper(), rank, "process rank", 1, 0));
271
272 // also register additional (non-standardized) user fields if any
273 for (auto&& field : fields_)
274 {
275 if (field.codim() == 0)
276 this->addCellData(field);
277 else if (field.codim() == dim)
278 this->addVertexData(field);
279 else
280 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
281 }
282 }
283
287 this->sequenceWriter().write(time, type);
288
292 this->writer().clear();
293
294 this->addedCellData_.clear();
295 this->addedVertexData_.clear();
296 }
297
299 virtual void writeNonConforming_(double time, Dune::VTK::OutputType type)
300 {
301 DUNE_THROW(Dune::NotImplemented, "Non-conforming VTK output");
302 }
303
305 template<class Vector>
306 std::size_t getNumberOfComponents_(const Vector& v)
307 {
308 if constexpr (Dune::IsIndexable<decltype(std::declval<Vector>()[0])>::value)
309 return v[0].size();
310 else
311 return 1;
312 }
313
314 const GridGeometry& gridGeometry_;
315 std::string name_;
316 const std::string paramGroup_;
317 Dune::VTK::DataMode dm_;
318 bool verbose_;
319 Dumux::Vtk::Precision precision_;
320
321 std::shared_ptr<Dune::VTKWriter<GridView>> writer_;
322 std::unique_ptr<Dune::VTKSequenceWriter<GridView>> sequenceWriter_;
323
324 std::vector<Field> fields_;
325
326 bool addProcessRank_ = true;
327};
328
341template<class GridVariables, class SolutionVector>
342class VtkOutputModule : public VtkOutputModuleBase<typename GridVariables::GridGeometry>
343{
345 using GridGeometry = typename GridVariables::GridGeometry;
346
347 using VV = typename GridVariables::VolumeVariables;
348 using Scalar = typename GridVariables::Scalar;
349
350 using GridView = typename GridGeometry::GridView;
351
352 enum {
353 dim = GridView::dimension,
354 dimWorld = GridView::dimensionworld
355 };
356
357 using Element = typename GridView::template Codim<0>::Entity;
358 using VolVarsVector = Dune::FieldVector<Scalar, dimWorld>;
359
360 static constexpr bool isBox = GridGeometry::discMethod == DiscretizationMethods::box;
361 static constexpr bool isDiamond = GridGeometry::discMethod == DiscretizationMethods::fcdiamond;
362 static constexpr bool isPQ1Bubble = GridGeometry::discMethod == DiscretizationMethods::pq1bubble;
363 static constexpr bool isPQ2 = GridGeometry::discMethod == DiscretizationMethods::pq2;
364
365 struct VolVarScalarDataInfo { std::function<Scalar(const VV&)> get; std::string name; Dumux::Vtk::Precision precision_; };
366 struct VolVarVectorDataInfo { std::function<VolVarsVector(const VV&)> get; std::string name; Dumux::Vtk::Precision precision_; };
367
369
370public:
372 using Field = Vtk::template Field<GridView>;
374 using VolumeVariables = VV;
375
376 VtkOutputModule(const GridVariables& gridVariables,
377 const SolutionVector& sol,
378 const std::string& name,
379 const std::string& paramGroup = "",
380 Dune::VTK::DataMode dm = Dune::VTK::conforming,
381 bool verbose = true)
383 , gridVariables_(gridVariables)
384 , sol_(sol)
385 , velocityOutput_(std::make_shared<VelocityOutputType>())
386 {
387 enableVelocityOutput_ = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddVelocity", false);
388 addProcessRank_ = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddProcessRank", true);
389 }
390
395
402 void addVelocityOutput(std::shared_ptr<VelocityOutputType> velocityOutput)
403 { velocityOutput_ = velocityOutput; }
404
408 void addVolumeVariable(std::function<Scalar(const VolumeVariables&)>&& f,
409 const std::string& name)
410 {
411 // data arrays in the vtk output have to be unique within cell and point data
412 // look if we have a field by the same name with the same codim, if so, replace it
413 for (auto i = 0UL; i < volVarScalarDataInfo_.size(); ++i)
414 {
415 if (volVarScalarDataInfo_[i].name == name)
416 {
417 volVarScalarDataInfo_[i] = VolVarScalarDataInfo{f, name, this->precision()};
418 std::cout << Fmt::format(
419 "VtkOutputModule: Replaced volume variable output \"{}\". "
420 "A field by the same name had already been registered previously.\n",
421 name
422 );
423 return;
424 }
425 }
426
427 // otherwise add it to the end of the fields
428 volVarScalarDataInfo_.push_back(VolVarScalarDataInfo{f, name, this->precision()});
429 }
430
435 template<class VVV = VolVarsVector, typename std::enable_if_t<(VVV::dimension > 1), int> = 0>
436 void addVolumeVariable(std::function<VolVarsVector(const VolumeVariables&)>&& f,
437 const std::string& name)
438 {
439 // data arrays in the vtk output have to be unique within cell and point data
440 // look if we have a field by the same name with the same codim, if so, replace it
441 for (auto i = 0UL; i < volVarVectorDataInfo_.size(); ++i)
442 {
443 if (volVarVectorDataInfo_[i].name == name)
444 {
445 volVarVectorDataInfo_[i] = VolVarVectorDataInfo{f, name, this->precision()};
446 std::cout << Fmt::format(
447 "VtkOutputModule: Replaced volume variable output \"{}\". "
448 "A field by the same name had already been registered previously.\n",
449 name
450 );
451 return;
452 }
453 }
454
455 // otherwise add it to the end of the fields
456 volVarVectorDataInfo_.push_back(VolVarVectorDataInfo{f, name, this->precision()});
457 }
458
459protected:
460 // some return functions for differing implementations to use
461 const auto& problem() const { return gridVariables_.curGridVolVars().problem(); }
462 const GridVariables& gridVariables() const { return gridVariables_; }
463 const GridGeometry& gridGeometry() const { return gridVariables_.gridGeometry(); }
464 const SolutionVector& sol() const { return sol_; }
465
466 const std::vector<VolVarScalarDataInfo>& volVarScalarDataInfo() const { return volVarScalarDataInfo_; }
467 const std::vector<VolVarVectorDataInfo>& volVarVectorDataInfo() const { return volVarVectorDataInfo_; }
468
470 const VelocityOutput& velocityOutput() const { return *velocityOutput_; }
471
472private:
473
475 void writeConforming_(double time, Dune::VTK::OutputType type) override
476 {
477 const Dune::VTK::DataMode dm = Dune::VTK::conforming;
481
482 // instantiate the velocity output
483 using VelocityVector = typename VelocityOutput::VelocityVector;
484 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
485
486 // process rank
487 std::vector<double> rank;
488
489 // volume variable data
490 std::vector<std::vector<Scalar>> volVarScalarData;
491 std::vector<std::vector<VolVarsVector>> volVarVectorData;
492
494 if (!volVarScalarDataInfo_.empty()
495 || !volVarVectorDataInfo_.empty()
496 || !this->fields().empty()
497 || velocityOutput_->enableOutput()
498 || addProcessRank_)
499 {
500 const auto numCells = gridGeometry().gridView().size(0);
501 const auto numDofs = numDofs_();
502
503 // get fields for all volume variables
504 if (!volVarScalarDataInfo_.empty())
505 volVarScalarData.resize(volVarScalarDataInfo_.size(), std::vector<Scalar>(numDofs));
506 if (!volVarVectorDataInfo_.empty())
507 volVarVectorData.resize(volVarVectorDataInfo_.size(), std::vector<VolVarsVector>(numDofs));
508
509 if (velocityOutput_->enableOutput())
510 {
511 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
512 {
513 if (velocityOutput_->fieldType() == VelocityOutput::FieldType::element)
514 velocity[phaseIdx].resize(numCells);
515 else if (velocityOutput_->fieldType() == VelocityOutput::FieldType::vertex)
516 velocity[phaseIdx].resize(numDofs);
517 else
518 {
519 if(isBox && dim == 1)
520 velocity[phaseIdx].resize(numCells);
521 else
522 velocity[phaseIdx].resize(numDofs);
523 }
524 }
525 }
526
527 // maybe allocate space for the process rank
528 if (addProcessRank_) rank.resize(numCells);
529
530 auto fvGeometry = localView(gridGeometry());
531 auto elemVolVars = localView(gridVariables_.curGridVolVars());
532 for (const auto& element : elements(gridGeometry().gridView()))
533 {
534 if (!velocityOutput_->enableOutput() &&
535 element.partitionType() != Dune::PartitionType::InteriorEntity)
536 {
537 continue;
538 }
539 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
540 // If velocity output is enabled we need to bind to the whole stencil
541 // otherwise element-local data is sufficient
542 if (velocityOutput_->enableOutput())
543 {
544 fvGeometry.bind(element);
545 elemVolVars.bind(element, fvGeometry, sol_);
546 }
547 else
548 {
549 fvGeometry.bindElement(element);
550 elemVolVars.bindElement(element, fvGeometry, sol_);
551 }
552
553 // velocity output
554 if (velocityOutput_->enableOutput())
555 {
556 const auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
557
558 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
559 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
560 }
561 else if (element.partitionType() != Dune::PartitionType::InteriorEntity)
562 {
563 continue;
564 }
565
566 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
567 {
568 for (const auto& scv : scvs(fvGeometry))
569 {
570 const auto dofIdxGlobal = scv.dofIndex();
571 const auto& volVars = elemVolVars[scv];
572
573 // get the scalar-valued data
574 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
575 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo_[i].get(volVars);
576
577 // get the vector-valued data
578 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
579 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo_[i].get(volVars);
580 }
581 }
582
584 if (addProcessRank_)
585 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
586 }
587
591
592 // volume variables if any
593 if constexpr (isBox || isPQ1Bubble || isPQ2)
594 {
595 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
596 this->addVertexData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarScalarData[i],
597 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/dim, dm, this->precision()) );
598 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
599 this->addVertexData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarVectorData[i],
600 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()) );
601
602 if constexpr (isPQ1Bubble)
603 {
604 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
605 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarScalarData[i],
606 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0,dm, this->precision()) );
607 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
608 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarVectorData[i],
609 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0,dm, this->precision()) );
610 }
611
612 }
613 else
614 {
615 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
616 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarScalarData[i],
617 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0,dm, this->precision()) );
618 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
619 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarVectorData[i],
620 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0,dm, this->precision()) );
621 }
622
623 // the velocity field
624 if (velocityOutput_->enableOutput())
625 {
626 if (velocityOutput_->fieldType() == VelocityOutput::FieldType::vertex
627 || ( (velocityOutput_->fieldType() == VelocityOutput::FieldType::automatic) && dim > 1 && isBox ))
628 {
629 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
630 this->addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
631 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
632 /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()) );
633 }
634 // cell-centered models
635 else
636 {
637 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
638 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
639 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
640 /*numComp*/dimWorld, /*codim*/0, dm, this->precision()) );
641 }
642 }
643
644 // the process rank
645 if (addProcessRank_)
646 this->addCellData(Field(gridGeometry().gridView(), gridGeometry().elementMapper(), rank, "process rank", 1, 0));
647
648 // also register additional (non-standardized) user fields if any
649 for (auto&& field : this->fields())
650 {
651 if (field.codim() == 0)
652 this->addCellData(field);
653 else if (field.codim() == dim)
654 this->addVertexData(field);
655 else
656 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
657 }
658 }
659
663 this->sequenceWriter().write(time, type);
664
668 this->writer().clear();
669
670 this->addedCellData_.clear();
671 this->addedVertexData_.clear();
672 }
673
675 void writeNonConforming_(double time, Dune::VTK::OutputType type) override
676 {
677 const Dune::VTK::DataMode dm = Dune::VTK::nonconforming;
678
679 // only supports finite-element-like discretization schemes
680 if(!isBox && !isDiamond)
681 DUNE_THROW(Dune::NotImplemented,
682 "Non-conforming output for discretization scheme " << GridGeometry::discMethod
683 );
684
688
689 // check the velocity output
690 if (enableVelocityOutput_ && !velocityOutput_->enableOutput())
691 std::cerr << "Warning! Velocity output was enabled in the input file"
692 << " but no velocity output policy was set for the VTK output module:"
693 << " There will be no velocity output."
694 << " Use the addVelocityOutput member function of the VTK output module." << std::endl;
695 using VelocityVector = typename VelocityOutput::VelocityVector;
696 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
697
698 // process rank
699 std::vector<double> rank;
700
701 // volume variable data (indexing: volvardata/element/localindex)
702 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
703 using VectorDataContainer = std::vector< std::vector<VolVarsVector> >;
704 std::vector< ScalarDataContainer > volVarScalarData;
705 std::vector< VectorDataContainer > volVarVectorData;
706
708 if (!volVarScalarDataInfo_.empty()
709 || !volVarVectorDataInfo_.empty()
710 || !this->fields().empty()
711 || velocityOutput_->enableOutput()
712 || addProcessRank_)
713 {
714 const auto numCells = gridGeometry().gridView().size(0);
715 const auto outputSize = numDofs_();
716
717 // get fields for all volume variables
718 if (!volVarScalarDataInfo_.empty())
719 volVarScalarData.resize(volVarScalarDataInfo_.size(), ScalarDataContainer(numCells));
720 if (!volVarVectorDataInfo_.empty())
721 volVarVectorData.resize(volVarVectorDataInfo_.size(), VectorDataContainer(numCells));
722
723 if (velocityOutput_->enableOutput())
724 {
725 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
726 {
727 if((isBox && dim == 1) || isDiamond)
728 velocity[phaseIdx].resize(numCells);
729 else
730 velocity[phaseIdx].resize(outputSize);
731 }
732 }
733
734 // maybe allocate space for the process rank
735 if (addProcessRank_) rank.resize(numCells);
736
737 // now we go element-local to extract values at local dof locations
738 auto fvGeometry = localView(gridGeometry());
739 auto elemVolVars = localView(gridVariables_.curGridVolVars());
740 for (const auto& element : elements(gridGeometry().gridView()))
741 {
742 if (!velocityOutput_->enableOutput() &&
743 element.partitionType() != Dune::PartitionType::InteriorEntity)
744 {
745 continue;
746 }
747 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
748 // If velocity output is enabled we need to bind to the whole stencil
749 // otherwise element-local data is sufficient
750 if (velocityOutput_->enableOutput())
751 {
752 fvGeometry.bind(element);
753 elemVolVars.bind(element, fvGeometry, sol_);
754 }
755 else
756 {
757 fvGeometry.bindElement(element);
758 elemVolVars.bindElement(element, fvGeometry, sol_);
759 }
760
761 const auto numLocalDofs = fvGeometry.numScv();
762 // resize element-local data containers
763 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
764 volVarScalarData[i][eIdxGlobal].resize(numLocalDofs);
765 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
766 volVarVectorData[i][eIdxGlobal].resize(numLocalDofs);
767
768 // velocity output
769 if (velocityOutput_->enableOutput())
770 {
771 const auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
772 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
773 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
774 }
775 else if (element.partitionType() != Dune::PartitionType::InteriorEntity)
776 {
777 continue;
778 }
779
780 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
781 {
782 for (const auto& scv : scvs(fvGeometry))
783 {
784 const auto& volVars = elemVolVars[scv];
785
786 // get the scalar-valued data
787 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
788 volVarScalarData[i][eIdxGlobal][scv.localDofIndex()] = volVarScalarDataInfo_[i].get(volVars);
789
790 // get the vector-valued data
791 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
792 volVarVectorData[i][eIdxGlobal][scv.localDofIndex()] = volVarVectorDataInfo_[i].get(volVars);
793 }
794 }
795
797 if (addProcessRank_)
798 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
799 }
800
804
805 // volume variables if any
806 static constexpr int dofLocCodim = isDiamond ? 1 : dim;
807 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
808 this->addVertexData(Field(
809 gridGeometry().gridView(), gridGeometry().elementMapper(),
810 volVarScalarData[i], volVarScalarDataInfo_[i].name,
811 /*numComp*/1, /*codim*/dofLocCodim, /*nonconforming*/dm, this->precision()
812 ));
813
814 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
815 this->addVertexData(Field(
816 gridGeometry().gridView(), gridGeometry().elementMapper(),
817 volVarVectorData[i], volVarVectorDataInfo_[i].name,
818 /*numComp*/dimWorld, /*codim*/dofLocCodim, /*nonconforming*/dm, this->precision()
819 ));
820
821 // the velocity field
822 if (velocityOutput_->enableOutput())
823 {
824 // node-wise velocities
825 if (dim > 1 && !isDiamond)
826 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
827 this->addVertexData(Field(
828 gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
829 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
830 /*numComp*/dimWorld, /*codim*/dofLocCodim, dm, this->precision()
831 ));
832
833 // cell-wise velocities
834 else
835 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
836 this->addCellData(Field(
837 gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
838 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
839 /*numComp*/dimWorld, /*codim*/0, dm, this->precision()
840 ));
841 }
842 }
843
844 // the process rank
845 if (addProcessRank_)
846 this->addCellData(Field(
847 gridGeometry().gridView(), gridGeometry().elementMapper(),
848 rank, "process rank", /*numComp*/1, /*codim*/0
849 ));
850
851 // also register additional (non-standardized) user fields if any
852 for (const auto& field : this->fields())
853 {
854 if (field.codim() == 0)
855 this->addCellData(field);
856 else if (field.codim() == dim || field.codim() == 1)
857 this->addVertexData(field);
858 else
859 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
860 }
861
865 this->sequenceWriter().write(time, type);
866
870 this->writer().clear();
871
872 this->addedCellData_.clear();
873 this->addedVertexData_.clear();
874 }
875
877 std::size_t numDofs_() const
878 {
879 // TODO this should actually always be dofMapper.size()
880 // maybe some discretizations needs special treatment (?)
881 if constexpr (isBox || isDiamond || isPQ1Bubble || isPQ2)
882 return gridGeometry().dofMapper().size();
883 else
884 return gridGeometry().elementMapper().size();
885 }
886
887 const GridVariables& gridVariables_;
888 const SolutionVector& sol_;
889
890 std::vector<VolVarScalarDataInfo> volVarScalarDataInfo_;
891 std::vector<VolVarVectorDataInfo> volVarVectorDataInfo_;
892
893 std::shared_ptr<VelocityOutput> velocityOutput_;
894 bool enableVelocityOutput_ = false;
895 bool addProcessRank_ = true;
896};
897
898} // end namespace Dumux
899
900#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:206
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:79
void addCellData(const Field &field)
Definition: io/vtkoutputmodule.hh:213
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Definition: io/vtkoutputmodule.hh:182
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:205
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:208
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:209
const std::string & name() const
Definition: io/vtkoutputmodule.hh:204
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:211
virtual ~VtkOutputModuleBase()=default
virtual void writeNonConforming_(double time, Dune::VTK::OutputType type)
Assembles the fields and adds them to the writer (nonconforming output)
Definition: io/vtkoutputmodule.hh:299
std::vector< std::string > addedCellData_
Definition: io/vtkoutputmodule.hh:233
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:203
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:55
std::vector< std::string > addedVertexData_
Definition: io/vtkoutputmodule.hh:234
virtual void writeConforming_(double time, Dune::VTK::OutputType type)
Assembles the fields and adds them to the writer (conforming output)
Definition: io/vtkoutputmodule.hh:238
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
void addVertexData(const Field &field)
Definition: io/vtkoutputmodule.hh:222
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:201
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:343
const auto & problem() const
Definition: io/vtkoutputmodule.hh:461
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:462
void addVolumeVariable(std::function< VolVarsVector(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:436
VV VolumeVariables
export type of the volume variables for the outputfields
Definition: io/vtkoutputmodule.hh:374
void addVolumeVariable(std::function< Scalar(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:408
void addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:402
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:470
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:372
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:464
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:463
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:376
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:466
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:467
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:162
constexpr PQ2 pq2
Definition: method.hh:157
constexpr Box box
Definition: method.hh:156
constexpr PQ1Bubble pq1bubble
Definition: method.hh:158
Definition: adapt.hh:17
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
TODO: docme!