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
32#include <dumux/common/concepts/variables_.hh>
34#include <dumux/io/format.hh>
36
39#include "velocityoutput.hh"
40
41namespace Dumux {
42
48template<class GridGeometry>
50{
51 using GridView = typename GridGeometry::GridView;
52 static constexpr int dim = GridView::dimension;
53
54public:
56 using Field = Vtk::template Field<GridView>;
57
58 VtkOutputModuleBase(const GridGeometry& gridGeometry,
59 const std::string& name,
60 const std::string& paramGroup = "",
61 Dune::VTK::DataMode dm = Dune::VTK::conforming,
62 bool verbose = true)
63 : gridGeometry_(gridGeometry)
64 , name_(name)
65 , paramGroup_(paramGroup)
66 , dm_(dm)
67 , verbose_(gridGeometry.gridView().comm().rank() == 0 && verbose)
68 {
69 const auto precisionString = getParamFromGroup<std::string>(paramGroup, "Vtk.Precision", "Float32");
70 precision_ = Dumux::Vtk::stringToPrecision(precisionString);
71 const auto coordPrecision = Dumux::Vtk::stringToPrecision(getParamFromGroup<std::string>(paramGroup, "Vtk.CoordPrecision", precisionString));
72 writer_ = std::make_shared<Dune::VTKWriter<GridView>>(gridGeometry.gridView(), dm, coordPrecision);
73 sequenceWriter_ = std::make_unique<Dune::VTKSequenceWriter<GridView>>(writer_, name);
74 addProcessRank_ = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddProcessRank", true);
75 }
76
77 virtual ~VtkOutputModuleBase() = default;
78
80 const std::string& paramGroup() const
81 { return paramGroup_; }
82
92 template<typename Vector>
93 void addField(const Vector& v,
94 const std::string& name,
96 { addField(v, name, this->precision(), fieldType); }
97
108 template<typename Vector>
109 void addField(const Vector& v,
110 const std::string& name,
111 Dumux::Vtk::Precision precision,
113 {
114 // Deduce the number of components from the given vector type
115 const auto nComp = getNumberOfComponents_(v);
116
117 const auto numElemDofs = gridGeometry().elementMapper().size();
118 const auto numVertexDofs = gridGeometry().vertexMapper().size();
119
120 // Automatically deduce the field type ...
121 if(fieldType == Vtk::FieldType::automatic)
122 {
123 if(numElemDofs == numVertexDofs)
124 DUNE_THROW(Dune::InvalidStateException, "Automatic deduction of FieldType failed. Please explicitly specify FieldType::element or FieldType::vertex.");
125
126 if(v.size() == numElemDofs)
127 fieldType = Vtk::FieldType::element;
128 else if(v.size() == numVertexDofs)
129 fieldType = Vtk::FieldType::vertex;
130 else
131 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
132 }
133 // ... or check if the user-specified type matches the size of v
134 else
135 {
136 if(fieldType == Vtk::FieldType::element)
137 if(v.size() != numElemDofs)
138 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
139
140 if(fieldType == Vtk::FieldType::vertex)
141 if(v.size() != numVertexDofs)
142 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
143 }
144
145 // add the appropriate field
146 if (fieldType == Vtk::FieldType::element)
147 addField(Field(gridGeometry_.gridView(), gridGeometry_.elementMapper(), v, name, nComp, 0, dm_, precision));
148 else
149 addField(Field(gridGeometry_.gridView(), gridGeometry_.vertexMapper(), v, name, nComp, dim, dm_, precision));
150 }
151
156 void addField(Field&& field)
157 {
158 // data arrays in the vtk output have to be unique within cell and point data
159 // look if we have a field by the same name with the same codim, if so, replace it
160 for (auto i = 0UL; i < fields_.size(); ++i)
161 {
162 if (fields_[i].name() == field.name() && fields_[i].codim() == field.codim())
163 {
164 fields_[i] = std::move(field);
165 std::cout << Fmt::format(
166 "VtkOutputModule: Replaced field \"{}\" (codim {}). "
167 "A field by the same name & codim had already been registered previously.\n",
168 field.name(), field.codim()
169 );
170 return;
171 }
172 }
173
174 // otherwise add it to the end of the fields
175 fields_.push_back(std::move(field));
176 }
177
183 void write(double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
184 {
185 Dune::Timer timer;
186
187 // write to file depending on data mode
188 if (dm_ == Dune::VTK::conforming)
189 writeConforming_(time, type);
190 else if (dm_ == Dune::VTK::nonconforming)
191 writeNonConforming_(time, type);
192 else
193 DUNE_THROW(Dune::NotImplemented, "Output for provided VTK data mode");
194
196 timer.stop();
197 if (verbose_)
198 std::cout << Fmt::format("Writing output for problem \"{}\". Took {:.2g} seconds.\n", name_, timer.elapsed());
199 }
200
201protected:
202 const GridGeometry& gridGeometry() const { return gridGeometry_; }
203
204 bool verbose() const { return verbose_; }
205 const std::string& name() const { return name_; }
206 Dune::VTK::DataMode dataMode() const { return dm_; }
207 Dumux::Vtk::Precision precision() const { return precision_; }
208
209 Dune::VTKWriter<GridView>& writer() { return *writer_; }
210 Dune::VTKSequenceWriter<GridView>& sequenceWriter() { return *sequenceWriter_; }
211
212 const std::vector<Field>& fields() const { return fields_; }
213
214 void addCellData(const Field& field)
215 {
216 if (std::ranges::find(addedCellData_, field.name()) == addedCellData_.end())
217 {
218 this->sequenceWriter().addCellData(field.get());
219 addedCellData_.push_back(field.name());
220 }
221 }
222
223 void addVertexData(const Field& field)
224 {
225 if (std::ranges::find(addedVertexData_, field.name()) == addedVertexData_.end())
226 {
227 this->sequenceWriter().addVertexData(field.get());
228 addedVertexData_.push_back(field.name());
229 }
230 }
231
232 // keep track of what has been already added because Dune::VTK::Writer
233 // potentially adds the same field twice which is not allowed in VTK/Paraview
234 std::vector<std::string> addedCellData_;
235 std::vector<std::string> addedVertexData_;
236
237private:
239 virtual void writeConforming_(double time, Dune::VTK::OutputType type)
240 {
244
245 // process rank
246 std::vector<int> rank;
247
249 if (!fields_.empty() || addProcessRank_)
250 {
251 const auto numCells = gridGeometry_.gridView().size(0);
252
253 // maybe allocate space for the process rank
254 if (addProcessRank_)
255 {
256 rank.resize(numCells);
257
258 for (const auto& element : elements(gridGeometry_.gridView(), Dune::Partitions::interior))
259 {
260 const auto eIdxGlobal = gridGeometry_.elementMapper().index(element);
261 rank[eIdxGlobal] = gridGeometry_.gridView().comm().rank();
262 }
263 }
264
268
269 // the process rank
270 if (addProcessRank_)
271 this->addCellData(Field(gridGeometry_.gridView(), gridGeometry_.elementMapper(), rank, "process rank", 1, 0));
272
273 // also register additional (non-standardized) user fields if any
274 for (auto&& field : fields_)
275 {
276 if (field.codim() == 0)
277 this->addCellData(field);
278 else if (field.codim() == dim)
279 this->addVertexData(field);
280 else
281 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
282 }
283 }
284
288 this->sequenceWriter().write(time, type);
289
293 this->writer().clear();
294
295 this->addedCellData_.clear();
296 this->addedVertexData_.clear();
297 }
298
300 virtual void writeNonConforming_(double time, Dune::VTK::OutputType type)
301 {
302 DUNE_THROW(Dune::NotImplemented, "Non-conforming VTK output");
303 }
304
306 template<class Vector>
307 std::size_t getNumberOfComponents_(const Vector& v)
308 {
309 if constexpr (Dune::IsIndexable<decltype(std::declval<Vector>()[0])>::value)
310 return v[0].size();
311 else
312 return 1;
313 }
314
315 const GridGeometry& gridGeometry_;
316 std::string name_;
317 const std::string paramGroup_;
318 Dune::VTK::DataMode dm_;
319 bool verbose_;
320 Dumux::Vtk::Precision precision_;
321
322 std::shared_ptr<Dune::VTKWriter<GridView>> writer_;
323 std::unique_ptr<Dune::VTKSequenceWriter<GridView>> sequenceWriter_;
324
325 std::vector<Field> fields_;
326
327 bool addProcessRank_ = true;
328};
329
342template<class GridVariables, class SolutionVector>
343class VtkOutputModule : public VtkOutputModuleBase<typename GridVariables::GridGeometry>
344{
346 using GridGeometry = typename GridVariables::GridGeometry;
347
348 using VV = Concept::Variables_t<GridVariables>;
349 using Scalar = typename GridVariables::Scalar;
350
351 using GridView = typename GridGeometry::GridView;
352
353 enum {
354 dim = GridView::dimension,
355 dimWorld = GridView::dimensionworld
356 };
357
358 using Element = typename GridView::template Codim<0>::Entity;
359 using VolVarsVector = Dune::FieldVector<Scalar, dimWorld>;
360
361 static constexpr bool isBox = GridGeometry::discMethod == DiscretizationMethods::box;
362 static constexpr bool isDiamond = GridGeometry::discMethod == DiscretizationMethods::fcdiamond;
363 static constexpr bool isPQ1Bubble = GridGeometry::discMethod == DiscretizationMethods::pq1bubble;
364 static constexpr bool isPQ2 = GridGeometry::discMethod == DiscretizationMethods::pq2;
365
366 struct VolVarScalarDataInfo { std::function<Scalar(const VV&)> get; std::string name; Dumux::Vtk::Precision precision_; };
367 struct VolVarVectorDataInfo { std::function<VolVarsVector(const VV&)> get; std::string name; Dumux::Vtk::Precision precision_; };
368
370
371public:
373 using Field = Vtk::template Field<GridView>;
375 using VolumeVariables = VV;
376
377 VtkOutputModule(const GridVariables& gridVariables,
378 const SolutionVector& sol,
379 const std::string& name,
380 const std::string& paramGroup = "",
381 Dune::VTK::DataMode dm = Dune::VTK::conforming,
382 bool verbose = true)
384 , gridVariables_(gridVariables)
385 , sol_(sol)
386 , velocityOutput_(std::make_shared<VelocityOutputType>())
387 {
388 enableVelocityOutput_ = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddVelocity", false);
389 addProcessRank_ = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddProcessRank", true);
390 }
391
396
403 void addVelocityOutput(std::shared_ptr<VelocityOutputType> velocityOutput)
404 { velocityOutput_ = velocityOutput; }
405
409 void addVolumeVariable(std::function<Scalar(const VolumeVariables&)>&& f,
410 const std::string& name)
411 {
412 // data arrays in the vtk output have to be unique within cell and point data
413 // look if we have a field by the same name with the same codim, if so, replace it
414 for (auto i = 0UL; i < volVarScalarDataInfo_.size(); ++i)
415 {
416 if (volVarScalarDataInfo_[i].name == name)
417 {
418 volVarScalarDataInfo_[i] = VolVarScalarDataInfo{f, name, this->precision()};
419 std::cout << Fmt::format(
420 "VtkOutputModule: Replaced volume variable output \"{}\". "
421 "A field by the same name had already been registered previously.\n",
422 name
423 );
424 return;
425 }
426 }
427
428 // otherwise add it to the end of the fields
429 volVarScalarDataInfo_.push_back(VolVarScalarDataInfo{f, name, this->precision()});
430 }
431
436 template<class VVV = VolVarsVector, typename std::enable_if_t<(VVV::dimension > 1), int> = 0>
437 void addVolumeVariable(std::function<VolVarsVector(const VolumeVariables&)>&& f,
438 const std::string& name)
439 {
440 // data arrays in the vtk output have to be unique within cell and point data
441 // look if we have a field by the same name with the same codim, if so, replace it
442 for (auto i = 0UL; i < volVarVectorDataInfo_.size(); ++i)
443 {
444 if (volVarVectorDataInfo_[i].name == name)
445 {
446 volVarVectorDataInfo_[i] = VolVarVectorDataInfo{f, name, this->precision()};
447 std::cout << Fmt::format(
448 "VtkOutputModule: Replaced volume variable output \"{}\". "
449 "A field by the same name had already been registered previously.\n",
450 name
451 );
452 return;
453 }
454 }
455
456 // otherwise add it to the end of the fields
457 volVarVectorDataInfo_.push_back(VolVarVectorDataInfo{f, name, this->precision()});
458 }
459
460protected:
461 // some return functions for differing implementations to use
462 const auto& problem() const { return curGridVariables_().problem(); }
463 const GridVariables& gridVariables() const { return gridVariables_; }
464 const GridGeometry& gridGeometry() const { return gridVariables_.gridGeometry(); }
465 const SolutionVector& sol() const { return sol_; }
466
467 const std::vector<VolVarScalarDataInfo>& volVarScalarDataInfo() const { return volVarScalarDataInfo_; }
468 const std::vector<VolVarVectorDataInfo>& volVarVectorDataInfo() const { return volVarVectorDataInfo_; }
469
471 const VelocityOutput& velocityOutput() const { return *velocityOutput_; }
472
473 const auto& curGridVariables_() const
474 {
475 if constexpr (Concept::FVGridVariables<GridVariables>)
476 return gridVariables_.curGridVolVars();
477 else
478 return gridVariables_.curGridVars();
479 }
480
481private:
482
484 void writeConforming_(double time, Dune::VTK::OutputType type) override
485 {
486 const Dune::VTK::DataMode dm = Dune::VTK::conforming;
490
491 // instantiate the velocity output
492 using VelocityVector = typename VelocityOutput::VelocityVector;
493 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
494
495 // process rank
496 std::vector<double> rank;
497
498 // volume variable data
499 std::vector<std::vector<Scalar>> volVarScalarData;
500 std::vector<std::vector<VolVarsVector>> volVarVectorData;
501
503 if (!volVarScalarDataInfo_.empty()
504 || !volVarVectorDataInfo_.empty()
505 || !this->fields().empty()
506 || velocityOutput_->enableOutput()
507 || addProcessRank_)
508 {
509 const auto numCells = gridGeometry().gridView().size(0);
510 const auto numDofs = numDofs_();
511
512 // get fields for all volume variables
513 if (!volVarScalarDataInfo_.empty())
514 volVarScalarData.resize(volVarScalarDataInfo_.size(), std::vector<Scalar>(numDofs));
515 if (!volVarVectorDataInfo_.empty())
516 volVarVectorData.resize(volVarVectorDataInfo_.size(), std::vector<VolVarsVector>(numDofs));
517
518 if (velocityOutput_->enableOutput())
519 {
520 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
521 {
522 if (velocityOutput_->fieldType() == VelocityOutput::FieldType::element)
523 velocity[phaseIdx].resize(numCells);
524 else if (velocityOutput_->fieldType() == VelocityOutput::FieldType::vertex)
525 velocity[phaseIdx].resize(numDofs);
526 else
527 {
528 if(isBox && dim == 1)
529 velocity[phaseIdx].resize(numCells);
530 else
531 velocity[phaseIdx].resize(numDofs);
532 }
533 }
534 }
535
536 // maybe allocate space for the process rank
537 if (addProcessRank_) rank.resize(numCells);
538
539 auto fvGeometry = localView(gridGeometry());
540 auto elemVolVars = localView(curGridVariables_());
541 for (const auto& element : elements(gridGeometry().gridView()))
542 {
543 if (!velocityOutput_->enableOutput() &&
544 element.partitionType() != Dune::PartitionType::InteriorEntity)
545 {
546 continue;
547 }
548 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
549 // If velocity output is enabled we need to bind to the whole stencil
550 // otherwise element-local data is sufficient
551 if (velocityOutput_->enableOutput())
552 {
553 fvGeometry.bind(element);
554 elemVolVars.bind(element, fvGeometry, sol_);
555 }
556 else
557 {
558 fvGeometry.bindElement(element);
559 elemVolVars.bindElement(element, fvGeometry, sol_);
560 }
561
562 // velocity output
563 if (velocityOutput_->enableOutput())
564 {
565 if constexpr (Concept::FVGridVariables<GridVariables>)
566 {
567 const auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
568
569 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
570 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
571 }
572 else
573 {
574 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
575 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, phaseIdx);
576 }
577 }
578 else if (element.partitionType() != Dune::PartitionType::InteriorEntity)
579 {
580 continue;
581 }
582
583 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
584 {
585 for (const auto& scv : scvs(fvGeometry))
586 {
587 const auto dofIdxGlobal = scv.dofIndex();
588 const auto& volVars = elemVolVars[scv];
589
590 // get the scalar-valued data
591 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
592 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo_[i].get(volVars);
593
594 // get the vector-valued data
595 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
596 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo_[i].get(volVars);
597 }
598 }
599
601 if (addProcessRank_)
602 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
603 }
604
608
609 // volume variables if any
610 if constexpr (isBox || isPQ1Bubble || isPQ2)
611 {
612 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
613 this->addVertexData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarScalarData[i],
614 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/dim, dm, this->precision()) );
615 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
616 this->addVertexData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarVectorData[i],
617 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()) );
618
619 if constexpr (isPQ1Bubble)
620 {
621 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
622 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarScalarData[i],
623 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0,dm, this->precision()) );
624 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
625 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarVectorData[i],
626 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0,dm, this->precision()) );
627 }
628
629 }
630 else
631 {
632 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
633 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarScalarData[i],
634 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0,dm, this->precision()) );
635 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
636 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarVectorData[i],
637 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0,dm, this->precision()) );
638 }
639
640 // the velocity field
641 if (velocityOutput_->enableOutput())
642 {
643 if (velocityOutput_->fieldType() == VelocityOutput::FieldType::vertex
644 || ( (velocityOutput_->fieldType() == VelocityOutput::FieldType::automatic) && dim > 1 && isBox ))
645 {
646 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
647 this->addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
648 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
649 /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()) );
650 }
651 // cell-centered models
652 else
653 {
654 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
655 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
656 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
657 /*numComp*/dimWorld, /*codim*/0, dm, this->precision()) );
658 }
659 }
660
661 // the process rank
662 if (addProcessRank_)
663 this->addCellData(Field(gridGeometry().gridView(), gridGeometry().elementMapper(), rank, "process rank", 1, 0));
664
665 // also register additional (non-standardized) user fields if any
666 for (auto&& field : this->fields())
667 {
668 if (field.codim() == 0)
669 this->addCellData(field);
670 else if (field.codim() == dim)
671 this->addVertexData(field);
672 else
673 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
674 }
675 }
676
680 this->sequenceWriter().write(time, type);
681
685 this->writer().clear();
686
687 this->addedCellData_.clear();
688 this->addedVertexData_.clear();
689 }
690
692 void writeNonConforming_(double time, Dune::VTK::OutputType type) override
693 {
694 const Dune::VTK::DataMode dm = Dune::VTK::nonconforming;
695
696 // only supports finite-element-like discretization schemes
697 if(!isBox && !isDiamond)
698 DUNE_THROW(Dune::NotImplemented,
699 "Non-conforming output for discretization scheme " << GridGeometry::discMethod
700 );
701
705
706 // check the velocity output
707 if (enableVelocityOutput_ && !velocityOutput_->enableOutput())
708 std::cerr << "Warning! Velocity output was enabled in the input file"
709 << " but no velocity output policy was set for the VTK output module:"
710 << " There will be no velocity output."
711 << " Use the addVelocityOutput member function of the VTK output module." << std::endl;
712 using VelocityVector = typename VelocityOutput::VelocityVector;
713 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
714
715 // process rank
716 std::vector<double> rank;
717
718 // volume variable data (indexing: volvardata/element/localindex)
719 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
720 using VectorDataContainer = std::vector< std::vector<VolVarsVector> >;
721 std::vector< ScalarDataContainer > volVarScalarData;
722 std::vector< VectorDataContainer > volVarVectorData;
723
725 if (!volVarScalarDataInfo_.empty()
726 || !volVarVectorDataInfo_.empty()
727 || !this->fields().empty()
728 || velocityOutput_->enableOutput()
729 || addProcessRank_)
730 {
731 const auto numCells = gridGeometry().gridView().size(0);
732 const auto outputSize = numDofs_();
733
734 // get fields for all volume variables
735 if (!volVarScalarDataInfo_.empty())
736 volVarScalarData.resize(volVarScalarDataInfo_.size(), ScalarDataContainer(numCells));
737 if (!volVarVectorDataInfo_.empty())
738 volVarVectorData.resize(volVarVectorDataInfo_.size(), VectorDataContainer(numCells));
739
740 if (velocityOutput_->enableOutput())
741 {
742 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
743 {
744 if((isBox && dim == 1) || isDiamond)
745 velocity[phaseIdx].resize(numCells);
746 else
747 velocity[phaseIdx].resize(outputSize);
748 }
749 }
750
751 // maybe allocate space for the process rank
752 if (addProcessRank_) rank.resize(numCells);
753
754 // now we go element-local to extract values at local dof locations
755 auto fvGeometry = localView(gridGeometry());
756 auto elemVolVars = localView(curGridVariables_());
757 for (const auto& element : elements(gridGeometry().gridView()))
758 {
759 if (!velocityOutput_->enableOutput() &&
760 element.partitionType() != Dune::PartitionType::InteriorEntity)
761 {
762 continue;
763 }
764 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
765 // If velocity output is enabled we need to bind to the whole stencil
766 // otherwise element-local data is sufficient
767 if (velocityOutput_->enableOutput())
768 {
769 fvGeometry.bind(element);
770 elemVolVars.bind(element, fvGeometry, sol_);
771 }
772 else
773 {
774 fvGeometry.bindElement(element);
775 elemVolVars.bindElement(element, fvGeometry, sol_);
776 }
777
778 const auto numLocalDofs = fvGeometry.numScv();
779 // resize element-local data containers
780 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
781 volVarScalarData[i][eIdxGlobal].resize(numLocalDofs);
782 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
783 volVarVectorData[i][eIdxGlobal].resize(numLocalDofs);
784
785 // velocity output
786 if (velocityOutput_->enableOutput())
787 {
788 if constexpr (Concept::FVGridVariables<GridVariables>)
789 {
790 const auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
791 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
792 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
793 }
794 else
795 {
796 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
797 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, phaseIdx);
798 }
799 }
800 else if (element.partitionType() != Dune::PartitionType::InteriorEntity)
801 {
802 continue;
803 }
804
805 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
806 {
807 for (const auto& scv : scvs(fvGeometry))
808 {
809 const auto& volVars = elemVolVars[scv];
810
811 // get the scalar-valued data
812 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
813 volVarScalarData[i][eIdxGlobal][scv.localDofIndex()] = volVarScalarDataInfo_[i].get(volVars);
814
815 // get the vector-valued data
816 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
817 volVarVectorData[i][eIdxGlobal][scv.localDofIndex()] = volVarVectorDataInfo_[i].get(volVars);
818 }
819 }
820
822 if (addProcessRank_)
823 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
824 }
825
829
830 // volume variables if any
831 static constexpr int dofLocCodim = isDiamond ? 1 : dim;
832 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
833 this->addVertexData(Field(
834 gridGeometry().gridView(), gridGeometry().elementMapper(),
835 volVarScalarData[i], volVarScalarDataInfo_[i].name,
836 /*numComp*/1, /*codim*/dofLocCodim, /*nonconforming*/dm, this->precision()
837 ));
838
839 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
840 this->addVertexData(Field(
841 gridGeometry().gridView(), gridGeometry().elementMapper(),
842 volVarVectorData[i], volVarVectorDataInfo_[i].name,
843 /*numComp*/dimWorld, /*codim*/dofLocCodim, /*nonconforming*/dm, this->precision()
844 ));
845
846 // the velocity field
847 if (velocityOutput_->enableOutput())
848 {
849 // node-wise velocities
850 if (dim > 1 && !isDiamond)
851 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
852 this->addVertexData(Field(
853 gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
854 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
855 /*numComp*/dimWorld, /*codim*/dofLocCodim, dm, this->precision()
856 ));
857
858 // cell-wise velocities
859 else
860 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
861 this->addCellData(Field(
862 gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
863 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
864 /*numComp*/dimWorld, /*codim*/0, dm, this->precision()
865 ));
866 }
867 }
868
869 // the process rank
870 if (addProcessRank_)
871 this->addCellData(Field(
872 gridGeometry().gridView(), gridGeometry().elementMapper(),
873 rank, "process rank", /*numComp*/1, /*codim*/0
874 ));
875
876 // also register additional (non-standardized) user fields if any
877 for (const auto& field : this->fields())
878 {
879 if (field.codim() == 0)
880 this->addCellData(field);
881 else if (field.codim() == dim || field.codim() == 1)
882 this->addVertexData(field);
883 else
884 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
885 }
886
890 this->sequenceWriter().write(time, type);
891
895 this->writer().clear();
896
897 this->addedCellData_.clear();
898 this->addedVertexData_.clear();
899 }
900
902 std::size_t numDofs_() const
903 {
904 // TODO this should actually always be dofMapper.size()
905 // maybe some discretizations needs special treatment (?)
906 if constexpr (isBox || isDiamond || isPQ1Bubble || isPQ2)
907 return gridGeometry().dofMapper().size();
908 else
909 return gridGeometry().elementMapper().size();
910 }
911
912 const GridVariables& gridVariables_;
913 const SolutionVector& sol_;
914
915 std::vector<VolVarScalarDataInfo> volVarScalarDataInfo_;
916 std::vector<VolVarVectorDataInfo> volVarVectorDataInfo_;
917
918 std::shared_ptr<VelocityOutput> velocityOutput_;
919 bool enableVelocityOutput_ = false;
920 bool addProcessRank_ = true;
921};
922
923} // end namespace Dumux
924
925#endif
Definition: io/velocityoutput.hh:31
std::vector< Dune::FieldVector< Scalar, dimWorld > > VelocityVector
Definition: io/velocityoutput.hh:40
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:50
Dumux::Vtk::Precision precision() const
Definition: io/vtkoutputmodule.hh:207
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:80
void addCellData(const Field &field)
Definition: io/vtkoutputmodule.hh:214
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Definition: io/vtkoutputmodule.hh:183
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:58
Dune::VTK::DataMode dataMode() const
Definition: io/vtkoutputmodule.hh:206
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:209
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:210
const std::string & name() const
Definition: io/vtkoutputmodule.hh:205
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:212
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:300
std::vector< std::string > addedCellData_
Definition: io/vtkoutputmodule.hh:234
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:109
bool verbose() const
Definition: io/vtkoutputmodule.hh:204
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:56
std::vector< std::string > addedVertexData_
Definition: io/vtkoutputmodule.hh:235
virtual void writeConforming_(double time, Dune::VTK::OutputType type)
Assembles the fields and adds them to the writer (conforming output)
Definition: io/vtkoutputmodule.hh:239
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:93
void addVertexData(const Field &field)
Definition: io/vtkoutputmodule.hh:223
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:202
void addField(Field &&field)
Add a scalar or vector valued vtk field.
Definition: io/vtkoutputmodule.hh:156
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:344
const auto & problem() const
Definition: io/vtkoutputmodule.hh:462
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:463
void addVolumeVariable(std::function< VolVarsVector(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:437
VV VolumeVariables
export type of the volume variables for the outputfields
Definition: io/vtkoutputmodule.hh:375
void addVolumeVariable(std::function< Scalar(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:409
void addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:403
const auto & curGridVariables_() const
Definition: io/vtkoutputmodule.hh:473
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:471
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:373
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:465
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:464
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:377
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:467
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:468
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!