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 "vtkmultiwriter.hh"
31#include "baseoutputmodule.hh"
32
35
36#include <opm/material/common/Valgrind.hpp>
37
38#include <dune/common/fvector.hh>
39
40#include <cstdio>
41#include <string_view>
42
43namespace Opm::Properties {
44
45namespace TTag {
46
47// create new type tag for the VTK multi-phase output
49
50} // namespace TTag
51
52// create the property tags needed for the multi phase module
53template<class TypeTag, class MyTypeTag>
55template<class TypeTag, class MyTypeTag>
57template<class TypeTag, class MyTypeTag>
59template<class TypeTag, class MyTypeTag>
61template<class TypeTag, class MyTypeTag>
63template<class TypeTag, class MyTypeTag>
65template<class TypeTag, class MyTypeTag>
67
68// set default values for what quantities to output
69template<class TypeTag>
70struct VtkWriteFractureSaturations<TypeTag, TTag::VtkDiscreteFracture> { static constexpr bool value = true; };
71template<class TypeTag>
72struct VtkWriteFractureMobilities<TypeTag, TTag::VtkDiscreteFracture> { static constexpr bool value = false; };
73template<class TypeTag>
74struct VtkWriteFractureRelativePermeabilities<TypeTag, TTag::VtkDiscreteFracture> { static constexpr bool value = true; };
75template<class TypeTag>
76struct VtkWriteFracturePorosity<TypeTag, TTag::VtkDiscreteFracture> { static constexpr bool value = true; };
77template<class TypeTag>
78struct VtkWriteFractureIntrinsicPermeabilities<TypeTag, TTag::VtkDiscreteFracture> { static constexpr bool value = false; };
79template<class TypeTag>
80struct VtkWriteFractureFilterVelocities<TypeTag, TTag::VtkDiscreteFracture> { static constexpr bool value = false; };
81template<class TypeTag>
82struct VtkWriteFractureVolumeFraction<TypeTag, TTag::VtkDiscreteFracture> { static constexpr bool value = true; };
83
84} // namespace Opm::Properties
85
86namespace Opm {
100template <class TypeTag>
102{
104
108
112
114
115 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
117
118 enum { dim = GridView::dimension };
119 enum { dimWorld = GridView::dimensionworld };
120 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
121
122 using ScalarBuffer = typename ParentType::ScalarBuffer;
123 using PhaseBuffer = typename ParentType::PhaseBuffer;
125
126public:
127 VtkDiscreteFractureModule(const Simulator& simulator)
128 : ParentType(simulator)
129 { }
130
134 static void registerParameters()
135 {
136 Parameters::registerParam<TypeTag, Properties::VtkWriteFractureSaturations>
137 ("Include the phase saturations in the VTK output files");
138 Parameters::registerParam<TypeTag, Properties::VtkWriteFractureMobilities>
139 ("Include the phase mobilities in the VTK output files");
140 Parameters::registerParam<TypeTag, Properties::VtkWriteFractureRelativePermeabilities>
141 ("Include the phase relative permeabilities in the "
142 "VTK output files");
143 Parameters::registerParam<TypeTag, Properties::VtkWriteFracturePorosity>
144 ("Include the porosity in the VTK output files");
145 Parameters::registerParam<TypeTag, Properties::VtkWriteFractureIntrinsicPermeabilities>
146 ("Include the intrinsic permeability in the VTK output files");
147 Parameters::registerParam<TypeTag, Properties::VtkWriteFractureFilterVelocities>
148 ("Include in the filter velocities of the phases in the VTK output files");
149 Parameters::registerParam<TypeTag, Properties::VtkWriteFractureVolumeFraction>
150 ("Add the fraction of the total volume which is "
151 "occupied by fractures in the VTK output");
152 }
153
159 {
160 if (saturationOutput_())
161 this->resizePhaseBuffer_(fractureSaturation_);
162 if (mobilityOutput_())
163 this->resizePhaseBuffer_(fractureMobility_);
164 if (relativePermeabilityOutput_())
165 this->resizePhaseBuffer_(fractureRelativePermeability_);
166
167 if (porosityOutput_())
168 this->resizeScalarBuffer_(fracturePorosity_);
169 if (intrinsicPermeabilityOutput_())
170 this->resizeScalarBuffer_(fractureIntrinsicPermeability_);
171 if (volumeFractionOutput_())
172 this->resizeScalarBuffer_(fractureVolumeFraction_);
173
174 if (velocityOutput_()) {
175 size_t nDof = this->simulator_.model().numGridDof();
176 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
177 fractureVelocity_[phaseIdx].resize(nDof);
178 for (unsigned dofIdx = 0; dofIdx < nDof; ++dofIdx) {
179 fractureVelocity_[phaseIdx][dofIdx].resize(dimWorld);
180 fractureVelocity_[phaseIdx][dofIdx] = 0.0;
181 }
182 }
183 this->resizePhaseBuffer_(fractureVelocityWeight_);
184 }
185 }
186
191 void processElement(const ElementContext& elemCtx)
192 {
193 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
194 return;
195
196 const auto& fractureMapper = elemCtx.simulator().vanguard().fractureMapper();
197
198 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
199 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
200 if (!fractureMapper.isFractureVertex(I))
201 continue;
202
203 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
204 const auto& fs = intQuants.fractureFluidState();
205
206 if (porosityOutput_()) {
207 Opm::Valgrind::CheckDefined(intQuants.fracturePorosity());
208 fracturePorosity_[I] = intQuants.fracturePorosity();
209 }
210 if (intrinsicPermeabilityOutput_()) {
211 const auto& K = intQuants.fractureIntrinsicPermeability();
212 fractureIntrinsicPermeability_[I] = K[0][0];
213 }
214
215 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
216 if (saturationOutput_()) {
217 Opm::Valgrind::CheckDefined(fs.saturation(phaseIdx));
218 fractureSaturation_[phaseIdx][I] = fs.saturation(phaseIdx);
219 }
220 if (mobilityOutput_()) {
221 Opm::Valgrind::CheckDefined(intQuants.fractureMobility(phaseIdx));
222 fractureMobility_[phaseIdx][I] = intQuants.fractureMobility(phaseIdx);
223 }
224 if (relativePermeabilityOutput_()) {
225 Opm::Valgrind::CheckDefined(intQuants.fractureRelativePermeability(phaseIdx));
226 fractureRelativePermeability_[phaseIdx][I] =
227 intQuants.fractureRelativePermeability(phaseIdx);
228 }
229 if (volumeFractionOutput_()) {
230 Opm::Valgrind::CheckDefined(intQuants.fractureVolume());
231 fractureVolumeFraction_[I] += intQuants.fractureVolume();
232 }
233 }
234 }
235
236 if (velocityOutput_()) {
237 // calculate velocities if requested by the simulator
238 for (unsigned scvfIdx = 0; scvfIdx < elemCtx.numInteriorFaces(/*timeIdx=*/0); ++ scvfIdx) {
239 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, /*timeIdx=*/0);
240
241 unsigned i = extQuants.interiorIndex();
242 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
243
244 unsigned j = extQuants.exteriorIndex();
245 unsigned J = elemCtx.globalSpaceIndex(j, /*timeIdx=*/0);
246
247 if (!fractureMapper.isFractureEdge(I, J))
248 continue;
249
250 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
251 Scalar weight =
252 std::max<Scalar>(1e-16, std::abs(extQuants.fractureVolumeFlux(phaseIdx)));
253 Opm::Valgrind::CheckDefined(extQuants.extrusionFactor());
254 assert(extQuants.extrusionFactor() > 0);
255 weight *= extQuants.extrusionFactor();
256
257 Dune::FieldVector<Scalar, dim> v(extQuants.fractureFilterVelocity(phaseIdx));
258 v *= weight;
259
260 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
261 fractureVelocity_[phaseIdx][I][dimIdx] += v[dimIdx];
262 fractureVelocity_[phaseIdx][J][dimIdx] += v[dimIdx];
263 }
264
265 fractureVelocityWeight_[phaseIdx][I] += weight;
266 fractureVelocityWeight_[phaseIdx][J] += weight;
267 }
268 }
269 }
270 }
271
276 {
277 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
278 if (!vtkWriter) {
279 return;
280 }
281
282 if (saturationOutput_())
283 this->commitPhaseBuffer_(baseWriter, "fractureSaturation_%s", fractureSaturation_);
284 if (mobilityOutput_())
285 this->commitPhaseBuffer_(baseWriter, "fractureMobility_%s", fractureMobility_);
286 if (relativePermeabilityOutput_())
287 this->commitPhaseBuffer_(baseWriter, "fractureRelativePerm_%s", fractureRelativePermeability_);
288
289 if (porosityOutput_())
290 this->commitScalarBuffer_(baseWriter, "fracturePorosity", fracturePorosity_);
291 if (intrinsicPermeabilityOutput_())
292 this->commitScalarBuffer_(baseWriter, "fractureIntrinsicPerm", fractureIntrinsicPermeability_);
293 if (volumeFractionOutput_()) {
294 // divide the fracture volume by the total volume of the finite volumes
295 for (unsigned I = 0; I < fractureVolumeFraction_.size(); ++I)
296 fractureVolumeFraction_[I] /= this->simulator_.model().dofTotalVolume(I);
297 this->commitScalarBuffer_(baseWriter, "fractureVolumeFraction", fractureVolumeFraction_);
298 }
299
300 if (velocityOutput_()) {
301 size_t nDof = this->simulator_.model().numGridDof();
302
303 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
304 // first, divide the velocity field by the
305 // respective finite volume's surface area
306 for (unsigned dofIdx = 0; dofIdx < nDof; ++dofIdx)
307 fractureVelocity_[phaseIdx][dofIdx] /=
308 std::max<Scalar>(1e-20, fractureVelocityWeight_[phaseIdx][dofIdx]);
309 // commit the phase velocity
310 char name[512];
311 snprintf(name, 512, "fractureFilterVelocity_%s", FluidSystem::phaseName(phaseIdx).data());
312
313 DiscBaseOutputModule::attachVectorDofData_(baseWriter, fractureVelocity_[phaseIdx], name);
314 }
315 }
316 }
317
318private:
319 static bool saturationOutput_()
320 {
321 static bool val = Parameters::get<TypeTag, Properties::VtkWriteFractureSaturations>();
322 return val;
323 }
324
325 static bool mobilityOutput_()
326 {
327 static bool val = Parameters::get<TypeTag, Properties::VtkWriteFractureMobilities>();
328 return val;
329 }
330
331 static bool relativePermeabilityOutput_()
332 {
333 static bool val = Parameters::get<TypeTag, Properties::VtkWriteFractureRelativePermeabilities>();
334 return val;
335 }
336
337 static bool porosityOutput_()
338 {
339 static bool val = Parameters::get<TypeTag, Properties::VtkWriteFracturePorosity>();
340 return val;
341 }
342
343 static bool intrinsicPermeabilityOutput_()
344 {
345 static bool val = Parameters::get<TypeTag, Properties::VtkWriteFractureIntrinsicPermeabilities>();
346 return val;
347 }
348
349 static bool volumeFractionOutput_()
350 {
351 static bool val = Parameters::get<TypeTag, Properties::VtkWriteFractureVolumeFraction>();
352 return val;
353 }
354
355 static bool velocityOutput_()
356 {
357 static bool val = Parameters::get<TypeTag, Properties::VtkWriteFractureFilterVelocities>();
358 return val;
359 }
360
361 PhaseBuffer fractureSaturation_;
362 PhaseBuffer fractureMobility_;
363 PhaseBuffer fractureRelativePermeability_;
364
365 ScalarBuffer fracturePorosity_;
366 ScalarBuffer fractureVolumeFraction_;
367 ScalarBuffer fractureIntrinsicPermeability_;
368
369 PhaseVectorBuffer fractureVelocity_;
370 PhaseBuffer fractureVelocityWeight_;
371
372 PhaseVectorBuffer potentialGradient_;
373 PhaseBuffer potentialWeight_;
374};
375
376} // namespace Opm
377
378#endif
The base class for writer modules.
Definition: baseoutputmodule.hh:73
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:91
void commitPhaseBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Add a phase-specific buffer to the result file.
Definition: baseoutputmodule.hh:419
void resizePhaseBuffer_(PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase-specific quantity.
Definition: baseoutputmodule.hh:246
const Simulator & simulator_
Definition: baseoutputmodule.hh:519
std::array< VectorBuffer, numPhases > PhaseVectorBuffer
Definition: baseoutputmodule.hh:100
std::array< ScalarBuffer, numPhases > PhaseBuffer
Definition: baseoutputmodule.hh:96
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:316
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a scalar quantity.
Definition: baseoutputmodule.hh:162
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:102
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkdiscretefracturemodule.hh:158
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkdiscretefracturemodule.hh:134
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkdiscretefracturemodule.hh:275
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkdiscretefracturemodule.hh:191
VtkDiscreteFractureModule(const Simulator &simulator)
Definition: vtkdiscretefracturemodule.hh:127
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Definition: blackoilmodel.hh:72
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:242
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: vtkdiscretefracturemodule.hh:48
a tag to mark properties as undefined
Definition: propertysystem.hh:40
Definition: vtkdiscretefracturemodule.hh:64
Definition: vtkdiscretefracturemodule.hh:62
Definition: vtkdiscretefracturemodule.hh:56
Definition: vtkdiscretefracturemodule.hh:60
Definition: vtkdiscretefracturemodule.hh:58
Definition: vtkdiscretefracturemodule.hh:54
Definition: vtkdiscretefracturemodule.hh:66