13#ifndef POROUSMEDIUMFLOW_BOXDFM_VTK_OUTPUT_MODULE_HH
14#define POROUSMEDIUMFLOW_BOXDFM_VTK_OUTPUT_MODULE_HH
18#include <dune/grid/common/gridfactory.hh>
19#include <dune/grid/common/mcmgmapper.hh>
42template<
class Gr
idVariables,
class SolutionVector,
class FractureGr
id>
46 using GridGeometry =
typename GridVariables::GridGeometry;
47 using VV =
typename GridVariables::VolumeVariables;
48 using FluidSystem =
typename VV::FluidSystem;
49 using Scalar =
typename GridVariables::Scalar;
51 using GridView =
typename GridGeometry::GridView;
52 using FractureGridView =
typename FractureGrid::LeafGridView;
53 using FractureMapper = Dune::MultipleCodimMultipleGeomTypeMapper<FractureGridView>;
57 dim = GridView::dimension,
58 dimWorld = GridView::dimensionworld
61 using GridIndexType =
typename GridView::IndexSet::IndexType;
62 using Element =
typename GridView::template Codim<0>::Entity;
63 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
68 static_assert(dim > 1,
"Box-Dfm output only works for dim > 1");
69 static_assert(FractureGrid::dimension == int(dim-1),
"Fracture grid must be of codimension one!");
70 static_assert(FractureGrid::dimensionworld == int(dimWorld),
"Fracture grid has to has the same coordinate dimension!");
71 static_assert(GridGeometry::discMethod ==
DiscretizationMethods::box,
"Box-Dfm output module can only be used with the box scheme!");
75 template<
class FractureGr
idAdapter >
77 const SolutionVector&
sol,
78 const std::string&
name,
79 const FractureGridAdapter& fractureGridAdapter,
81 Dune::VTK::DataMode dm = Dune::VTK::conforming,
86 initializeFracture_(fractureGridAdapter);
98 void write(
double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
104 if (dm == Dune::VTK::conforming)
105 writeConforming_(time, type);
106 else if (dm == Dune::VTK::nonconforming)
107 writeNonConforming_(time, type);
109 DUNE_THROW(Dune::NotImplemented,
"Output for provided vtk data mode");
114 std::cout <<
"Writing output for problem \"" << this->
name() <<
"\". Took " << timer.elapsed() <<
" seconds." << std::endl;
119 void writeConforming_(
double time, Dune::VTK::OutputType type)
126 std::vector<typename ParentType::VelocityOutput::VelocityVector> velocity;
129 static bool addProcessRank = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank");
130 std::vector<double> rank;
133 std::vector<std::vector<Scalar>> volVarScalarData;
134 std::vector<std::vector<Scalar>> volVarScalarDataFracture;
135 std::vector<std::vector<GlobalPosition>> volVarVectorData;
136 std::vector<std::vector<GlobalPosition>> volVarVectorDataFracture;
140 const auto& fractureGridView = fractureGrid_->leafGridView();
147 || !this->fields().empty()
148 || this->velocityOutput().enableOutput()
151 const auto numCells = gridView.size(0);
152 const auto numDofs = gridView.size(dim);
153 const auto numFractureVert = fractureGridView.size(FractureGridView::dimension);
159 volVarScalarDataFracture.resize(
volVarScalarDataInfo.size(), std::vector<Scalar>(numFractureVert));
164 volVarVectorDataFracture.resize(
volVarVectorDataInfo.size(), std::vector<GlobalPosition>(numFractureVert));
169 velocity[phaseIdx].resize(numDofs);
172 if (addProcessRank) rank.resize(numCells);
177 for (
const auto& element : elements(gridView, Dune::Partitions::interior))
179 const auto eIdxGlobal = this->
gridGeometry().elementMapper().index(element);
185 fvGeometry.bind(element);
186 elemVolVars.bind(element, fvGeometry, this->
sol());
190 fvGeometry.bindElement(element);
191 elemVolVars.bindElement(element, fvGeometry, this->
sol());
196 for (
auto&& scv : scvs(fvGeometry))
198 const auto dofIdxGlobal = scv.dofIndex();
199 const auto& volVars = elemVolVars[scv];
201 if (!scv.isOnFracture())
211 volVarScalarDataFracture[i][vertexToFractureVertexIdx_[dofIdxGlobal]] =
volVarScalarDataInfo[i].get(volVars);
213 volVarVectorDataFracture[i][vertexToFractureVertexIdx_[dofIdxGlobal]] =
volVarVectorDataInfo[i].get(volVars);
221 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
229 rank[eIdxGlobal] =
static_cast<double>(gridView.comm().rank());
247 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureVertexMapper_, volVarVectorDataFracture[i],
256 "velocity_" + std::string(this->
velocityOutput().phaseName(phaseIdx)) +
" (m/s)",
257 dimWorld, dim).get() );
263 "process rank", 1, 0).get() );
266 for (
auto&& field : this->
fields())
268 if (field.codim() == 0)
270 else if (field.codim() == dim)
273 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
281 fractureSequenceWriter_->write(time, type);
287 fractureWriter_->clear();
291 void writeNonConforming_(
double time, Dune::VTK::OutputType type)
298 std::vector<typename ParentType::VelocityOutput::VelocityVector> velocity;
301 static bool addProcessRank = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank");
302 std::vector<double> rank;
305 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
306 using VectorDataContainer = std::vector< std::vector<GlobalPosition> >;
307 std::vector< ScalarDataContainer > volVarScalarData;
308 std::vector< ScalarDataContainer > volVarScalarDataFracture;
309 std::vector< VectorDataContainer > volVarVectorData;
310 std::vector< VectorDataContainer > volVarVectorDataFracture;
314 const auto& fractureGridView = fractureGrid_->leafGridView();
321 || !this->fields().empty()
322 || this->velocityOutput().enableOutput()
325 const auto numCells = gridView.size(0);
326 const auto numDofs = gridView.size(dim);
327 const auto numFractureCells = fractureGridView.size(0);
333 volVarScalarDataFracture.resize(
volVarScalarDataInfo.size(), ScalarDataContainer(numFractureCells));
338 volVarVectorDataFracture.resize(
volVarVectorDataInfo.size(), VectorDataContainer(numFractureCells));
343 velocity[phaseIdx].resize(numDofs);
346 if (addProcessRank) rank.resize(numCells);
348 for (
const auto& element : elements(gridView, Dune::Partitions::interior))
350 const auto eIdxGlobal = this->
gridGeometry().elementMapper().index(element);
358 volVarScalarData[i][eIdxGlobal].resize(
numCorners);
360 volVarVectorData[i][eIdxGlobal].resize(
numCorners);
366 fvGeometry.bind(element);
367 elemVolVars.bind(element, fvGeometry, this->
sol());
371 fvGeometry.bindElement(element);
372 elemVolVars.bindElement(element, fvGeometry, this->
sol());
377 for (
auto&& scv : scvs(fvGeometry))
379 const auto& volVars = elemVolVars[scv];
381 if (!scv.isOnFracture())
390 const auto fIdx = scv.facetIndexInElement();
391 const auto& localMap = fractureElementMap_[eIdxGlobal];
392 const auto fracEIdx = std::find_if(localMap.begin(), localMap.end(), [fIdx] (
const auto& p) { return p.first == fIdx; })->second;
404 const auto elemFluxVarsCache =
localView(this->
gridVariables().gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
412 rank[eIdxGlobal] =
static_cast<double>(gridView.comm().rank());
425 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureElementMapper_, volVarScalarDataFracture[i],
435 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureElementMapper_, volVarVectorDataFracture[i],
445 "velocity_" + std::string(this->
velocityOutput().phaseName(phaseIdx)) +
" (m/s)",
446 dimWorld, dim).get() );
452 "process rank", 1, 0).get() );
455 for (
auto&& field : this->
fields())
457 if (field.codim() == 0)
459 else if (field.codim() == dim)
462 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
470 fractureSequenceWriter_->write(time, type);
476 fractureWriter_->clear();
480 template<
class FractureGr
idAdapter >
481 void initializeFracture_(
const FractureGridAdapter& fractureGridAdapter)
485 Dune::GridFactory<FractureGrid> gridFactory;
488 std::size_t fracVertexCount = 0;
489 vertexToFractureVertexIdx_.resize(gridView.size(dim));
490 for (
const auto& v : vertices(gridView))
492 if (fractureGridAdapter.isOnFacetGrid(v))
494 gridFactory.insertVertex(v.geometry().center());
495 vertexToFractureVertexIdx_[
gridGeometry.vertexMapper().index(v)] = fracVertexCount++;
500 std::size_t fractureElementCount = 0;
501 fractureElementMap_.resize(gridView.size(0));
502 std::set< std::pair<GridIndexType, unsigned int> > handledFacets;
503 for (
const auto& element : elements(gridView))
505 const auto eIdxGlobal =
gridGeometry.elementMapper().index(element);
506 const auto refElement = referenceElement(element);
508 for (
const auto& is : intersections(gridView, element))
511 const auto& isGeometry = is.geometry();
513 const auto indexInInside = is.indexInInside();
515 std::vector<GridIndexType> isVertexIndices(
numCorners);
517 isVertexIndices[i] =
gridGeometry.vertexMapper().subIndex(element,
518 refElement.subEntity(indexInInside, 1, i, dim),
522 bool insertFacet =
false;
523 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
529 const auto outsideEIdx =
gridGeometry.elementMapper().index(is.outside());
530 const auto idxInOutside = is.indexInOutside();
531 const auto pair = std::make_pair(outsideEIdx, idxInOutside);
532 if (handledFacets.count( pair ) != 0)
537 const auto& outsideMap = fractureElementMap_[outsideEIdx];
538 auto it = std::find_if(outsideMap.begin(), outsideMap.end(), [idxInOutside] (
const auto& p) { return p.first == idxInOutside; });
539 fractureElementMap_[eIdxGlobal].push_back( std::make_pair(indexInInside, it->second) );
547 std::for_each( isVertexIndices.begin(),
548 isVertexIndices.end(),
549 [&] (
auto& idx) { idx = this->vertexToFractureVertexIdx_[idx]; } );
552 gridFactory.insertElement(isGeometry.type(), isVertexIndices);
555 handledFacets.insert( std::make_pair(eIdxGlobal, indexInInside) );
556 fractureElementMap_[eIdxGlobal].push_back( std::make_pair(indexInInside, fractureElementCount) );
557 fractureElementCount++;
563 fractureGrid_ = std::shared_ptr<FractureGrid>(gridFactory.createGrid());
566 const auto& fractureGridView = fractureGrid_->leafGridView();
567 fractureVertexMapper_ = std::make_unique<FractureMapper>(fractureGridView, Dune::mcmgVertexLayout());
568 fractureElementMapper_ = std::make_unique<FractureMapper>(fractureGridView, Dune::mcmgElementLayout());
571 std::vector<GridIndexType> insToVertexIdx(fractureGridView.size(FractureGridView::dimension));
572 std::vector<GridIndexType> insToElemIdx(fractureGridView.size(0));
573 for (
const auto& v : vertices(fractureGridView)) insToVertexIdx[ gridFactory.insertionIndex(v) ] = fractureVertexMapper_->index(v);
574 for (
const auto& e : elements(fractureGridView)) insToElemIdx[ gridFactory.insertionIndex(e) ] = fractureElementMapper_->index(e);
577 for (GridIndexType dofIdx = 0; dofIdx < gridView.size(GridView::dimension); ++dofIdx)
579 vertexToFractureVertexIdx_[dofIdx] = insToVertexIdx[ vertexToFractureVertexIdx_[dofIdx] ];
582 for (
auto& elemLocalMap : fractureElementMap_)
583 for (
auto& dataPair : elemLocalMap)
584 dataPair.second = insToElemIdx[ dataPair.second ];
587 fractureWriter_ = std::make_shared< Dune::VTKWriter<FractureGridView> >(fractureGridView, this->
dataMode());
588 fractureSequenceWriter_ = std::make_unique< Dune::VTKSequenceWriter<FractureGridView> >(fractureWriter_, this->
name() +
"_fracture");
591 std::shared_ptr<FractureGrid> fractureGrid_;
593 std::unique_ptr<FractureMapper> fractureVertexMapper_;
594 std::unique_ptr<FractureMapper> fractureElementMapper_;
596 std::shared_ptr<Dune::VTKWriter<FractureGridView>> fractureWriter_;
597 std::unique_ptr< Dune::VTKSequenceWriter<FractureGridView> > fractureSequenceWriter_;
600 std::vector<GridIndexType> vertexToFractureVertexIdx_;
603 std::vector< std::vector<std::pair<GridIndexType, unsigned int>> > fractureElementMap_;
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: porousmediumflow/boxdfm/vtkoutputmodule.hh:44
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Writing data.
Definition: porousmediumflow/boxdfm/vtkoutputmodule.hh:98
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:76
virtual int numFluidPhases() const
returns the number of phases
Definition: io/velocityoutput.hh:66
virtual void calculateVelocity(VelocityVector &velocity, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVarsCache &elemFluxVarsCache, int phaseIdx) const
Definition: io/velocityoutput.hh:70
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:49
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:79
Dune::VTK::DataMode dataMode() const
Definition: io/vtkoutputmodule.hh:188
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:191
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:192
const std::string & name() const
Definition: io/vtkoutputmodule.hh:187
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:194
bool verbose() const
Definition: io/vtkoutputmodule.hh:186
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:300
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:384
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:392
Vtk::template Field< GridView > Field
the type of Field that can be added to this writer
Definition: io/vtkoutputmodule.hh:328
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:386
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:385
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:388
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:389
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
A VTK output module to simplify writing dumux simulation data to VTK format.
constexpr Box box
Definition: method.hh:147
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220