3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/boxdfm/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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef POROUSMEDIUMFLOW_BOXDFM_VTK_OUTPUT_MODULE_HH
26#define POROUSMEDIUMFLOW_BOXDFM_VTK_OUTPUT_MODULE_HH
27
28#include <set>
29
30#include <dune/grid/common/gridfactory.hh>
31#include <dune/grid/common/mcmgmapper.hh>
32
34
35namespace Dumux {
36
54template<class GridVariables, class SolutionVector, class FractureGrid>
55class BoxDfmVtkOutputModule : public VtkOutputModule<GridVariables, SolutionVector>
56{
58 using GridGeometry = typename GridVariables::GridGeometry;
59 using VV = typename GridVariables::VolumeVariables;
60 using FluidSystem = typename VV::FluidSystem;
61 using Scalar = typename GridVariables::Scalar;
62
63 using GridView = typename GridGeometry::GridView;
64 using FractureGridView = typename FractureGrid::LeafGridView;
65 using FractureMapper = Dune::MultipleCodimMultipleGeomTypeMapper<FractureGridView>;
66
67 enum
68 {
69 dim = GridView::dimension,
70 dimWorld = GridView::dimensionworld
71 };
72
73 using GridIndexType = typename GridView::IndexSet::IndexType;
74 using Element = typename GridView::template Codim<0>::Entity;
75 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
76
77 using Field = Vtk::template Field<GridView>;
78 using FractureField = Vtk::template Field<FractureGridView>;
79
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!");
84public:
85
87 template< class FractureGridAdapter >
89 const SolutionVector& sol,
90 const std::string& name,
91 const FractureGridAdapter& fractureGridAdapter,
92 const std::string& paramGroup = "",
93 Dune::VTK::DataMode dm = Dune::VTK::conforming,
94 bool verbose = true)
96 {
97 // create the fracture grid and all objects needed on it
98 initializeFracture_(fractureGridAdapter);
99 }
100
104
110 void write(double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
111 {
112 Dune::Timer timer;
113
114 // write to file depending on data mode
115 const auto dm = this->dataMode();
116 if (dm == Dune::VTK::conforming)
117 writeConforming_(time, type);
118 else if (dm == Dune::VTK::nonconforming)
119 writeNonConforming_(time, type);
120 else
121 DUNE_THROW(Dune::NotImplemented, "Output for provided vtk data mode");
122
124 timer.stop();
125 if (this->verbose())
126 std::cout << "Writing output for problem \"" << this->name() << "\". Took " << timer.elapsed() << " seconds." << std::endl;
127 }
128
129private:
131 void writeConforming_(double time, Dune::VTK::OutputType type)
132 {
136
137 // instatiate the velocity output
138 std::vector<typename ParentType::VelocityOutput::VelocityVector> velocity;
139
140 // process rank
141 static bool addProcessRank = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddProcessRank");
142 std::vector<double> rank;
143
144 // volume variable data
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;
149
150 // some references for convenience
151 const auto& gridView = this->gridGeometry().gridView();
152 const auto& fractureGridView = fractureGrid_->leafGridView();
153 const auto& volVarScalarDataInfo = this->volVarScalarDataInfo();
154 const auto& volVarVectorDataInfo = this->volVarVectorDataInfo();
155
157 if (!volVarScalarDataInfo.empty()
158 || !volVarVectorDataInfo.empty()
159 || !this->fields().empty()
160 || this->velocityOutput().enableOutput()
161 || addProcessRank)
162 {
163 const auto numCells = gridView.size(0);
164 const auto numDofs = gridView.size(dim);
165 const auto numFractureVert = fractureGridView.size(FractureGridView::dimension);
166
167 // get fields for all volume variables
168 if (!this->volVarScalarDataInfo().empty())
169 {
170 volVarScalarData.resize(volVarScalarDataInfo.size(), std::vector<Scalar>(numDofs));
171 volVarScalarDataFracture.resize(volVarScalarDataInfo.size(), std::vector<Scalar>(numFractureVert));
172 }
173 if (!this->volVarVectorDataInfo().empty())
174 {
175 volVarVectorData.resize(volVarVectorDataInfo.size(), std::vector<GlobalPosition>(numDofs));
176 volVarVectorDataFracture.resize(volVarVectorDataInfo.size(), std::vector<GlobalPosition>(numFractureVert));
177 }
178
179 if (this->velocityOutput().enableOutput())
180 for (int phaseIdx = 0; phaseIdx < this->velocityOutput().numFluidPhases(); ++phaseIdx)
181 velocity[phaseIdx].resize(numDofs);
182
183 // maybe allocate space for the process rank
184 if (addProcessRank) rank.resize(numCells);
185
186 auto fvGeometry = localView(this->gridGeometry());
187 auto elemVolVars = localView(this->gridVariables().curGridVolVars());
188 auto elemFluxVarsCache = localView(this->gridVariables().gridFluxVarsCache());
189 for (const auto& element : elements(gridView, Dune::Partitions::interior))
190 {
191 const auto eIdxGlobal = this->gridGeometry().elementMapper().index(element);
192
193 // If velocity output is enabled we need to bind to the whole stencil
194 // otherwise element-local data is sufficient
195 if (this->velocityOutput().enableOutput())
196 {
197 fvGeometry.bind(element);
198 elemVolVars.bind(element, fvGeometry, this->sol());
199 }
200 else
201 {
202 fvGeometry.bindElement(element);
203 elemVolVars.bindElement(element, fvGeometry, this->sol());
204 }
205
206 if (!volVarScalarDataInfo.empty() || !volVarVectorDataInfo.empty())
207 {
208 for (auto&& scv : scvs(fvGeometry))
209 {
210 const auto dofIdxGlobal = scv.dofIndex();
211 const auto& volVars = elemVolVars[scv];
212
213 if (!scv.isOnFracture())
214 {
215 for (std::size_t i = 0; i < volVarScalarDataInfo.size(); ++i)
216 volVarScalarData[i][dofIdxGlobal] = volVarScalarDataInfo[i].get(volVars);
217 for (std::size_t i = 0; i < volVarVectorDataInfo.size(); ++i)
218 volVarVectorData[i][dofIdxGlobal] = volVarVectorDataInfo[i].get(volVars);
219 }
220 else
221 {
222 for (std::size_t i = 0; i < volVarScalarDataInfo.size(); ++i)
223 volVarScalarDataFracture[i][vertexToFractureVertexIdx_[dofIdxGlobal]] = volVarScalarDataInfo[i].get(volVars);
224 for (std::size_t i = 0; i < volVarVectorDataInfo.size(); ++i)
225 volVarVectorDataFracture[i][vertexToFractureVertexIdx_[dofIdxGlobal]] = volVarVectorDataInfo[i].get(volVars);
226 }
227 }
228 }
229
230 // velocity output
231 if (this->velocityOutput().enableOutput())
232 {
233 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
234
235 for (int phaseIdx = 0; phaseIdx < this->velocityOutput().numFluidPhases(); ++phaseIdx)
236 this->velocityOutput().calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
237 }
238
240 if (addProcessRank)
241 rank[eIdxGlobal] = static_cast<double>(gridView.comm().rank());
242 }
243
247
248 // volume variables if any
249 for (std::size_t i = 0; i < volVarScalarDataInfo.size(); ++i)
250 {
251 this->sequenceWriter().addVertexData(volVarScalarData[i], volVarScalarDataInfo[i].name);
252 fractureSequenceWriter_->addVertexData(volVarScalarDataFracture[i], volVarScalarDataInfo[i].name);
253 }
254
255 for (std::size_t i = 0; i < volVarVectorDataInfo.size(); ++i)
256 {
257 this->sequenceWriter().addVertexData( Field(gridView, this->gridGeometry().vertexMapper(), volVarVectorData[i],
258 volVarVectorDataInfo[i].name, /*numComp*/dimWorld, /*codim*/dim).get() );
259 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureVertexMapper_, volVarVectorDataFracture[i],
260 volVarVectorDataInfo[i].name, /*numComp*/dimWorld, /*codim*/dim-1).get() );
261 }
262
263 // the velocity field
264 if (this->velocityOutput().enableOutput())
265 {
266 for (int phaseIdx = 0; phaseIdx < this->velocityOutput().numFluidPhases(); ++phaseIdx)
267 this->sequenceWriter().addVertexData( Field(gridView, this->gridGeometry().vertexMapper(), velocity[phaseIdx],
268 "velocity_" + std::string(this->velocityOutput().phaseName(phaseIdx)) + " (m/s)",
269 /*numComp*/dimWorld, /*codim*/dim).get() );
270 }
271
272 // the process rank
273 if (addProcessRank)
274 this->sequenceWriter().addCellData( Field(gridView, this->gridGeometry().elementMapper(), rank,
275 "process rank", /*numComp*/1, /*codim*/0).get() );
276
277 // also register additional (non-standardized) user fields if any (only on matrix grid)
278 for (auto&& field : this->fields())
279 {
280 if (field.codim() == 0)
281 this->sequenceWriter().addCellData(field.get());
282 else if (field.codim() == dim)
283 this->sequenceWriter().addVertexData(field.get());
284 else
285 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
286 }
287 }
288
292 this->sequenceWriter().write(time, type);
293 fractureSequenceWriter_->write(time, type);
294
298 this->writer().clear();
299 fractureWriter_->clear();
300 }
301
303 void writeNonConforming_(double time, Dune::VTK::OutputType type)
304 {
308
309 // instatiate the velocity output
310 std::vector<typename ParentType::VelocityOutput::VelocityVector> velocity;
311
312 // process rank
313 static bool addProcessRank = getParamFromGroup<bool>(this->paramGroup(), "Vtk.AddProcessRank");
314 std::vector<double> rank;
315
316 // volume variable data (indexing: volvardata/element/localcorner)
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;
323
324 // some references for convenience
325 const auto& gridView = this->gridGeometry().gridView();
326 const auto& fractureGridView = fractureGrid_->leafGridView();
327 const auto& volVarScalarDataInfo = this->volVarScalarDataInfo();
328 const auto& volVarVectorDataInfo = this->volVarVectorDataInfo();
329
331 if (!volVarScalarDataInfo.empty()
332 || !volVarVectorDataInfo.empty()
333 || !this->fields().empty()
334 || this->velocityOutput().enableOutput()
335 || addProcessRank)
336 {
337 const auto numCells = gridView.size(0);
338 const auto numDofs = gridView.size(dim);
339 const auto numFractureCells = fractureGridView.size(0);
340
341 // get fields for all volume variables
342 if (!this->volVarScalarDataInfo().empty())
343 {
344 volVarScalarData.resize(volVarScalarDataInfo.size(), ScalarDataContainer(numCells));
345 volVarScalarDataFracture.resize(volVarScalarDataInfo.size(), ScalarDataContainer(numFractureCells));
346 }
347 if (!this->volVarVectorDataInfo().empty())
348 {
349 volVarVectorData.resize(volVarVectorDataInfo.size(), VectorDataContainer(numCells));
350 volVarVectorDataFracture.resize(volVarVectorDataInfo.size(), VectorDataContainer(numFractureCells));
351 }
352
353 if (this->velocityOutput().enableOutput())
354 for (int phaseIdx = 0; phaseIdx < this->velocityOutput().numFluidPhases(); ++phaseIdx)
355 velocity[phaseIdx].resize(numDofs);
356
357 // maybe allocate space for the process rank
358 if (addProcessRank) rank.resize(numCells);
359
360 for (const auto& element : elements(gridView, Dune::Partitions::interior))
361 {
362 const auto eIdxGlobal = this->gridGeometry().elementMapper().index(element);
363 const auto numCorners = element.subEntities(dim);
364
365 auto fvGeometry = localView(this->gridGeometry());
366 auto elemVolVars = localView(this->gridVariables().curGridVolVars());
367
368 // resize element-local data containers (for bulk grid)
369 for (std::size_t i = 0; i < volVarScalarDataInfo.size(); ++i)
370 volVarScalarData[i][eIdxGlobal].resize(numCorners);
371 for (std::size_t i = 0; i < volVarVectorDataInfo.size(); ++i)
372 volVarVectorData[i][eIdxGlobal].resize(numCorners);
373
374 // If velocity output is enabled we need to bind to the whole stencil
375 // otherwise element-local data is sufficient
376 if (this->velocityOutput().enableOutput())
377 {
378 fvGeometry.bind(element);
379 elemVolVars.bind(element, fvGeometry, this->sol());
380 }
381 else
382 {
383 fvGeometry.bindElement(element);
384 elemVolVars.bindElement(element, fvGeometry, this->sol());
385 }
386
387 if (!volVarScalarDataInfo.empty() || !volVarVectorDataInfo.empty())
388 {
389 for (auto&& scv : scvs(fvGeometry))
390 {
391 const auto& volVars = elemVolVars[scv];
392
393 if (!scv.isOnFracture())
394 {
395 for (std::size_t i = 0; i < volVarScalarDataInfo.size(); ++i)
396 volVarScalarData[i][eIdxGlobal][scv.localDofIndex()] = volVarScalarDataInfo[i].get(volVars);
397 for (std::size_t i = 0; i < volVarVectorDataInfo.size(); ++i)
398 volVarVectorData[i][eIdxGlobal][scv.localDofIndex()] = volVarVectorDataInfo[i].get(volVars);
399 }
400 else
401 {
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;
405 for (std::size_t i = 0; i < volVarScalarDataInfo.size(); ++i)
406 volVarScalarDataFracture[i][fracEIdx].push_back(volVarScalarDataInfo[i].get(volVars));
407 for (std::size_t i = 0; i < volVarVectorDataInfo.size(); ++i)
408 volVarVectorDataFracture[i][fracEIdx].push_back(volVarVectorDataInfo[i].get(volVars));
409 }
410 }
411 }
412
413 // velocity output
414 if (this->velocityOutput().enableOutput())
415 {
416 const auto elemFluxVarsCache = localView(this->gridVariables().gridFluxVarsCache()).bind(element, fvGeometry, elemVolVars);
417
418 for (int phaseIdx = 0; phaseIdx < this->velocityOutput().numFluidPhases(); ++phaseIdx)
419 this->velocityOutput().calculateVelocity(velocity[phaseIdx], element, fvGeometry, elemVolVars, elemFluxVarsCache, phaseIdx);
420 }
421
423 if (addProcessRank)
424 rank[eIdxGlobal] = static_cast<double>(gridView.comm().rank());
425 }
426
430
431 // volume variables if any
432 for (std::size_t i = 0; i < volVarScalarDataInfo.size(); ++i)
433 {
434 this->sequenceWriter().addVertexData( Field(gridView, this->gridGeometry().elementMapper(), volVarScalarData[i],
435 volVarScalarDataInfo[i].name, /*numComp*/1, /*codim*/dim,
436 /*nonconforming*/this->dataMode()).get() );
437 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureElementMapper_, volVarScalarDataFracture[i],
438 volVarScalarDataInfo[i].name, /*numComp*/1, /*codim*/dim-1,
439 /*nonconforming*/this->dataMode()).get() );
440 }
441
442 for (std::size_t i = 0; i < volVarVectorDataInfo.size(); ++i)
443 {
444 this->sequenceWriter().addVertexData( Field(gridView, this->gridGeometry().elementMapper(), volVarVectorData[i],
445 volVarVectorDataInfo[i].name, /*numComp*/dimWorld, /*codim*/dim,
446 /*nonconforming*/this->dataMode()).get() );
447 fractureSequenceWriter_->addVertexData( FractureField(fractureGridView, *fractureElementMapper_, volVarVectorDataFracture[i],
448 volVarVectorDataInfo[i].name, /*numComp*/dimWorld, /*codim*/dim-1,
449 /*nonconforming*/this->dataMode()).get() );
450 }
451
452 // the velocity field
453 if (this->velocityOutput().enableOutput())
454 {
455 for (int phaseIdx = 0; phaseIdx < this->velocityOutput().numFluidPhases(); ++phaseIdx)
456 this->sequenceWriter().addVertexData( Field(gridView, this->gridGeometry().vertexMapper(), velocity[phaseIdx],
457 "velocity_" + std::string(this->velocityOutput().phaseName(phaseIdx)) + " (m/s)",
458 /*numComp*/dimWorld, /*codim*/dim).get() );
459 }
460
461 // the process rank
462 if (addProcessRank)
463 this->sequenceWriter().addCellData( Field(gridView, this->gridGeometry().elementMapper(), rank,
464 "process rank", /*numComp*/1, /*codim*/0).get() );
465
466 // also register additional (non-standardized) user fields if any (only on matrix grid)
467 for (auto&& field : this->fields())
468 {
469 if (field.codim() == 0)
470 this->sequenceWriter().addCellData(field.get());
471 else if (field.codim() == dim)
472 this->sequenceWriter().addVertexData(field.get());
473 else
474 DUNE_THROW(Dune::RangeError, "Cannot add wrongly sized vtk scalar field!");
475 }
476 }
477
481 this->sequenceWriter().write(time, type);
482 fractureSequenceWriter_->write(time, type);
483
487 this->writer().clear();
488 fractureWriter_->clear();
489 }
490
492 template< class FractureGridAdapter >
493 void initializeFracture_(const FractureGridAdapter& fractureGridAdapter)
494 {
495 const auto& gridGeometry = this->gridGeometry();
496 const auto& gridView = gridGeometry.gridView();
497 Dune::GridFactory<FractureGrid> gridFactory;
498
499 // insert fracture vertices
500 std::size_t fracVertexCount = 0;
501 vertexToFractureVertexIdx_.resize(gridView.size(dim));
502 for (const auto& v : vertices(gridView))
503 {
504 if (fractureGridAdapter.isOnFacetGrid(v))
505 {
506 gridFactory.insertVertex(v.geometry().center());
507 vertexToFractureVertexIdx_[gridGeometry.vertexMapper().index(v)] = fracVertexCount++;
508 }
509 }
510
511 // insert fracture elements
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))
516 {
517 const auto eIdxGlobal = gridGeometry.elementMapper().index(element);
518 const auto refElement = referenceElement(element);
519
520 for (const auto& is : intersections(gridView, element))
521 {
522 // obtain all vertex indices on this intersection
523 const auto& isGeometry = is.geometry();
524 const auto numCorners = isGeometry.corners();
525 const auto indexInInside = is.indexInInside();
526
527 std::vector<GridIndexType> isVertexIndices(numCorners);
528 for (unsigned int i = 0; i < numCorners; ++i)
529 isVertexIndices[i] = gridGeometry.vertexMapper().subIndex(element,
530 refElement.subEntity(indexInInside, 1, i, dim),
531 dim);
532
533 // determine if this is a fracture facet & if it has to be inserted
534 bool insertFacet = false;
535 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
536 {
537 insertFacet = true;
538 if (!is.boundary())
539 {
540 // only proceed if facet has not been handled yet
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)
545 {
546 insertFacet = false;
547
548 // obtain the fracture grid elem idx from outside map and insert to map
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) );
552 }
553 }
554 }
555
556 if (insertFacet)
557 {
558 // transform intersection vertex indices to frac grid indices
559 std::for_each( isVertexIndices.begin(),
560 isVertexIndices.end(),
561 [&] (auto& idx) { idx = this->vertexToFractureVertexIdx_[idx]; } );
562
563 // insert the element
564 gridFactory.insertElement(isGeometry.type(), isVertexIndices);
565
566 // insert to set of handled facets
567 handledFacets.insert( std::make_pair(eIdxGlobal, indexInInside) );
568 fractureElementMap_[eIdxGlobal].push_back( std::make_pair(indexInInside, fractureElementCount) );
569 fractureElementCount++;
570 }
571 }
572 }
573
574 // make grid and get grid view
575 fractureGrid_ = std::shared_ptr<FractureGrid>(gridFactory.createGrid());
576
577 // update fracture mappers
578 const auto& fractureGridView = fractureGrid_->leafGridView();
579 fractureVertexMapper_ = std::make_unique<FractureMapper>(fractureGridView, Dune::mcmgVertexLayout());
580 fractureElementMapper_ = std::make_unique<FractureMapper>(fractureGridView, Dune::mcmgElementLayout());
581
582 // obtain map fracture insertion indices -> fracture grid indices
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);
587
588 // update vertex index map
589 for (GridIndexType dofIdx = 0; dofIdx < gridView.size(GridView::dimension); ++dofIdx)
590 if (gridGeometry.dofOnFracture(dofIdx))
591 vertexToFractureVertexIdx_[dofIdx] = insToVertexIdx[ vertexToFractureVertexIdx_[dofIdx] ];
592
593 // update fracture element map
594 for (auto& elemLocalMap : fractureElementMap_)
595 for (auto& dataPair : elemLocalMap)
596 dataPair.second = insToElemIdx[ dataPair.second ];
597
598 // instantiate writers for the fracture
599 fractureWriter_ = std::make_shared< Dune::VTKWriter<FractureGridView> >(fractureGridView, this->dataMode());
600 fractureSequenceWriter_ = std::make_unique< Dune::VTKSequenceWriter<FractureGridView> >(fractureWriter_, this->name() + "_fracture");
601 }
602
603 std::shared_ptr<FractureGrid> fractureGrid_;
604
605 std::unique_ptr<FractureMapper> fractureVertexMapper_;
606 std::unique_ptr<FractureMapper> fractureElementMapper_;
607
608 std::shared_ptr<Dune::VTKWriter<FractureGridView>> fractureWriter_;
609 std::unique_ptr< Dune::VTKSequenceWriter<FractureGridView> > fractureSequenceWriter_;
610
611 // maps to a bulk grid vertex the vertex index within the fracture grid
612 std::vector<GridIndexType> vertexToFractureVertexIdx_;
613
614 // maps to the local facet indices of an element the corresponding fracture element indices
615 std::vector< std::vector<std::pair<GridIndexType, unsigned int>> > fractureElementMap_;
616};
617
618} // end namespace Dumux
619
620#endif
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Definition: adapt.hh:29
constexpr Box box
Definition: method.hh:139
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:97
Dune::VTK::DataMode dataMode() const
Definition: io/vtkoutputmodule.hh:197
Dune::VTKWriter< GridView > & writer()
Definition: io/vtkoutputmodule.hh:200
Dune::VTKSequenceWriter< GridView > & sequenceWriter()
Definition: io/vtkoutputmodule.hh:201
const std::string & name() const
Definition: io/vtkoutputmodule.hh:196
const std::vector< Field > & fields() const
Definition: io/vtkoutputmodule.hh:203
bool verbose() const
Definition: io/vtkoutputmodule.hh:195
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:309
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:391
const VelocityOutput & velocityOutput() const
Definition: io/vtkoutputmodule.hh:399
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:393
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:392
const std::vector< VolVarScalarDataInfo > & volVarScalarDataInfo() const
Definition: io/vtkoutputmodule.hh:395
const std::vector< VolVarVectorDataInfo > & volVarVectorDataInfo() const
Definition: io/vtkoutputmodule.hh:396
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 &paramGroup="", 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.