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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef EWOMS_VTK_DISCRETE_FRACTURE_MODULE_HH
28#define EWOMS_VTK_DISCRETE_FRACTURE_MODULE_HH
29
30#include <dune/common/fvector.hh>
31
32#include <opm/material/common/Valgrind.hpp>
33
35
38
41
42#include <cstdio>
43
44namespace Opm::Parameters {
45
46// set default values for what quantities to output
47struct VtkWriteFractureSaturations { static constexpr bool value = true; };
48struct VtkWriteFractureMobilities { static constexpr bool value = false; };
49struct VtkWriteFractureRelativePermeabilities { static constexpr bool value = true; };
50struct VtkWriteFracturePorosity { static constexpr bool value = true; };
51struct VtkWriteFractureIntrinsicPermeabilities { static constexpr bool value = false; };
52struct VtkWriteFractureFilterVelocities { static constexpr bool value = false; };
53struct VtkWriteFractureVolumeFraction { static constexpr bool value = true; };
54
55} // namespace Opm::Parameters
56
57namespace Opm {
71template <class TypeTag>
73{
75
79
83
85
86 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
88
89 enum { dim = GridView::dimension };
90 enum { dimWorld = GridView::dimensionworld };
91 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
92
94 using PhaseBuffer = typename ParentType::PhaseBuffer;
96
97public:
98 VtkDiscreteFractureModule(const Simulator& simulator)
99 : ParentType(simulator)
100 { }
101
105 static void registerParameters()
106 {
107 Parameters::Register<Parameters::VtkWriteFractureSaturations>
108 ("Include the phase saturations in the VTK output files");
109 Parameters::Register<Parameters::VtkWriteFractureMobilities>
110 ("Include the phase mobilities in the VTK output files");
111 Parameters::Register<Parameters::VtkWriteFractureRelativePermeabilities>
112 ("Include the phase relative permeabilities in the "
113 "VTK output files");
114 Parameters::Register<Parameters::VtkWriteFracturePorosity>
115 ("Include the porosity in the VTK output files");
116 Parameters::Register<Parameters::VtkWriteFractureIntrinsicPermeabilities>
117 ("Include the intrinsic permeability in the VTK output files");
118 Parameters::Register<Parameters::VtkWriteFractureFilterVelocities>
119 ("Include in the filter velocities of the phases in the VTK output files");
120 Parameters::Register<Parameters::VtkWriteFractureVolumeFraction>
121 ("Add the fraction of the total volume which is "
122 "occupied by fractures in the VTK output");
123 }
124
130 {
131 if (saturationOutput_())
132 this->resizePhaseBuffer_(fractureSaturation_);
133 if (mobilityOutput_())
134 this->resizePhaseBuffer_(fractureMobility_);
135 if (relativePermeabilityOutput_())
136 this->resizePhaseBuffer_(fractureRelativePermeability_);
137
138 if (porosityOutput_())
139 this->resizeScalarBuffer_(fracturePorosity_);
140 if (intrinsicPermeabilityOutput_())
141 this->resizeScalarBuffer_(fractureIntrinsicPermeability_);
142 if (volumeFractionOutput_())
143 this->resizeScalarBuffer_(fractureVolumeFraction_);
144
145 if (velocityOutput_()) {
146 size_t nDof = this->simulator_.model().numGridDof();
147 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
148 fractureVelocity_[phaseIdx].resize(nDof);
149 for (unsigned dofIdx = 0; dofIdx < nDof; ++dofIdx) {
150 fractureVelocity_[phaseIdx][dofIdx].resize(dimWorld);
151 fractureVelocity_[phaseIdx][dofIdx] = 0.0;
152 }
153 }
154 this->resizePhaseBuffer_(fractureVelocityWeight_);
155 }
156 }
157
162 void processElement(const ElementContext& elemCtx)
163 {
164 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
165 return;
166 }
167
168 const auto& fractureMapper = elemCtx.simulator().vanguard().fractureMapper();
169
170 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
171 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
172 if (!fractureMapper.isFractureVertex(I))
173 continue;
174
175 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
176 const auto& fs = intQuants.fractureFluidState();
177
178 if (porosityOutput_()) {
179 Opm::Valgrind::CheckDefined(intQuants.fracturePorosity());
180 fracturePorosity_[I] = intQuants.fracturePorosity();
181 }
182 if (intrinsicPermeabilityOutput_()) {
183 const auto& K = intQuants.fractureIntrinsicPermeability();
184 fractureIntrinsicPermeability_[I] = K[0][0];
185 }
186
187 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
188 if (saturationOutput_()) {
189 Opm::Valgrind::CheckDefined(fs.saturation(phaseIdx));
190 fractureSaturation_[phaseIdx][I] = fs.saturation(phaseIdx);
191 }
192 if (mobilityOutput_()) {
193 Opm::Valgrind::CheckDefined(intQuants.fractureMobility(phaseIdx));
194 fractureMobility_[phaseIdx][I] = intQuants.fractureMobility(phaseIdx);
195 }
196 if (relativePermeabilityOutput_()) {
197 Opm::Valgrind::CheckDefined(intQuants.fractureRelativePermeability(phaseIdx));
198 fractureRelativePermeability_[phaseIdx][I] =
199 intQuants.fractureRelativePermeability(phaseIdx);
200 }
201 if (volumeFractionOutput_()) {
202 Opm::Valgrind::CheckDefined(intQuants.fractureVolume());
203 fractureVolumeFraction_[I] += intQuants.fractureVolume();
204 }
205 }
206 }
207
208 if (velocityOutput_()) {
209 // calculate velocities if requested by the simulator
210 for (unsigned scvfIdx = 0; scvfIdx < elemCtx.numInteriorFaces(/*timeIdx=*/0); ++ scvfIdx) {
211 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, /*timeIdx=*/0);
212
213 unsigned i = extQuants.interiorIndex();
214 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
215
216 unsigned j = extQuants.exteriorIndex();
217 unsigned J = elemCtx.globalSpaceIndex(j, /*timeIdx=*/0);
218
219 if (!fractureMapper.isFractureEdge(I, J))
220 continue;
221
222 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
223 Scalar weight =
224 std::max<Scalar>(1e-16, std::abs(extQuants.fractureVolumeFlux(phaseIdx)));
225 Opm::Valgrind::CheckDefined(extQuants.extrusionFactor());
226 assert(extQuants.extrusionFactor() > 0);
227 weight *= extQuants.extrusionFactor();
228
229 Dune::FieldVector<Scalar, dim> v(extQuants.fractureFilterVelocity(phaseIdx));
230 v *= weight;
231
232 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
233 fractureVelocity_[phaseIdx][I][dimIdx] += v[dimIdx];
234 fractureVelocity_[phaseIdx][J][dimIdx] += v[dimIdx];
235 }
236
237 fractureVelocityWeight_[phaseIdx][I] += weight;
238 fractureVelocityWeight_[phaseIdx][J] += weight;
239 }
240 }
241 }
242 }
243
248 {
249 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
250 if (!vtkWriter) {
251 return;
252 }
253
254 if (saturationOutput_())
255 this->commitPhaseBuffer_(baseWriter, "fractureSaturation_%s", fractureSaturation_);
256 if (mobilityOutput_())
257 this->commitPhaseBuffer_(baseWriter, "fractureMobility_%s", fractureMobility_);
258 if (relativePermeabilityOutput_())
259 this->commitPhaseBuffer_(baseWriter, "fractureRelativePerm_%s", fractureRelativePermeability_);
260
261 if (porosityOutput_())
262 this->commitScalarBuffer_(baseWriter, "fracturePorosity", fracturePorosity_);
263 if (intrinsicPermeabilityOutput_())
264 this->commitScalarBuffer_(baseWriter, "fractureIntrinsicPerm", fractureIntrinsicPermeability_);
265 if (volumeFractionOutput_()) {
266 // divide the fracture volume by the total volume of the finite volumes
267 for (unsigned I = 0; I < fractureVolumeFraction_.size(); ++I)
268 fractureVolumeFraction_[I] /= this->simulator_.model().dofTotalVolume(I);
269 this->commitScalarBuffer_(baseWriter, "fractureVolumeFraction", fractureVolumeFraction_);
270 }
271
272 if (velocityOutput_()) {
273 size_t nDof = this->simulator_.model().numGridDof();
274
275 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
276 // first, divide the velocity field by the
277 // respective finite volume's surface area
278 for (unsigned dofIdx = 0; dofIdx < nDof; ++dofIdx)
279 fractureVelocity_[phaseIdx][dofIdx] /=
280 std::max<Scalar>(1e-20, fractureVelocityWeight_[phaseIdx][dofIdx]);
281 // commit the phase velocity
282 char name[512];
283 snprintf(name, 512, "fractureFilterVelocity_%s", FluidSystem::phaseName(phaseIdx).data());
284
285 DiscBaseOutputModule::attachVectorDofData_(baseWriter, fractureVelocity_[phaseIdx], name);
286 }
287 }
288 }
289
290private:
291 static bool saturationOutput_()
292 {
293 static bool val = Parameters::Get<Parameters::VtkWriteFractureSaturations>();
294 return val;
295 }
296
297 static bool mobilityOutput_()
298 {
299 static bool val = Parameters::Get<Parameters::VtkWriteFractureMobilities>();
300 return val;
301 }
302
303 static bool relativePermeabilityOutput_()
304 {
305 static bool val = Parameters::Get<Parameters::VtkWriteFractureRelativePermeabilities>();
306 return val;
307 }
308
309 static bool porosityOutput_()
310 {
311 static bool val = Parameters::Get<Parameters::VtkWriteFracturePorosity>();
312 return val;
313 }
314
315 static bool intrinsicPermeabilityOutput_()
316 {
317 static bool val = Parameters::Get<Parameters::VtkWriteFractureIntrinsicPermeabilities>();
318 return val;
319 }
320
321 static bool volumeFractionOutput_()
322 {
323 static bool val = Parameters::Get<Parameters::VtkWriteFractureVolumeFraction>();
324 return val;
325 }
326
327 static bool velocityOutput_()
328 {
329 static bool val = Parameters::Get<Parameters::VtkWriteFractureFilterVelocities>();
330 return val;
331 }
332
333 PhaseBuffer fractureSaturation_;
334 PhaseBuffer fractureMobility_;
335 PhaseBuffer fractureRelativePermeability_;
336
337 ScalarBuffer fracturePorosity_;
338 ScalarBuffer fractureVolumeFraction_;
339 ScalarBuffer fractureIntrinsicPermeability_;
340
341 PhaseVectorBuffer fractureVelocity_;
342 PhaseBuffer fractureVelocityWeight_;
343
344 PhaseVectorBuffer potentialGradient_;
345 PhaseBuffer potentialWeight_;
346};
347
348} // namespace Opm
349
350#endif
The base class for writer modules.
Definition: baseoutputmodule.hh:67
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:85
void commitPhaseBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Add a phase-specific buffer to the result file.
Definition: baseoutputmodule.hh:345
void resizePhaseBuffer_(PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase-specific quantity.
Definition: baseoutputmodule.hh:201
const Simulator & simulator_
Definition: baseoutputmodule.hh:434
std::array< VectorBuffer, numPhases > PhaseVectorBuffer
Definition: baseoutputmodule.hh:94
std::array< ScalarBuffer, numPhases > PhaseBuffer
Definition: baseoutputmodule.hh:90
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:244
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a scalar quantity.
Definition: baseoutputmodule.hh:156
The base class for all output writers.
Definition: baseoutputwriter.hh:44
VTK output module for quantities which make sense for all models which deal with discrete fractures i...
Definition: vtkdiscretefracturemodule.hh:73
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkdiscretefracturemodule.hh:129
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkdiscretefracturemodule.hh:105
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkdiscretefracturemodule.hh:247
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkdiscretefracturemodule.hh:162
VtkDiscreteFractureModule(const Simulator &simulator)
Definition: vtkdiscretefracturemodule.hh:98
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilnewtonmethodparameters.hh:31
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: vtkdiscretefracturemodule.hh:52
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:52
Definition: vtkdiscretefracturemodule.hh:51
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:51
Definition: vtkdiscretefracturemodule.hh:48
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:48
Definition: vtkdiscretefracturemodule.hh:50
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:50
Definition: vtkdiscretefracturemodule.hh:49
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:49
Definition: vtkdiscretefracturemodule.hh:47
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:47
Definition: vtkdiscretefracturemodule.hh:53
static constexpr bool value
Definition: vtkdiscretefracturemodule.hh:53