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;
77 using Field = Vtk::template Field<GridView>;
78 using FractureField = Vtk::template Field<FractureGridView>;
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 ==
DiscretizationMethod::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);
186 for (
const auto& element : elements(gridView, Dune::Partitions::interior))
188 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);
234 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
242 rank[eIdxGlobal] =
static_cast<double>(gridView.comm().rank());
260 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureVertexMapper_, volVarVectorDataFracture[i],
269 "velocity_" + std::string(this->
velocityOutput().phaseName(phaseIdx)) +
" (m/s)",
270 dimWorld, dim).get() );
276 "process rank", 1, 0).get() );
279 for (
auto&& field : this->
fields())
281 if (field.codim() == 0)
283 else if (field.codim() == dim)
286 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
294 fractureSequenceWriter_->write(time, type);
300 fractureWriter_->clear();
304 void writeNonConforming_(
double time, Dune::VTK::OutputType type)
311 std::vector<typename ParentType::VelocityOutput::VelocityVector> velocity;
314 static bool addProcessRank = getParamFromGroup<bool>(this->
paramGroup(),
"Vtk.AddProcessRank");
315 std::vector<double> rank;
318 using ScalarDataContainer = std::vector< std::vector<Scalar> >;
319 using VectorDataContainer = std::vector< std::vector<GlobalPosition> >;
320 std::vector< ScalarDataContainer > volVarScalarData;
321 std::vector< ScalarDataContainer > volVarScalarDataFracture;
322 std::vector< VectorDataContainer > volVarVectorData;
323 std::vector< VectorDataContainer > volVarVectorDataFracture;
327 const auto& fractureGridView = fractureGrid_->leafGridView();
334 || !this->fields().empty()
335 || this->velocityOutput().enableOutput()
338 const auto numCells = gridView.size(0);
339 const auto numDofs = gridView.size(dim);
340 const auto numFractureCells = fractureGridView.size(0);
346 volVarScalarDataFracture.resize(
volVarScalarDataInfo.size(), ScalarDataContainer(numFractureCells));
351 volVarVectorDataFracture.resize(
volVarVectorDataInfo.size(), VectorDataContainer(numFractureCells));
356 velocity[phaseIdx].resize(numDofs);
359 if (addProcessRank) rank.resize(numCells);
361 for (
const auto& element : elements(gridView, Dune::Partitions::interior))
363 const auto eIdxGlobal = this->
gridGeometry().elementMapper().index(element);
371 volVarScalarData[i][eIdxGlobal].resize(
numCorners);
373 volVarVectorData[i][eIdxGlobal].resize(
numCorners);
379 fvGeometry.bind(element);
380 elemVolVars.bind(element, fvGeometry, this->
sol());
384 fvGeometry.bindElement(element);
385 elemVolVars.bindElement(element, fvGeometry, this->
sol());
390 for (
auto&& scv : scvs(fvGeometry))
392 const auto& volVars = elemVolVars[scv];
394 if (!scv.isOnFracture())
403 const auto fIdx = scv.facetIndexInElement();
404 const auto& localMap = fractureElementMap_[eIdxGlobal];
405 const auto fracEIdx = std::find_if(localMap.begin(), localMap.end(), [fIdx] (
const auto& p) { return p.first == fIdx; })->second;
418 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
426 rank[eIdxGlobal] =
static_cast<double>(gridView.comm().rank());
439 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureElementMapper_, volVarScalarDataFracture[i],
449 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureElementMapper_, volVarVectorDataFracture[i],
459 "velocity_" + std::string(this->
velocityOutput().phaseName(phaseIdx)) +
" (m/s)",
460 dimWorld, dim).get() );
466 "process rank", 1, 0).get() );
469 for (
auto&& field : this->
fields())
471 if (field.codim() == 0)
473 else if (field.codim() == dim)
476 DUNE_THROW(Dune::RangeError,
"Cannot add wrongly sized vtk scalar field!");
484 fractureSequenceWriter_->write(time, type);
490 fractureWriter_->clear();
494 template<
class FractureGr
idAdapter >
495 void initializeFracture_(
const FractureGridAdapter& fractureGridAdapter)
499 Dune::GridFactory<FractureGrid> gridFactory;
502 std::size_t fracVertexCount = 0;
503 vertexToFractureVertexIdx_.resize(gridView.size(dim));
504 for (
const auto& v : vertices(gridView))
506 if (fractureGridAdapter.isOnFacetGrid(v))
508 gridFactory.insertVertex(v.geometry().center());
509 vertexToFractureVertexIdx_[
gridGeometry.vertexMapper().index(v)] = fracVertexCount++;
514 std::size_t fractureElementCount = 0;
515 fractureElementMap_.resize(gridView.size(0));
516 std::set< std::pair<GridIndexType, unsigned int> > handledFacets;
517 for (
const auto& element : elements(gridView))
519 const auto eIdxGlobal =
gridGeometry.elementMapper().index(element);
520 const auto refElement = referenceElement(element);
522 for (
const auto& is : intersections(gridView, element))
525 const auto& isGeometry = is.geometry();
527 const auto indexInInside = is.indexInInside();
529 std::vector<GridIndexType> isVertexIndices(
numCorners);
531 isVertexIndices[i] =
gridGeometry.vertexMapper().subIndex(element,
532 refElement.subEntity(indexInInside, 1, i, dim),
536 bool insertFacet =
false;
537 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
543 const auto outsideEIdx =
gridGeometry.elementMapper().index(is.outside());
544 const auto idxInOutside = is.indexInOutside();
545 const auto pair = std::make_pair(outsideEIdx, idxInOutside);
546 if (handledFacets.count( pair ) != 0)
551 const auto& outsideMap = fractureElementMap_[outsideEIdx];
552 auto it = std::find_if(outsideMap.begin(), outsideMap.end(), [idxInOutside] (
const auto& p) { return p.first == idxInOutside; });
553 fractureElementMap_[eIdxGlobal].push_back( std::make_pair(indexInInside, it->second) );
561 std::for_each( isVertexIndices.begin(),
562 isVertexIndices.end(),
563 [&] (
auto& idx) { idx = this->vertexToFractureVertexIdx_[idx]; } );
566 gridFactory.insertElement(isGeometry.type(), isVertexIndices);
569 handledFacets.insert( std::make_pair(eIdxGlobal, indexInInside) );
570 fractureElementMap_[eIdxGlobal].push_back( std::make_pair(indexInInside, fractureElementCount) );
571 fractureElementCount++;
577 fractureGrid_ = std::shared_ptr<FractureGrid>(gridFactory.createGrid());
580 const auto& fractureGridView = fractureGrid_->leafGridView();
581 fractureVertexMapper_ = std::make_unique<FractureMapper>(fractureGridView, Dune::mcmgVertexLayout());
582 fractureElementMapper_ = std::make_unique<FractureMapper>(fractureGridView, Dune::mcmgElementLayout());
585 std::vector<GridIndexType> insToVertexIdx(fractureGridView.size(FractureGridView::dimension));
586 std::vector<GridIndexType> insToElemIdx(fractureGridView.size(0));
587 for (
const auto& v : vertices(fractureGridView)) insToVertexIdx[ gridFactory.insertionIndex(v) ] = fractureVertexMapper_->index(v);
588 for (
const auto& e : elements(fractureGridView)) insToElemIdx[ gridFactory.insertionIndex(e) ] = fractureElementMapper_->index(e);
591 for (GridIndexType dofIdx = 0; dofIdx < gridView.size(GridView::dimension); ++dofIdx)
593 vertexToFractureVertexIdx_[dofIdx] = insToVertexIdx[ vertexToFractureVertexIdx_[dofIdx] ];
596 for (
auto& elemLocalMap : fractureElementMap_)
597 for (
auto& dataPair : elemLocalMap)
598 dataPair.second = insToElemIdx[ dataPair.second ];
601 fractureWriter_ = std::make_shared< Dune::VTKWriter<FractureGridView> >(fractureGridView, this->
dataMode());
602 fractureSequenceWriter_ = std::make_unique< Dune::VTKSequenceWriter<FractureGridView> >(fractureWriter_, this->
name() +
"_fracture");
605 std::shared_ptr<FractureGrid> fractureGrid_;
607 std::unique_ptr<FractureMapper> fractureVertexMapper_;
608 std::unique_ptr<FractureMapper> fractureElementMapper_;
610 std::shared_ptr<Dune::VTKWriter<FractureGridView>> fractureWriter_;
611 std::unique_ptr< Dune::VTKSequenceWriter<FractureGridView> > fractureSequenceWriter_;
614 std::vector<GridIndexType> vertexToFractureVertexIdx_;
617 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
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:215
virtual int numFluidPhases() const
returns the number of phases
Definition: io/velocityoutput.hh:67
virtual void calculateVelocity(VelocityVector &velocity, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVarsCache &elemFluxVarsCache, int phaseIdx) const
Definition: io/velocityoutput.hh:71
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:96
Dune::VTK::DataMode dataMode() const
Definition: io/vtkoutputmodule.hh:196
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:199
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:200
const std::string & name() const
Definition: io/vtkoutputmodule.hh:195
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:202
bool verbose() const
Definition: io/vtkoutputmodule.hh:194
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:307
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:386
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:394
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:388
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:387
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:390
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:391
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.