vtkdiscretefracturemodule.hpp
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 OPM_VTK_DISCRETE_FRACTURE_MODULE_HPP
28#define OPM_VTK_DISCRETE_FRACTURE_MODULE_HPP
29
30#include <dune/common/fvector.hh>
31
32#include <opm/material/common/Valgrind.hpp>
33
35
39
42
43#include <cstdio>
44
45namespace Opm {
46
60template <class TypeTag>
62{
64
68
72
74
75 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
77
78 enum { dim = GridView::dimension };
79 enum { dimWorld = GridView::dimensionworld };
80 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
81
83 using PhaseBuffer = typename ParentType::PhaseBuffer;
85
86public:
87 VtkDiscreteFractureModule(const Simulator& simulator)
88 : ParentType(simulator)
89 {
90 params_.read();
91 }
92
96 static void registerParameters()
97 {
99 }
100
106 {
107 if (params_.saturationOutput_) {
108 this->resizePhaseBuffer_(fractureSaturation_);
109 }
110 if (params_.mobilityOutput_) {
111 this->resizePhaseBuffer_(fractureMobility_);
112 }
113 if (params_.relativePermeabilityOutput_) {
114 this->resizePhaseBuffer_(fractureRelativePermeability_);
115 }
116
117 if (params_.porosityOutput_) {
118 this->resizeScalarBuffer_(fracturePorosity_);
119 }
120 if (params_.intrinsicPermeabilityOutput_) {
121 this->resizeScalarBuffer_(fractureIntrinsicPermeability_);
122 }
123 if (params_.volumeFractionOutput_) {
124 this->resizeScalarBuffer_(fractureVolumeFraction_);
125 }
126
127 if (params_.velocityOutput_) {
128 size_t nDof = this->simulator_.model().numGridDof();
129 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
130 fractureVelocity_[phaseIdx].resize(nDof);
131 for (unsigned dofIdx = 0; dofIdx < nDof; ++dofIdx) {
132 fractureVelocity_[phaseIdx][dofIdx].resize(dimWorld);
133 fractureVelocity_[phaseIdx][dofIdx] = 0.0;
134 }
135 }
136 this->resizePhaseBuffer_(fractureVelocityWeight_);
137 }
138 }
139
144 void processElement(const ElementContext& elemCtx)
145 {
146 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
147 return;
148 }
149
150 const auto& fractureMapper = elemCtx.simulator().vanguard().fractureMapper();
151
152 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
153 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
154 if (!fractureMapper.isFractureVertex(I)) {
155 continue;
156 }
157
158 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
159 const auto& fs = intQuants.fractureFluidState();
160
161 if (params_.porosityOutput_) {
162 Opm::Valgrind::CheckDefined(intQuants.fracturePorosity());
163 fracturePorosity_[I] = intQuants.fracturePorosity();
164 }
165 if (params_.intrinsicPermeabilityOutput_) {
166 const auto& K = intQuants.fractureIntrinsicPermeability();
167 fractureIntrinsicPermeability_[I] = K[0][0];
168 }
169
170 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
171 if (params_.saturationOutput_) {
172 Opm::Valgrind::CheckDefined(fs.saturation(phaseIdx));
173 fractureSaturation_[phaseIdx][I] = fs.saturation(phaseIdx);
174 }
175 if (params_.mobilityOutput_) {
176 Opm::Valgrind::CheckDefined(intQuants.fractureMobility(phaseIdx));
177 fractureMobility_[phaseIdx][I] = intQuants.fractureMobility(phaseIdx);
178 }
179 if (params_.relativePermeabilityOutput_) {
180 Opm::Valgrind::CheckDefined(intQuants.fractureRelativePermeability(phaseIdx));
181 fractureRelativePermeability_[phaseIdx][I] =
182 intQuants.fractureRelativePermeability(phaseIdx);
183 }
184 if (params_.volumeFractionOutput_) {
185 Opm::Valgrind::CheckDefined(intQuants.fractureVolume());
186 fractureVolumeFraction_[I] += intQuants.fractureVolume();
187 }
188 }
189 }
190
191 if (params_.velocityOutput_) {
192 // calculate velocities if requested by the simulator
193 for (unsigned scvfIdx = 0; scvfIdx < elemCtx.numInteriorFaces(/*timeIdx=*/0); ++ scvfIdx) {
194 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, /*timeIdx=*/0);
195
196 unsigned i = extQuants.interiorIndex();
197 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
198
199 unsigned j = extQuants.exteriorIndex();
200 unsigned J = elemCtx.globalSpaceIndex(j, /*timeIdx=*/0);
201
202 if (!fractureMapper.isFractureEdge(I, J)) {
203 continue;
204 }
205
206 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
207 Scalar weight =
208 std::max<Scalar>(1e-16, std::abs(extQuants.fractureVolumeFlux(phaseIdx)));
209 Opm::Valgrind::CheckDefined(extQuants.extrusionFactor());
210 assert(extQuants.extrusionFactor() > 0);
211 weight *= extQuants.extrusionFactor();
212
213 Dune::FieldVector<Scalar, dim> v(extQuants.fractureFilterVelocity(phaseIdx));
214 v *= weight;
215
216 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
217 fractureVelocity_[phaseIdx][I][dimIdx] += v[dimIdx];
218 fractureVelocity_[phaseIdx][J][dimIdx] += v[dimIdx];
219 }
220
221 fractureVelocityWeight_[phaseIdx][I] += weight;
222 fractureVelocityWeight_[phaseIdx][J] += weight;
223 }
224 }
225 }
226 }
227
232 {
233 VtkMultiWriter* vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
234 if (!vtkWriter) {
235 return;
236 }
237
238 if (params_.saturationOutput_) {
239 this->commitPhaseBuffer_(baseWriter, "fractureSaturation_%s", fractureSaturation_);
240 }
241 if (params_.mobilityOutput_) {
242 this->commitPhaseBuffer_(baseWriter, "fractureMobility_%s", fractureMobility_);
243 }
244 if (params_.relativePermeabilityOutput_) {
245 this->commitPhaseBuffer_(baseWriter, "fractureRelativePerm_%s", fractureRelativePermeability_);
246 }
247
248 if (params_.porosityOutput_) {
249 this->commitScalarBuffer_(baseWriter, "fracturePorosity", fracturePorosity_);
250 }
251 if (params_.intrinsicPermeabilityOutput_) {
252 this->commitScalarBuffer_(baseWriter, "fractureIntrinsicPerm", fractureIntrinsicPermeability_);
253 }
254 if (params_.volumeFractionOutput_) {
255 // divide the fracture volume by the total volume of the finite volumes
256 for (unsigned I = 0; I < fractureVolumeFraction_.size(); ++I) {
257 fractureVolumeFraction_[I] /= this->simulator_.model().dofTotalVolume(I);
258 }
259 this->commitScalarBuffer_(baseWriter, "fractureVolumeFraction", fractureVolumeFraction_);
260 }
261
262 if (params_.velocityOutput_) {
263 size_t nDof = this->simulator_.model().numGridDof();
264
265 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
266 // first, divide the velocity field by the
267 // respective finite volume's surface area
268 for (unsigned dofIdx = 0; dofIdx < nDof; ++dofIdx) {
269 fractureVelocity_[phaseIdx][dofIdx] /=
270 std::max<Scalar>(1e-20, fractureVelocityWeight_[phaseIdx][dofIdx]);
271 }
272 // commit the phase velocity
273 char name[512];
274 snprintf(name, 512, "fractureFilterVelocity_%s", FluidSystem::phaseName(phaseIdx).data());
275
276 DiscBaseOutputModule::attachVectorDofData_(baseWriter, fractureVelocity_[phaseIdx], name);
277 }
278 }
279 }
280
281private:
283 PhaseBuffer fractureSaturation_{};
284 PhaseBuffer fractureMobility_{};
285 PhaseBuffer fractureRelativePermeability_{};
286
287 ScalarBuffer fracturePorosity_{};
288 ScalarBuffer fractureVolumeFraction_{};
289 ScalarBuffer fractureIntrinsicPermeability_{};
290
291 PhaseVectorBuffer fractureVelocity_{};
292 PhaseBuffer fractureVelocityWeight_{};
293
294 PhaseVectorBuffer potentialGradient_{};
295 PhaseBuffer potentialWeight_{};
296};
297
298} // namespace Opm
299
300#endif // OPM_VTK_DISCRETE_FRACTURE_MODULE_HPP
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.hpp:62
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkdiscretefracturemodule.hpp:105
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkdiscretefracturemodule.hpp:96
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkdiscretefracturemodule.hpp:231
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkdiscretefracturemodule.hpp:144
VtkDiscreteFractureModule(const Simulator &simulator)
Definition: vtkdiscretefracturemodule.hpp:87
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Declare the properties used by the infrastructure code of the finite volume discretizations.
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.
Struct holding the parameters for VtkDiscreteFractureModule.
Definition: vtkdiscretefractureparams.hpp:49
void read()
Reads the parameter values from the parameter system.
bool velocityOutput_
Definition: vtkdiscretefractureparams.hpp:62
bool relativePermeabilityOutput_
Definition: vtkdiscretefractureparams.hpp:58
bool porosityOutput_
Definition: vtkdiscretefractureparams.hpp:59
static void registerParameters()
Registers the parameters in parameter system.
bool intrinsicPermeabilityOutput_
Definition: vtkdiscretefractureparams.hpp:60
bool volumeFractionOutput_
Definition: vtkdiscretefractureparams.hpp:61
bool saturationOutput_
Definition: vtkdiscretefractureparams.hpp:56
bool mobilityOutput_
Definition: vtkdiscretefractureparams.hpp:57