version 3.11-dev
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
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
364 struct VolVarScalarDataInfo { std::function<Scalar(const VV&)> get; std::string name; Dumux::Vtk::Precision precision_; };
365 struct VolVarVectorDataInfo { std::function<VolVarsVector(const VV&)> get; std::string name; Dumux::Vtk::Precision precision_; };
366
368
369public:
371 using Field = Vtk::template Field<GridView>;
373 using VolumeVariables = VV;
374
375 VtkOutputModule(const GridVariables& gridVariables,
376 const SolutionVector& sol,
377 const std::string& name,
378 const std::string& paramGroup = "",
379 Dune::VTK::DataMode dm = Dune::VTK::conforming,
380 bool verbose = true)
382 , gridVariables_(gridVariables)
383 , sol_(sol)
384 , velocityOutput_(std::make_shared<VelocityOutputType>())
385 {
386 enableVelocityOutput_ = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddVelocity", false);
387 addProcessRank_ = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddProcessRank", true);
388 }
389
394
401 void addVelocityOutput(std::shared_ptr<VelocityOutputType> velocityOutput)
402 { velocityOutput_ = velocityOutput; }
403
407 void addVolumeVariable(std::function<Scalar(const VolumeVariables&)>&& f,
408 const std::string& name)
409 {
410 // data arrays in the vtk output have to be unique within cell and point data
411 // look if we have a field by the same name with the same codim, if so, replace it
412 for (auto i = 0UL; i < volVarScalarDataInfo_.size(); ++i)
413 {
414 if (volVarScalarDataInfo_[i].name == name)
415 {
416 volVarScalarDataInfo_[i] = VolVarScalarDataInfo{f, name, this->precision()};
417 std::cout << Fmt::format(
418 "VtkOutputModule: Replaced volume variable output \"{}\". "
419 "A field by the same name had already been registered previously.\n",
420 name
421 );
422 return;
423 }
424 }
425
426 // otherwise add it to the end of the fields
427 volVarScalarDataInfo_.push_back(VolVarScalarDataInfo{f, name, this->precision()});
428 }
429
434 template<class VVV = VolVarsVector, typename std::enable_if_t<(VVV::dimension > 1), int> = 0>
435 void addVolumeVariable(std::function<VolVarsVector(const VolumeVariables&)>&& f,
436 const std::string& name)
437 {
438 // data arrays in the vtk output have to be unique within cell and point data
439 // look if we have a field by the same name with the same codim, if so, replace it
440 for (auto i = 0UL; i < volVarVectorDataInfo_.size(); ++i)
441 {
442 if (volVarVectorDataInfo_[i].name == name)
443 {
444 volVarVectorDataInfo_[i] = VolVarVectorDataInfo{f, name, this->precision()};
445 std::cout << Fmt::format(
446 "VtkOutputModule: Replaced volume variable output \"{}\". "
447 "A field by the same name had already been registered previously.\n",
448 name
449 );
450 return;
451 }
452 }
453
454 // otherwise add it to the end of the fields
455 volVarVectorDataInfo_.push_back(VolVarVectorDataInfo{f, name, this->precision()});
456 }
457
458protected:
459 // some return functions for differing implementations to use
460 const auto& problem() const { return gridVariables_.curGridVolVars().problem(); }
461 const GridVariables& gridVariables() const { return gridVariables_; }
462 const GridGeometry& gridGeometry() const { return gridVariables_.gridGeometry(); }
463 const SolutionVector& sol() const { return sol_; }
464
465 const std::vector<VolVarScalarDataInfo>& volVarScalarDataInfo() const { return volVarScalarDataInfo_; }
466 const std::vector<VolVarVectorDataInfo>& volVarVectorDataInfo() const { return volVarVectorDataInfo_; }
467
469 const VelocityOutput& velocityOutput() const { return *velocityOutput_; }
470
471private:
472
474 void writeConforming_(double time, Dune::VTK::OutputType type) override
475 {
476 const Dune::VTK::DataMode dm = Dune::VTK::conforming;
480
481 // instantiate the velocity output
482 using VelocityVector = typename VelocityOutput::VelocityVector;
483 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
484
485 // process rank
486 std::vector<double> rank;
487
488 // volume variable data
489 std::vector<std::vector<Scalar>> volVarScalarData;
490 std::vector<std::vector<VolVarsVector>> volVarVectorData;
491
493 if (!volVarScalarDataInfo_.empty()
494 || !volVarVectorDataInfo_.empty()
495 || !this->fields().empty()
496 || velocityOutput_->enableOutput()
497 || addProcessRank_)
498 {
499 const auto numCells = gridGeometry().gridView().size(0);
500 const auto numDofs = numDofs_();
501
502 // get fields for all volume variables
503 if (!volVarScalarDataInfo_.empty())
504 volVarScalarData.resize(volVarScalarDataInfo_.size(), std::vector<Scalar>(numDofs));
505 if (!volVarVectorDataInfo_.empty())
506 volVarVectorData.resize(volVarVectorDataInfo_.size(), std::vector<VolVarsVector>(numDofs));
507
508 if (velocityOutput_->enableOutput())
509 {
510 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
511 {
512 if (velocityOutput_->fieldType() == VelocityOutput::FieldType::element)
513 velocity[phaseIdx].resize(numCells);
514 else if (velocityOutput_->fieldType() == VelocityOutput::FieldType::vertex)
515 velocity[phaseIdx].resize(numDofs);
516 else
517 {
518 if(isBox && dim == 1)
519 velocity[phaseIdx].resize(numCells);
520 else
521 velocity[phaseIdx].resize(numDofs);
522 }
523 }
524 }
525
526 // maybe allocate space for the process rank
527 if (addProcessRank_) rank.resize(numCells);
528
529 auto fvGeometry = localView(gridGeometry());
530 auto elemVolVars = localView(gridVariables_.curGridVolVars());
531 for (const auto& element : elements(gridGeometry().gridView(), Dune::Partitions::interior))
532 {
533 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
534 // If velocity output is enabled we need to bind to the whole stencil
535 // otherwise element-local data is sufficient
536 if (velocityOutput_->enableOutput())
537 {
538 fvGeometry.bind(element);
539 elemVolVars.bind(element, fvGeometry, sol_);
540 }
541 else
542 {
543 fvGeometry.bindElement(element);
544 elemVolVars.bindElement(element, fvGeometry, sol_);
545 }
546
547 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
548 {
549 for (const auto& scv : scvs(fvGeometry))
550 {
551 const auto dofIdxGlobal = scv.dofIndex();
552 const auto& volVars = elemVolVars[scv];
553
554 // get the scalar-valued data
555 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
556 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo_[i].get(volVars);
557
558 // get the vector-valued data
559 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
560 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo_[i].get(volVars);
561 }
562 }
563
564 // velocity output
565 if (velocityOutput_->enableOutput())
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
574 if (addProcessRank_)
575 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
576 }
577
581
582 // volume variables if any
583 if constexpr (isBox || isPQ1Bubble)
584 {
585 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
586 this->addVertexData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarScalarData[i],
587 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/dim, dm, this->precision()) );
588 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
589 this->addVertexData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarVectorData[i],
590 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()) );
591
592 if constexpr (isPQ1Bubble)
593 {
594 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
595 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarScalarData[i],
596 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0,dm, this->precision()) );
597 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
598 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().dofMapper(), volVarVectorData[i],
599 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0,dm, this->precision()) );
600 }
601
602 }
603 else
604 {
605 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
606 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarScalarData[i],
607 volVarScalarDataInfo_[i].name, /*numComp*/1, /*codim*/0,dm, this->precision()) );
608 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
609 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), volVarVectorData[i],
610 volVarVectorDataInfo_[i].name, /*numComp*/dimWorld, /*codim*/0,dm, this->precision()) );
611 }
612
613 // the velocity field
614 if (velocityOutput_->enableOutput())
615 {
616 if (velocityOutput_->fieldType() == VelocityOutput::FieldType::vertex
617 || ( (velocityOutput_->fieldType() == VelocityOutput::FieldType::automatic) && dim > 1 && isBox ))
618 {
619 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
620 this->addVertexData( Field(gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
621 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
622 /*numComp*/dimWorld, /*codim*/dim, dm, this->precision()) );
623 }
624 // cell-centered models
625 else
626 {
627 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
628 this->addCellData( Field(gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
629 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
630 /*numComp*/dimWorld, /*codim*/0, dm, this->precision()) );
631 }
632 }
633
634 // the process rank
635 if (addProcessRank_)
636 this->addCellData(Field(gridGeometry().gridView(), gridGeometry().elementMapper(), rank, "process rank", 1, 0));
637
638 // also register additional (non-standardized) user fields if any
639 for (auto&& field : this->fields())
640 {
641 if (field.codim() == 0)
642 this->addCellData(field);
643 else if (field.codim() == dim)
644 this->addVertexData(field);
645 else
646 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
647 }
648 }
649
653 this->sequenceWriter().write(time, type);
654
658 this->writer().clear();
659
660 this->addedCellData_.clear();
661 this->addedVertexData_.clear();
662 }
663
665 void writeNonConforming_(double time, Dune::VTK::OutputType type) override
666 {
667 const Dune::VTK::DataMode dm = Dune::VTK::nonconforming;
668
669 // only supports finite-element-like discretization schemes
670 if(!isBox && !isDiamond)
671 DUNE_THROW(Dune::NotImplemented,
672 "Non-conforming output for discretization scheme " << GridGeometry::discMethod
673 );
674
678
679 // check the velocity output
680 if (enableVelocityOutput_ && !velocityOutput_->enableOutput())
681 std::cerr << "Warning! Velocity output was enabled in the input file"
682 << " but no velocity output policy was set for the VTK output module:"
683 << " There will be no velocity output."
684 << " Use the addVelocityOutput member function of the VTK output module." << std::endl;
685 using VelocityVector = typename VelocityOutput::VelocityVector;
686 std::vector<VelocityVector> velocity(velocityOutput_->numFluidPhases());
687
688 // process rank
689 std::vector<double> rank;
690
691 // volume variable data (indexing: volvardata/element/localindex)
692 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
693 using VectorDataContainer = std::vector< std::vector<VolVarsVector> >;
694 std::vector< ScalarDataContainer > volVarScalarData;
695 std::vector< VectorDataContainer > volVarVectorData;
696
698 if (!volVarScalarDataInfo_.empty()
699 || !volVarVectorDataInfo_.empty()
700 || !this->fields().empty()
701 || velocityOutput_->enableOutput()
702 || addProcessRank_)
703 {
704 const auto numCells = gridGeometry().gridView().size(0);
705 const auto outputSize = numDofs_();
706
707 // get fields for all volume variables
708 if (!volVarScalarDataInfo_.empty())
709 volVarScalarData.resize(volVarScalarDataInfo_.size(), ScalarDataContainer(numCells));
710 if (!volVarVectorDataInfo_.empty())
711 volVarVectorData.resize(volVarVectorDataInfo_.size(), VectorDataContainer(numCells));
712
713 if (velocityOutput_->enableOutput())
714 {
715 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
716 {
717 if((isBox && dim == 1) || isDiamond)
718 velocity[phaseIdx].resize(numCells);
719 else
720 velocity[phaseIdx].resize(outputSize);
721 }
722 }
723
724 // maybe allocate space for the process rank
725 if (addProcessRank_) rank.resize(numCells);
726
727 // now we go element-local to extract values at local dof locations
728 auto fvGeometry = localView(gridGeometry());
729 auto elemVolVars = localView(gridVariables_.curGridVolVars());
730 for (const auto& element : elements(gridGeometry().gridView(), Dune::Partitions::interior))
731 {
732 const auto eIdxGlobal = gridGeometry().elementMapper().index(element);
733 // If velocity output is enabled we need to bind to the whole stencil
734 // otherwise element-local data is sufficient
735 if (velocityOutput_->enableOutput())
736 {
737 fvGeometry.bind(element);
738 elemVolVars.bind(element, fvGeometry, sol_);
739 }
740 else
741 {
742 fvGeometry.bindElement(element);
743 elemVolVars.bindElement(element, fvGeometry, sol_);
744 }
745
746 const auto numLocalDofs = fvGeometry.numScv();
747 // resize element-local data containers
748 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
749 volVarScalarData[i][eIdxGlobal].resize(numLocalDofs);
750 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
751 volVarVectorData[i][eIdxGlobal].resize(numLocalDofs);
752
753 if (!volVarScalarDataInfo_.empty() || !volVarVectorDataInfo_.empty())
754 {
755 for (const auto& scv : scvs(fvGeometry))
756 {
757 const auto& volVars = elemVolVars[scv];
758
759 // get the scalar-valued data
760 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
761 volVarScalarData[i][eIdxGlobal][scv.localDofIndex()] = volVarScalarDataInfo_[i].get(volVars);
762
763 // get the vector-valued data
764 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
765 volVarVectorData[i][eIdxGlobal][scv.localDofIndex()] = volVarVectorDataInfo_[i].get(volVars);
766 }
767 }
768
769 // velocity output
770 if (velocityOutput_->enableOutput())
771 {
772 const auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
773 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
774 velocityOutput_->calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
775 }
776
778 if (addProcessRank_)
779 rank[eIdxGlobal] = static_cast<double>(gridGeometry().gridView().comm().rank());
780 }
781
785
786 // volume variables if any
787 static constexpr int dofLocCodim = isDiamond ? 1 : dim;
788 for (std::size_t i = 0; i < volVarScalarDataInfo_.size(); ++i)
789 this->addVertexData(Field(
790 gridGeometry().gridView(), gridGeometry().elementMapper(),
791 volVarScalarData[i], volVarScalarDataInfo_[i].name,
792 /*numComp*/1, /*codim*/dofLocCodim, /*nonconforming*/dm, this->precision()
793 ));
794
795 for (std::size_t i = 0; i < volVarVectorDataInfo_.size(); ++i)
796 this->addVertexData(Field(
797 gridGeometry().gridView(), gridGeometry().elementMapper(),
798 volVarVectorData[i], volVarVectorDataInfo_[i].name,
799 /*numComp*/dimWorld, /*codim*/dofLocCodim, /*nonconforming*/dm, this->precision()
800 ));
801
802 // the velocity field
803 if (velocityOutput_->enableOutput())
804 {
805 // node-wise velocities
806 if (dim > 1 && !isDiamond)
807 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
808 this->addVertexData(Field(
809 gridGeometry().gridView(), gridGeometry().vertexMapper(), velocity[phaseIdx],
810 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
811 /*numComp*/dimWorld, /*codim*/dofLocCodim, dm, this->precision()
812 ));
813
814 // cell-wise velocities
815 else
816 for (int phaseIdx = 0; phaseIdx < velocityOutput_->numFluidPhases(); ++phaseIdx)
817 this->addCellData(Field(
818 gridGeometry().gridView(), gridGeometry().elementMapper(), velocity[phaseIdx],
819 "velocity_" + velocityOutput_->phaseName(phaseIdx) + " (m/s)",
820 /*numComp*/dimWorld, /*codim*/0, dm, this->precision()
821 ));
822 }
823 }
824
825 // the process rank
826 if (addProcessRank_)
827 this->addCellData(Field(
828 gridGeometry().gridView(), gridGeometry().elementMapper(),
829 rank, "process rank", /*numComp*/1, /*codim*/0
830 ));
831
832 // also register additional (non-standardized) user fields if any
833 for (const auto& field : this->fields())
834 {
835 if (field.codim() == 0)
836 this->addCellData(field);
837 else if (field.codim() == dim || field.codim() == 1)
838 this->addVertexData(field);
839 else
840 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
841 }
842
846 this->sequenceWriter().write(time, type);
847
851 this->writer().clear();
852
853 this->addedCellData_.clear();
854 this->addedVertexData_.clear();
855 }
856
858 std::size_t numDofs_() const
859 {
860 // TODO this should actually always be dofMapper.size()
861 // maybe some discretizations needs special treatment (?)
862 if constexpr (isBox || isDiamond || isPQ1Bubble)
863 return gridGeometry().dofMapper().size();
864 else
865 return gridGeometry().elementMapper().size();
866 }
867
868 const GridVariables& gridVariables_;
869 const SolutionVector& sol_;
870
871 std::vector<VolVarScalarDataInfo> volVarScalarDataInfo_;
872 std::vector<VolVarVectorDataInfo> volVarVectorDataInfo_;
873
874 std::shared_ptr<VelocityOutput> velocityOutput_;
875 bool enableVelocityOutput_ = false;
876 bool addProcessRank_ = true;
877};
878
879} // end namespace Dumux
880
881#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
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
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:460
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:461
void addVolumeVariable(std::function< VolVarsVector(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:435
VV VolumeVariables
export type of the volume variables for the outputfields
Definition: io/vtkoutputmodule.hh:373
void addVolumeVariable(std::function< Scalar(const VolumeVariables &)> &&f, const std::string &name)
Definition: io/vtkoutputmodule.hh:407
void addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:401
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:469
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:371
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:463
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:462
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:375
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:465
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:466
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!