25#ifndef POROUSMEDIUMFLOW_BOXDFM_VTK_OUTPUT_MODULE_HH
26#define POROUSMEDIUMFLOW_BOXDFM_VTK_OUTPUT_MODULE_HH
30#include <dune/grid/common/gridfactory.hh>
31#include <dune/grid/common/mcmgmapper.hh>
54template<
class Gr
idVariables,
class SolutionVector,
class FractureGr
id>
58 using GridGeometry =
typename GridVariables::GridGeometry;
59 using VV =
typename GridVariables::VolumeVariables;
60 using FluidSystem =
typename VV::FluidSystem;
61 using Scalar =
typename GridVariables::Scalar;
63 using GridView =
typename GridGeometry::GridView;
64 using FractureGridView =
typename FractureGrid::LeafGridView;
65 using FractureMapper = Dune::MultipleCodimMultipleGeomTypeMapper<FractureGridView>;
69 dim = GridView::dimension,
70 dimWorld = GridView::dimensionworld
73 using GridIndexType =
typename GridView::IndexSet::IndexType;
74 using Element =
typename GridView::template Codim<0>::Entity;
75 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
80 static_assert(dim > 1,
"Box-Dfm output only works for dim > 1");
81 static_assert(FractureGrid::dimension == int(dim-1),
"Fracture grid must be of codimension one!");
82 static_assert(FractureGrid::dimensionworld == int(dimWorld),
"Fracture grid has to has the same coordinate dimension!");
83 static_assert(GridGeometry::discMethod ==
DiscretizationMethods::box,
"Box-Dfm output module can only be used with the box scheme!");
87 template<
class FractureGr
idAdapter >
89 const SolutionVector&
sol,
90 const std::string&
name,
91 const FractureGridAdapter& fractureGridAdapter,
93 Dune::VTK::DataMode dm = Dune::VTK::conforming,
98 initializeFracture_(fractureGridAdapter);
110 void write(
double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
116 if (dm == Dune::VTK::conforming)
117 writeConforming_(time, type);
118 else if (dm == Dune::VTK::nonconforming)
119 writeNonConforming_(time, type);
121 DUNE_THROW(Dune::NotImplemented,
"Output for provided vtk data mode");
126 std::cout <<
"Writing output for problem \"" << this->
name() <<
"\". Took " << timer.elapsed() <<
" seconds." << std::endl;
131 void writeConforming_(
double time, Dune::VTK::OutputType type)
138 std::vector<typename ParentType::VelocityOutput::VelocityVector> velocity;
141 static bool addProcessRank = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank");
142 std::vector<double> rank;
145 std::vector<std::vector<Scalar>> volVarScalarData;
146 std::vector<std::vector<Scalar>> volVarScalarDataFracture;
147 std::vector<std::vector<GlobalPosition>> volVarVectorData;
148 std::vector<std::vector<GlobalPosition>> volVarVectorDataFracture;
152 const auto& fractureGridView = fractureGrid_->leafGridView();
159 || !this->fields().empty()
160 || this->velocityOutput().enableOutput()
163 const auto numCells = gridView.size(0);
164 const auto numDofs = gridView.size(dim);
165 const auto numFractureVert = fractureGridView.size(FractureGridView::dimension);
171 volVarScalarDataFracture.resize(
volVarScalarDataInfo.size(), std::vector<Scalar>(numFractureVert));
176 volVarVectorDataFracture.resize(
volVarVectorDataInfo.size(), std::vector<GlobalPosition>(numFractureVert));
181 velocity[phaseIdx].resize(numDofs);
184 if (addProcessRank) rank.resize(numCells);
189 for (
const auto& element : elements(gridView, Dune::Partitions::interior))
191 const auto eIdxGlobal = this->
gridGeometry().elementMapper().index(element);
197 fvGeometry.bind(element);
198 elemVolVars.bind(element, fvGeometry, this->
sol());
202 fvGeometry.bindElement(element);
203 elemVolVars.bindElement(element, fvGeometry, this->
sol());
208 for (
auto&& scv : scvs(fvGeometry))
210 const auto dofIdxGlobal = scv.dofIndex();
211 const auto& volVars = elemVolVars[scv];
213 if (!scv.isOnFracture())
223 volVarScalarDataFracture[i][vertexToFractureVertexIdx_[dofIdxGlobal]] =
volVarScalarDataInfo[i].get(volVars);
225 volVarVectorDataFracture[i][vertexToFractureVertexIdx_[dofIdxGlobal]] =
volVarVectorDataInfo[i].get(volVars);
233 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
241 rank[eIdxGlobal] =
static_cast<double>(gridView.comm().rank());
259 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureVertexMapper_, volVarVectorDataFracture[i],
268 "velocity_" + std::string(this->
velocityOutput().phaseName(phaseIdx)) +
" (m/s)",
269 dimWorld, dim).get() );
275 "process rank", 1, 0).get() );
278 for (
auto&& field : this->
fields())
280 if (field.codim() == 0)
282 else if (field.codim() == dim)
285 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
293 fractureSequenceWriter_->write(time, type);
299 fractureWriter_->clear();
303 void writeNonConforming_(
double time, Dune::VTK::OutputType type)
310 std::vector<typename ParentType::VelocityOutput::VelocityVector> velocity;
313 static bool addProcessRank = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank");
314 std::vector<double> rank;
317 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
318 using VectorDataContainer = std::vector< std::vector<GlobalPosition> >;
319 std::vector< ScalarDataContainer > volVarScalarData;
320 std::vector< ScalarDataContainer > volVarScalarDataFracture;
321 std::vector< VectorDataContainer > volVarVectorData;
322 std::vector< VectorDataContainer > volVarVectorDataFracture;
326 const auto& fractureGridView = fractureGrid_->leafGridView();
333 || !this->fields().empty()
334 || this->velocityOutput().enableOutput()
337 const auto numCells = gridView.size(0);
338 const auto numDofs = gridView.size(dim);
339 const auto numFractureCells = fractureGridView.size(0);
345 volVarScalarDataFracture.resize(
volVarScalarDataInfo.size(), ScalarDataContainer(numFractureCells));
350 volVarVectorDataFracture.resize(
volVarVectorDataInfo.size(), VectorDataContainer(numFractureCells));
355 velocity[phaseIdx].resize(numDofs);
358 if (addProcessRank) rank.resize(numCells);
360 for (
const auto& element : elements(gridView, Dune::Partitions::interior))
362 const auto eIdxGlobal = this->
gridGeometry().elementMapper().index(element);
370 volVarScalarData[i][eIdxGlobal].resize(
numCorners);
372 volVarVectorData[i][eIdxGlobal].resize(
numCorners);
378 fvGeometry.bind(element);
379 elemVolVars.bind(element, fvGeometry, this->
sol());
383 fvGeometry.bindElement(element);
384 elemVolVars.bindElement(element, fvGeometry, this->
sol());
389 for (
auto&& scv : scvs(fvGeometry))
391 const auto& volVars = elemVolVars[scv];
393 if (!scv.isOnFracture())
402 const auto fIdx = scv.facetIndexInElement();
403 const auto& localMap = fractureElementMap_[eIdxGlobal];
404 const auto fracEIdx = std::find_if(localMap.begin(), localMap.end(), [fIdx] (
const auto& p) { return p.first == fIdx; })->second;
416 const auto elemFluxVarsCache =
localView(this->
gridVariables().gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
424 rank[eIdxGlobal] =
static_cast<double>(gridView.comm().rank());
437 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureElementMapper_, volVarScalarDataFracture[i],
447 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureElementMapper_, volVarVectorDataFracture[i],
457 "velocity_" + std::string(this->
velocityOutput().phaseName(phaseIdx)) +
" (m/s)",
458 dimWorld, dim).get() );
464 "process rank", 1, 0).get() );
467 for (
auto&& field : this->
fields())
469 if (field.codim() == 0)
471 else if (field.codim() == dim)
474 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
482 fractureSequenceWriter_->write(time, type);
488 fractureWriter_->clear();
492 template<
class FractureGr
idAdapter >
493 void initializeFracture_(
const FractureGridAdapter& fractureGridAdapter)
497 Dune::GridFactory<FractureGrid> gridFactory;
500 std::size_t fracVertexCount = 0;
501 vertexToFractureVertexIdx_.resize(gridView.size(dim));
502 for (
const auto& v : vertices(gridView))
504 if (fractureGridAdapter.isOnFacetGrid(v))
506 gridFactory.insertVertex(v.geometry().center());
507 vertexToFractureVertexIdx_[
gridGeometry.vertexMapper().index(v)] = fracVertexCount++;
512 std::size_t fractureElementCount = 0;
513 fractureElementMap_.resize(gridView.size(0));
514 std::set< std::pair<GridIndexType, unsigned int> > handledFacets;
515 for (
const auto& element : elements(gridView))
517 const auto eIdxGlobal =
gridGeometry.elementMapper().index(element);
518 const auto refElement = referenceElement(element);
520 for (
const auto& is : intersections(gridView, element))
523 const auto& isGeometry = is.geometry();
525 const auto indexInInside = is.indexInInside();
527 std::vector<GridIndexType> isVertexIndices(
numCorners);
529 isVertexIndices[i] =
gridGeometry.vertexMapper().subIndex(element,
530 refElement.subEntity(indexInInside, 1, i, dim),
534 bool insertFacet =
false;
535 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
541 const auto outsideEIdx =
gridGeometry.elementMapper().index(is.outside());
542 const auto idxInOutside = is.indexInOutside();
543 const auto pair = std::make_pair(outsideEIdx, idxInOutside);
544 if (handledFacets.count( pair ) != 0)
549 const auto& outsideMap = fractureElementMap_[outsideEIdx];
550 auto it = std::find_if(outsideMap.begin(), outsideMap.end(), [idxInOutside] (
const auto& p) { return p.first == idxInOutside; });
551 fractureElementMap_[eIdxGlobal].push_back( std::make_pair(indexInInside, it->second) );
559 std::for_each( isVertexIndices.begin(),
560 isVertexIndices.end(),
561 [&] (
auto& idx) { idx = this->vertexToFractureVertexIdx_[idx]; } );
564 gridFactory.insertElement(isGeometry.type(), isVertexIndices);
567 handledFacets.insert( std::make_pair(eIdxGlobal, indexInInside) );
568 fractureElementMap_[eIdxGlobal].push_back( std::make_pair(indexInInside, fractureElementCount) );
569 fractureElementCount++;
575 fractureGrid_ = std::shared_ptr<FractureGrid>(gridFactory.createGrid());
578 const auto& fractureGridView = fractureGrid_->leafGridView();
579 fractureVertexMapper_ = std::make_unique<FractureMapper>(fractureGridView, Dune::mcmgVertexLayout());
580 fractureElementMapper_ = std::make_unique<FractureMapper>(fractureGridView, Dune::mcmgElementLayout());
583 std::vector<GridIndexType> insToVertexIdx(fractureGridView.size(FractureGridView::dimension));
584 std::vector<GridIndexType> insToElemIdx(fractureGridView.size(0));
585 for (
const auto& v : vertices(fractureGridView)) insToVertexIdx[ gridFactory.insertionIndex(v) ] = fractureVertexMapper_->index(v);
586 for (
const auto& e : elements(fractureGridView)) insToElemIdx[ gridFactory.insertionIndex(e) ] = fractureElementMapper_->index(e);
589 for (GridIndexType dofIdx = 0; dofIdx < gridView.size(GridView::dimension); ++dofIdx)
591 vertexToFractureVertexIdx_[dofIdx] = insToVertexIdx[ vertexToFractureVertexIdx_[dofIdx] ];
594 for (
auto& elemLocalMap : fractureElementMap_)
595 for (
auto& dataPair : elemLocalMap)
596 dataPair.second = insToElemIdx[ dataPair.second ];
599 fractureWriter_ = std::make_shared< Dune::VTKWriter<FractureGridView> >(fractureGridView, this->
dataMode());
600 fractureSequenceWriter_ = std::make_unique< Dune::VTKSequenceWriter<FractureGridView> >(fractureWriter_, this->
name() +
"_fracture");
603 std::shared_ptr<FractureGrid> fractureGrid_;
605 std::unique_ptr<FractureMapper> fractureVertexMapper_;
606 std::unique_ptr<FractureMapper> fractureElementMapper_;
608 std::shared_ptr<Dune::VTKWriter<FractureGridView>> fractureWriter_;
609 std::unique_ptr< Dune::VTKSequenceWriter<FractureGridView> > fractureSequenceWriter_;
612 std::vector<GridIndexType> vertexToFractureVertexIdx_;
615 std::vector< std::vector<std::pair<GridIndexType, unsigned int>> > fractureElementMap_;
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
constexpr Box box
Definition: method.hh:136
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:232
virtual int numFluidPhases() const
returns the number of phases
Definition: io/velocityoutput.hh:78
virtual void calculateVelocity(VelocityVector &velocity, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVarsCache &elemFluxVarsCache, int phaseIdx) const
Definition: io/velocityoutput.hh:82
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:61
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:91
Dune::VTK::DataMode dataMode() const
Definition: io/vtkoutputmodule.hh:200
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:203
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:204
const std::string & name() const
Definition: io/vtkoutputmodule.hh:199
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:206
bool verbose() const
Definition: io/vtkoutputmodule.hh:198
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:312
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:396
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:404
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:340
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:398
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:397
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:400
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:401
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: porousmediumflow/boxdfm/vtkoutputmodule.hh:56
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Writing data.
Definition: porousmediumflow/boxdfm/vtkoutputmodule.hh:110
BoxDfmVtkOutputModule(const GridVariables &gridVariables, const SolutionVector &sol, const std::string &name, const FractureGridAdapter &fractureGridAdapter, const std::string ¶mGroup="", Dune::VTK::DataMode dm=Dune::VTK::conforming, bool verbose=true)
The constructor.
Definition: porousmediumflow/boxdfm/vtkoutputmodule.hh:88
A VTK output module to simplify writing dumux simulation data to VTK format.