vtkdiscretefracturemodule.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  Copyright (C) 2011-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef EWOMS_VTK_DISCRETE_FRACTURE_MODULE_HH
26 #define EWOMS_VTK_DISCRETE_FRACTURE_MODULE_HH
27 
28 #include "vtkmultiwriter.hh"
29 #include "baseoutputmodule.hh"
30 
33 
34 #include <dune/common/fvector.hh>
35 
36 #include <cstdio>
37 
38 namespace Ewoms {
39 namespace Properties {
40 // create new type tag for the VTK multi-phase output
41 NEW_TYPE_TAG(VtkDiscreteFracture);
42 
43 // create the property tags needed for the multi phase module
44 NEW_PROP_TAG(GridManager);
45 NEW_PROP_TAG(VtkWriteFractureSaturations);
46 NEW_PROP_TAG(VtkWriteFractureMobilities);
47 NEW_PROP_TAG(VtkWriteFractureRelativePermeabilities);
48 NEW_PROP_TAG(VtkWriteFracturePorosity);
49 NEW_PROP_TAG(VtkWriteFractureIntrinsicPermeabilities);
50 NEW_PROP_TAG(VtkWriteFractureFilterVelocities);
51 NEW_PROP_TAG(VtkWriteFractureVolumeFraction);
52 NEW_PROP_TAG(VtkOutputFormat);
53 NEW_PROP_TAG(DiscBaseOutputModule);
54 
55 // set default values for what quantities to output
56 SET_BOOL_PROP(VtkDiscreteFracture, VtkWriteFractureSaturations, true);
57 SET_BOOL_PROP(VtkDiscreteFracture, VtkWriteFractureMobilities, false);
58 SET_BOOL_PROP(VtkDiscreteFracture, VtkWriteFractureRelativePermeabilities, true);
59 SET_BOOL_PROP(VtkDiscreteFracture, VtkWriteFracturePorosity, true);
60 SET_BOOL_PROP(VtkDiscreteFracture, VtkWriteFractureIntrinsicPermeabilities, false);
61 SET_BOOL_PROP(VtkDiscreteFracture, VtkWriteFractureFilterVelocities, false);
62 SET_BOOL_PROP(VtkDiscreteFracture, VtkWriteFractureVolumeFraction, true);
63 } // namespace Properties
64 
78 template <class TypeTag>
80 {
82 
83  typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
84  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
85  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
86 
87  typedef typename GET_PROP_TYPE(TypeTag, GridManager) GridManager;
88  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
89  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
90 
91  typedef typename GET_PROP_TYPE(TypeTag, DiscBaseOutputModule) DiscBaseOutputModule;
92 
93  static const int vtkFormat = GET_PROP_VALUE(TypeTag, VtkOutputFormat);
95 
96  enum { dim = GridView::dimension };
97  enum { dimWorld = GridView::dimensionworld };
98  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
99 
100  typedef typename ParentType::ScalarBuffer ScalarBuffer;
101  typedef typename ParentType::PhaseBuffer PhaseBuffer;
103 
104 public:
105  VtkDiscreteFractureModule(const Simulator &simulator)
106  : ParentType(simulator)
107  { }
108 
112  static void registerParameters()
113  {
114  EWOMS_REGISTER_PARAM(TypeTag, bool, VtkWriteFractureSaturations,
115  "Include the phase saturations in the VTK output files");
116  EWOMS_REGISTER_PARAM(TypeTag, bool, VtkWriteFractureMobilities,
117  "Include the phase mobilities in the VTK output files");
118  EWOMS_REGISTER_PARAM(TypeTag, bool, VtkWriteFractureRelativePermeabilities,
119  "Include the phase relative permeabilities in the "
120  "VTK output files");
121  EWOMS_REGISTER_PARAM(TypeTag, bool, VtkWriteFracturePorosity,
122  "Include the porosity in the VTK output files");
123  EWOMS_REGISTER_PARAM(TypeTag, bool, VtkWriteFractureIntrinsicPermeabilities,
124  "Include the intrinsic permeability in the VTK output files");
125  EWOMS_REGISTER_PARAM(TypeTag, bool, VtkWriteFractureFilterVelocities,
126  "Include in the filter velocities of the phases "
127  "the VTK output files");
128  EWOMS_REGISTER_PARAM(TypeTag, bool, VtkWriteFractureVolumeFraction,
129  "Add the fraction of the total volume which is "
130  "occupied by fractures in the VTK output");
131  }
132 
138  {
139  if (saturationOutput_())
140  this->resizePhaseBuffer_(fractureSaturation_);
141  if (mobilityOutput_())
142  this->resizePhaseBuffer_(fractureMobility_);
143  if (relativePermeabilityOutput_())
144  this->resizePhaseBuffer_(fractureRelativePermeability_);
145 
146  if (porosityOutput_())
147  this->resizeScalarBuffer_(fracturePorosity_);
148  if (intrinsicPermeabilityOutput_())
149  this->resizeScalarBuffer_(fractureIntrinsicPermeability_);
150  if (volumeFractionOutput_())
151  this->resizeScalarBuffer_(fractureVolumeFraction_);
152 
153  if (velocityOutput_()) {
154  int nDof = this->simulator_.model().numGridDof();
155  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
156  fractureVelocity_[phaseIdx].resize(nDof);
157  for (int dofIdx = 0; dofIdx < nDof; ++dofIdx) {
158  fractureVelocity_[phaseIdx][dofIdx].resize(dimWorld);
159  fractureVelocity_[phaseIdx][dofIdx] = 0.0;
160  }
161  }
162  this->resizePhaseBuffer_(fractureVelocityWeight_);
163  }
164  }
165 
170  void processElement(const ElementContext &elemCtx)
171  {
172  const auto &fractureMapper = elemCtx.simulator().gridManager().fractureMapper();
173 
174  for (int i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
175  int I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
176  if (!fractureMapper.isFractureVertex(I))
177  continue;
178 
179  const auto &intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
180  const auto &fs = intQuants.fractureFluidState();
181 
182  if (porosityOutput_())
183  fracturePorosity_[I] = intQuants.fracturePorosity();
184  if (intrinsicPermeabilityOutput_()) {
185  const auto &K = intQuants.fractureIntrinsicPermeability();
186  fractureIntrinsicPermeability_[I] = K[0][0];
187  }
188 
189  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
190  if (saturationOutput_())
191  fractureSaturation_[phaseIdx][I] = fs.saturation(phaseIdx);
192  if (mobilityOutput_())
193  fractureMobility_[phaseIdx][I] = intQuants.fractureMobility(phaseIdx);
194  if (relativePermeabilityOutput_())
195  fractureRelativePermeability_[phaseIdx][I] =
196  intQuants.fractureRelativePermeability(phaseIdx);
197  if (volumeFractionOutput_())
198  fractureVolumeFraction_[I] += intQuants.fractureVolume();
199  }
200  }
201 
202  if (velocityOutput_()) {
203  // calculate velocities if requested by the simulator
204  for (int scvfIdx = 0; scvfIdx < elemCtx.numInteriorFaces(/*timeIdx=*/0); ++ scvfIdx) {
205  const auto &extQuants = elemCtx.extensiveQuantities(scvfIdx, /*timeIdx=*/0);
206 
207  int i = extQuants.interiorIndex();
208  int I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
209 
210  int j = extQuants.exteriorIndex();
211  int J = elemCtx.globalSpaceIndex(j, /*timeIdx=*/0);
212 
213  if (!fractureMapper.isFractureEdge(I, J))
214  continue;
215 
216  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
217  Scalar weight = std::max(1e-16, std::abs(extQuants.fractureVolumeFlux(phaseIdx)));
218  Valgrind::CheckDefined(extQuants.extrusionFactor());
219  assert(extQuants.extrusionFactor() > 0);
220  weight *= extQuants.extrusionFactor();
221 
222  Dune::FieldVector<Scalar, dim> v(extQuants.fractureFilterVelocity(phaseIdx));
223  v *= weight;
224 
225  for (int dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
226  fractureVelocity_[phaseIdx][I][dimIdx] += v[dimIdx];
227  fractureVelocity_[phaseIdx][J][dimIdx] += v[dimIdx];
228  }
229 
230  fractureVelocityWeight_[phaseIdx][I] += weight;
231  fractureVelocityWeight_[phaseIdx][J] += weight;
232  }
233  }
234  }
235  }
236 
240  void commitBuffers(BaseOutputWriter &baseWriter)
241  {
242  VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
243  if (!vtkWriter) {
244  return;
245  }
246 
247  if (saturationOutput_())
248  this->commitPhaseBuffer_(baseWriter, "fractureSaturation_%s", fractureSaturation_);
249  if (mobilityOutput_())
250  this->commitPhaseBuffer_(baseWriter, "fractureMobility_%s", fractureMobility_);
251  if (relativePermeabilityOutput_())
252  this->commitPhaseBuffer_(baseWriter, "fractureRelativePerm_%s", fractureRelativePermeability_);
253 
254  if (porosityOutput_())
255  this->commitScalarBuffer_(baseWriter, "fracturePorosity", fracturePorosity_);
256  if (intrinsicPermeabilityOutput_())
257  this->commitScalarBuffer_(baseWriter, "fractureIntrinsicPerm", fractureIntrinsicPermeability_);
258  if (volumeFractionOutput_()) {
259  // divide the fracture volume by the total volume of the finite volumes
260  for (unsigned I = 0; I < fractureVolumeFraction_.size(); ++I)
261  fractureVolumeFraction_[I] /= this->simulator_.model().dofTotalVolume(I);
262  this->commitScalarBuffer_(baseWriter, "fractureVolumeFraction", fractureVolumeFraction_);
263  }
264 
265  if (velocityOutput_()) {
266  int nDof = this->simulator_.model().numGridDof();
267 
268  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
269  // first, divide the velocity field by the
270  // respective finite volume's surface area
271  for (int dofIdx = 0; dofIdx < nDof; ++dofIdx)
272  fractureVelocity_[phaseIdx][dofIdx] /=
273  std::max<Scalar>(1e-20, fractureVelocityWeight_[phaseIdx][dofIdx]);
274  // commit the phase velocity
275  char name[512];
276  snprintf(name, 512, "fractureFilterVelocity_%s", FluidSystem::phaseName(phaseIdx));
277 
278  DiscBaseOutputModule::attachVectorDofData_(baseWriter, fractureVelocity_[phaseIdx], name);
279  }
280  }
281  }
282 
283 private:
284  static bool saturationOutput_()
285  { return EWOMS_GET_PARAM(TypeTag, bool, VtkWriteFractureSaturations); }
286 
287  static bool mobilityOutput_()
288  { return EWOMS_GET_PARAM(TypeTag, bool, VtkWriteFractureMobilities); }
289 
290  static bool relativePermeabilityOutput_()
291  { return EWOMS_GET_PARAM(TypeTag, bool, VtkWriteFractureRelativePermeabilities); }
292 
293  static bool porosityOutput_()
294  { return EWOMS_GET_PARAM(TypeTag, bool, VtkWriteFracturePorosity); }
295 
296  static bool intrinsicPermeabilityOutput_()
297  { return EWOMS_GET_PARAM(TypeTag, bool, VtkWriteFractureIntrinsicPermeabilities); }
298 
299  static bool volumeFractionOutput_()
300  { return EWOMS_GET_PARAM(TypeTag, bool, VtkWriteFractureVolumeFraction); }
301 
302  static bool velocityOutput_()
303  { return EWOMS_GET_PARAM(TypeTag, bool, VtkWriteFractureFilterVelocities); }
304 
305  PhaseBuffer fractureSaturation_;
306  PhaseBuffer fractureMobility_;
307  PhaseBuffer fractureRelativePermeability_;
308 
309  ScalarBuffer fracturePorosity_;
310  ScalarBuffer fractureVolumeFraction_;
311  ScalarBuffer fractureIntrinsicPermeability_;
312 
313  PhaseVectorBuffer fractureVelocity_;
314  PhaseBuffer fractureVelocityWeight_;
315 
316  PhaseVectorBuffer potentialGradient_;
317  PhaseBuffer potentialWeight_;
318 };
319 
320 } // namespace Ewoms
321 
322 #endif
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a scalar quantity.
Definition: baseoutputmodule.hh:148
void resizePhaseBuffer_(PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase-specific quantity.
Definition: baseoutputmodule.hh:213
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkdiscretefracturemodule.hh:137
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkdiscretefracturemodule.hh:112
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkdiscretefracturemodule.hh:240
std::array< VectorBuffer, numPhases > PhaseVectorBuffer
Definition: baseoutputmodule.hh:99
NEW_PROP_TAG(Grid)
The type of the DUNE grid.
This file provides the infrastructure to retrieve run-time parameters.
The base class for all output writers.
Definition: baseoutputwriter.hh:41
Simplifies writing multi-file VTK datasets.
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:73
Definition: baseauxiliarymodule.hh:35
#define EWOMS_REGISTER_PARAM(TypeTag, ParamType, ParamName, Description)
Register a run-time parameter.
Definition: parametersystem.hh:64
Model & model()
Return the physical model used in the simulation.
Definition: simulator.hh:176
VtkDiscreteFractureModule(const Simulator &simulator)
Definition: vtkdiscretefracturemodule.hh:105
std::array< ScalarBuffer, numPhases > PhaseBuffer
Definition: baseoutputmodule.hh:95
void commitScalarBuffer_(BaseOutputWriter &baseWriter, const char *name, ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Add a buffer containing scalar quantities to the result file.
Definition: baseoutputmodule.hh:283
VTK output module for quantities which make sense for all models which deal with discrete fractures i...
Definition: vtkdiscretefracturemodule.hh:79
NEW_TYPE_TAG(AuxModule)
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element...
Definition: vtkdiscretefracturemodule.hh:170
Provides the magic behind the eWoms property system.
The base class for writer modules.
void commitPhaseBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Add a phase-specific buffer to the result file.
Definition: baseoutputmodule.hh:386
SET_BOOL_PROP(FvBaseDiscretization, EnableVtkOutput, true)
Enable the VTK output by default.
The base class for writer modules.
Definition: baseoutputmodule.hh:71
const Simulator & simulator_
Definition: baseoutputmodule.hh:487
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:61
#define EWOMS_GET_PARAM(TypeTag, ParamType, ParamName)
Retrieve a runtime parameter.
Definition: parametersystem.hh:95
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:90