vtkdiffusionmodule.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 */
26 #ifndef EWOMS_VTK_DIFFUSION_MODULE_HH
27 #define EWOMS_VTK_DIFFUSION_MODULE_HH
28 
29 #include "vtkmultiwriter.hh"
30 #include "baseoutputmodule.hh"
31 
34 
35 namespace Ewoms {
36 namespace Properties {
37 // create new type tag for the VTK output of the quantities for molecular
38 // diffusion
39 NEW_TYPE_TAG(VtkDiffusion);
40 
41 // create the property tags needed for the diffusion module
42 NEW_PROP_TAG(VtkWriteTortuosities);
43 NEW_PROP_TAG(VtkWriteDiffusionCoefficients);
44 NEW_PROP_TAG(VtkWriteEffectiveDiffusionCoefficients);
45 NEW_PROP_TAG(VtkOutputFormat);
46 
47 // set default values for what quantities to output
48 SET_BOOL_PROP(VtkDiffusion, VtkWriteTortuosities, false);
49 SET_BOOL_PROP(VtkDiffusion, VtkWriteDiffusionCoefficients, false);
50 SET_BOOL_PROP(VtkDiffusion, VtkWriteEffectiveDiffusionCoefficients, false);
51 } // namespace Properties
52 } // namespace Ewoms
53 
54 namespace Ewoms {
66 template <class TypeTag>
67 class VtkDiffusionModule : public BaseOutputModule<TypeTag>
68 {
70 
71  typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
72  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
73 
74  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
75 
77  typedef typename ParentType::PhaseBuffer PhaseBuffer;
78 
79  static const int vtkFormat = GET_PROP_VALUE(TypeTag, VtkOutputFormat);
81 
82  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
83  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
84 
85 public:
86  VtkDiffusionModule(const Simulator &simulator)
87  : ParentType(simulator)
88  { }
89 
93  static void registerParameters()
94  {
95  EWOMS_REGISTER_PARAM(TypeTag, bool, VtkWriteTortuosities,
96  "Include the tortuosity for each phase in the VTK "
97  "output files");
98  EWOMS_REGISTER_PARAM(TypeTag, bool, VtkWriteDiffusionCoefficients,
99  "Include the molecular diffusion coefficients in "
100  "the VTK output files");
101  EWOMS_REGISTER_PARAM(TypeTag, bool,
102  VtkWriteEffectiveDiffusionCoefficients,
103  "Include the effective molecular diffusion "
104  "coefficients the medium in the VTK output files");
105  }
106 
112  {
113  if (tortuosityOutput_())
114  this->resizePhaseBuffer_(tortuosity_);
115  if (diffusionCoefficientOutput_())
116  this->resizePhaseComponentBuffer_(diffusionCoefficient_);
117  if (effectiveDiffusionCoefficientOutput_())
118  this->resizePhaseComponentBuffer_(effectiveDiffusionCoefficient_);
119  }
120 
125  void processElement(const ElementContext &elemCtx)
126  {
127  for (int i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
128  int I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
129  const auto &intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
130 
131  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
132  if (tortuosityOutput_())
133  tortuosity_[phaseIdx][I] = intQuants.tortuosity(phaseIdx);
134  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
135  if (diffusionCoefficientOutput_())
136  diffusionCoefficient_[phaseIdx][compIdx][I] =
137  intQuants.diffusionCoefficient(phaseIdx, compIdx);
138  if (effectiveDiffusionCoefficientOutput_())
139  effectiveDiffusionCoefficient_[phaseIdx][compIdx][I] =
140  intQuants.effectiveDiffusionCoefficient(phaseIdx, compIdx);
141  }
142  }
143  }
144  }
145 
149  void commitBuffers(BaseOutputWriter &baseWriter)
150  {
151  VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
152  if (!vtkWriter) {
153  return;
154  }
155 
156  if (tortuosityOutput_())
157  this->commitPhaseBuffer_(baseWriter, "tortuosity", tortuosity_);
158  if (diffusionCoefficientOutput_())
159  this->commitPhaseComponentBuffer_(baseWriter, "diffusionCoefficient",
160  diffusionCoefficient_);
161  if (effectiveDiffusionCoefficientOutput_())
162  this->commitPhaseComponentBuffer_(baseWriter,
163  "effectiveDiffusionCoefficient",
164  effectiveDiffusionCoefficient_);
165  }
166 
167 private:
168  static bool tortuosityOutput_()
169  { return EWOMS_GET_PARAM(TypeTag, bool, VtkWriteTortuosities); }
170 
171  static bool diffusionCoefficientOutput_()
172  { return EWOMS_GET_PARAM(TypeTag, bool, VtkWriteDiffusionCoefficients); }
173 
174  static bool effectiveDiffusionCoefficientOutput_()
175  {
176  return EWOMS_GET_PARAM(TypeTag, bool,
177  VtkWriteEffectiveDiffusionCoefficients);
178  }
179 
180  PhaseBuffer tortuosity_;
181  PhaseComponentBuffer diffusionCoefficient_;
182  PhaseComponentBuffer effectiveDiffusionCoefficient_;
183 };
184 
185 } // namespace Ewoms
186 
187 #endif
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quanties relevant for an element.
Definition: vtkdiffusionmodule.hh:125
VtkDiffusionModule(const Simulator &simulator)
Definition: vtkdiffusionmodule.hh:86
void resizePhaseBuffer_(PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase-specific quantity.
Definition: baseoutputmodule.hh:213
void resizePhaseComponentBuffer_(PhaseComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase and component specific buffer.
Definition: baseoutputmodule.hh:259
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
std::array< std::array< ScalarBuffer, numComponents >, numPhases > PhaseComponentBuffer
Definition: baseoutputmodule.hh:97
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hh:93
void commitPhaseComponentBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Add a phase and component specific quantities to the output.
Definition: baseoutputmodule.hh:432
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
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkdiffusionmodule.hh:111
std::array< ScalarBuffer, numPhases > PhaseBuffer
Definition: baseoutputmodule.hh:95
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkdiffusionmodule.hh:149
NEW_TYPE_TAG(AuxModule)
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.
VTK output module for quantities which make sense for models which incorperate molecular diffusion...
Definition: vtkdiffusionmodule.hh:67
The base class for writer modules.
Definition: baseoutputmodule.hh:71
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