vtkmultiphasemodule.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_MULTI_PHASE_MODULE_HH
28#define EWOMS_VTK_MULTI_PHASE_MODULE_HH
29
30#include <dune/common/fvector.hh>
31
32#include <opm/material/common/MathToolbox.hpp>
33#include <opm/material/common/Valgrind.hpp>
34
36
39
42
43#include <cstdio>
44
45namespace Opm::Parameters {
46
47// set default values for what quantities to output
48struct VtkWriteExtrusionFactor { static constexpr bool value = false; };
49struct VtkWritePressures { static constexpr bool value = true; };
50struct VtkWriteDensities { static constexpr bool value = true; };
51struct VtkWriteSaturations { static constexpr bool value = true; };
52struct VtkWriteMobilities { static constexpr bool value = false; };
53struct VtkWriteRelativePermeabilities { static constexpr bool value = true; };
54struct VtkWriteViscosities { static constexpr bool value = false; };
55struct VtkWriteAverageMolarMasses { static constexpr bool value = false; };
56struct VtkWritePorosity { static constexpr bool value = true; };
57struct VtkWriteIntrinsicPermeabilities { static constexpr bool value = false; };
58struct VtkWritePotentialGradients { static constexpr bool value = false; };
59struct VtkWriteFilterVelocities { static constexpr bool value = false; };
60
61} // namespace Opm::Parameters
62
63namespace Opm {
64
83template<class TypeTag>
85{
87
91
95
96 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
98
99 enum { dimWorld = GridView::dimensionworld };
100 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
101
102 using ScalarBuffer = typename ParentType::ScalarBuffer;
103 using VectorBuffer = typename ParentType::VectorBuffer;
104 using TensorBuffer = typename ParentType::TensorBuffer;
105 using PhaseBuffer = typename ParentType::PhaseBuffer;
106
107 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
108
109 using PhaseVectorBuffer = std::array<VectorBuffer, numPhases>;
110
111public:
112 VtkMultiPhaseModule(const Simulator& simulator)
113 : ParentType(simulator)
114 {}
115
119 static void registerParameters()
120 {
121 Parameters::Register<Parameters::VtkWriteExtrusionFactor>
122 ("Include the extrusion factor of the degrees of freedom into the VTK output files");
123 Parameters::Register<Parameters::VtkWritePressures>
124 ("Include the phase pressures in the VTK output files");
125 Parameters::Register<Parameters::VtkWriteDensities>
126 ("Include the phase densities in the VTK output files");
127 Parameters::Register<Parameters::VtkWriteSaturations>
128 ("Include the phase saturations in the VTK output files");
129 Parameters::Register<Parameters::VtkWriteMobilities>
130 ("Include the phase mobilities in the VTK output files");
131 Parameters::Register<Parameters::VtkWriteRelativePermeabilities>
132 ("Include the phase relative permeabilities in the VTK output files");
133 Parameters::Register<Parameters::VtkWriteViscosities>
134 ("Include component phase viscosities in the VTK output files");
135 Parameters::Register<Parameters::VtkWriteAverageMolarMasses>
136 ("Include the average phase mass in the VTK output files");
137 Parameters::Register<Parameters::VtkWritePorosity>
138 ("Include the porosity in the VTK output files");
139 Parameters::Register<Parameters::VtkWriteIntrinsicPermeabilities>
140 ("Include the intrinsic permeability in the VTK output files");
141 Parameters::Register<Parameters::VtkWriteFilterVelocities>
142 ("Include in the filter velocities of the phases the VTK output files");
143 Parameters::Register<Parameters::VtkWritePotentialGradients>
144 ("Include the phase pressure potential gradients in the VTK output files");
145 }
146
152 {
153 if (extrusionFactorOutput_()) this->resizeScalarBuffer_(extrusionFactor_);
154 if (pressureOutput_()) this->resizePhaseBuffer_(pressure_);
155 if (densityOutput_()) this->resizePhaseBuffer_(density_);
156 if (saturationOutput_()) this->resizePhaseBuffer_(saturation_);
157 if (mobilityOutput_()) this->resizePhaseBuffer_(mobility_);
158 if (relativePermeabilityOutput_()) this->resizePhaseBuffer_(relativePermeability_);
159 if (viscosityOutput_()) this->resizePhaseBuffer_(viscosity_);
160 if (averageMolarMassOutput_()) this->resizePhaseBuffer_(averageMolarMass_);
161
162 if (porosityOutput_()) this->resizeScalarBuffer_(porosity_);
163 if (intrinsicPermeabilityOutput_()) this->resizeTensorBuffer_(intrinsicPermeability_);
164
165 if (velocityOutput_()) {
166 size_t nDof = this->simulator_.model().numGridDof();
167 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
168 velocity_[phaseIdx].resize(nDof);
169 for (unsigned dofIdx = 0; dofIdx < nDof; ++ dofIdx) {
170 velocity_[phaseIdx][dofIdx].resize(dimWorld);
171 velocity_[phaseIdx][dofIdx] = 0.0;
172 }
173 }
174 this->resizePhaseBuffer_(velocityWeight_);
175 }
176
177 if (potentialGradientOutput_()) {
178 size_t nDof = this->simulator_.model().numGridDof();
179 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
180 potentialGradient_[phaseIdx].resize(nDof);
181 for (unsigned dofIdx = 0; dofIdx < nDof; ++ dofIdx) {
182 potentialGradient_[phaseIdx][dofIdx].resize(dimWorld);
183 potentialGradient_[phaseIdx][dofIdx] = 0.0;
184 }
185 }
186
187 this->resizePhaseBuffer_(potentialWeight_);
188 }
189 }
190
195 void processElement(const ElementContext& elemCtx)
196 {
197 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
198 return;
199 }
200
201 const auto& problem = elemCtx.problem();
202 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
203 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
204 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
205 const auto& fs = intQuants.fluidState();
206
207 if (extrusionFactorOutput_()) extrusionFactor_[I] = intQuants.extrusionFactor();
208 if (porosityOutput_()) porosity_[I] = getValue(intQuants.porosity());
209
210 if (intrinsicPermeabilityOutput_()) {
211 const auto& K = problem.intrinsicPermeability(elemCtx, i, /*timeIdx=*/0);
212 for (unsigned rowIdx = 0; rowIdx < K.rows; ++rowIdx)
213 for (unsigned colIdx = 0; colIdx < K.cols; ++colIdx)
214 intrinsicPermeability_[I][rowIdx][colIdx] = K[rowIdx][colIdx];
215 }
216
217 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
218 if (!FluidSystem::phaseIsActive(phaseIdx)) {
219 continue;
220 }
221 if (pressureOutput_())
222 pressure_[phaseIdx][I] = getValue(fs.pressure(phaseIdx));
223 if (densityOutput_())
224 density_[phaseIdx][I] = getValue(fs.density(phaseIdx));
225 if (saturationOutput_())
226 saturation_[phaseIdx][I] = getValue(fs.saturation(phaseIdx));
227 if (mobilityOutput_())
228 mobility_[phaseIdx][I] = getValue(intQuants.mobility(phaseIdx));
229 if (relativePermeabilityOutput_())
230 relativePermeability_[phaseIdx][I] = getValue(intQuants.relativePermeability(phaseIdx));
231 if (viscosityOutput_())
232 viscosity_[phaseIdx][I] = getValue(fs.viscosity(phaseIdx));
233 if (averageMolarMassOutput_())
234 averageMolarMass_[phaseIdx][I] = getValue(fs.averageMolarMass(phaseIdx));
235 }
236 }
237
238 if (potentialGradientOutput_()) {
239 // calculate velocities if requested
240 for (unsigned faceIdx = 0; faceIdx < elemCtx.numInteriorFaces(/*timeIdx=*/0); ++ faceIdx) {
241 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, /*timeIdx=*/0);
242
243 unsigned i = extQuants.interiorIndex();
244 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
245
246 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
247 Scalar weight = extQuants.extrusionFactor();
248
249 potentialWeight_[phaseIdx][I] += weight;
250
251 const auto& inputPGrad = extQuants.potentialGrad(phaseIdx);
252 DimVector pGrad;
253 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
254 pGrad[dimIdx] = getValue(inputPGrad[dimIdx])*weight;
255 potentialGradient_[phaseIdx][I] += pGrad;
256 } // end for all phases
257 } // end for all faces
258 }
259
260 if (velocityOutput_()) {
261 // calculate velocities if requested
262 for (unsigned faceIdx = 0; faceIdx < elemCtx.numInteriorFaces(/*timeIdx=*/0); ++ faceIdx) {
263 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, /*timeIdx=*/0);
264
265 unsigned i = extQuants.interiorIndex();
266 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
267
268 unsigned j = extQuants.exteriorIndex();
269 unsigned J = elemCtx.globalSpaceIndex(j, /*timeIdx=*/0);
270
271 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
272 Scalar weight = std::max<Scalar>(1e-16,
273 std::abs(getValue(extQuants.volumeFlux(phaseIdx))));
274 Valgrind::CheckDefined(extQuants.extrusionFactor());
275 assert(extQuants.extrusionFactor() > 0);
276 weight *= extQuants.extrusionFactor();
277
278 const auto& inputV = extQuants.filterVelocity(phaseIdx);
279 DimVector v;
280 for (unsigned k = 0; k < dimWorld; ++k)
281 v[k] = getValue(inputV[k]);
282 if (v.two_norm() > 1e-20)
283 weight /= v.two_norm();
284 v *= weight;
285
286 velocity_[phaseIdx][I] += v;
287 velocity_[phaseIdx][J] += v;
288
289 velocityWeight_[phaseIdx][I] += weight;
290 velocityWeight_[phaseIdx][J] += weight;
291 } // end for all phases
292 } // end for all faces
293 }
294 }
295
300 {
301 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
302 if (!vtkWriter)
303 return;
304
305 if (extrusionFactorOutput_())
306 this->commitScalarBuffer_(baseWriter, "extrusionFactor", extrusionFactor_);
307 if (pressureOutput_())
308 this->commitPhaseBuffer_(baseWriter, "pressure_%s", pressure_);
309 if (densityOutput_())
310 this->commitPhaseBuffer_(baseWriter, "density_%s", density_);
311 if (saturationOutput_())
312 this->commitPhaseBuffer_(baseWriter, "saturation_%s", saturation_);
313 if (mobilityOutput_())
314 this->commitPhaseBuffer_(baseWriter, "mobility_%s", mobility_);
315 if (relativePermeabilityOutput_())
316 this->commitPhaseBuffer_(baseWriter, "relativePerm_%s", relativePermeability_);
317 if (viscosityOutput_())
318 this->commitPhaseBuffer_(baseWriter, "viscosity_%s", viscosity_);
319 if (averageMolarMassOutput_())
320 this->commitPhaseBuffer_(baseWriter, "averageMolarMass_%s", averageMolarMass_);
321
322 if (porosityOutput_())
323 this->commitScalarBuffer_(baseWriter, "porosity", porosity_);
324 if (intrinsicPermeabilityOutput_())
325 this->commitTensorBuffer_(baseWriter, "intrinsicPerm", intrinsicPermeability_);
326
327 if (velocityOutput_()) {
328 size_t numDof = this->simulator_.model().numGridDof();
329
330 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
331 // first, divide the velocity field by the
332 // respective finite volume's surface area
333 for (unsigned i = 0; i < numDof; ++i)
334 velocity_[phaseIdx][i] /= velocityWeight_[phaseIdx][i];
335 // commit the phase velocity
336 char name[512];
337 snprintf(name, 512, "filterVelocity_%s", FluidSystem::phaseName(phaseIdx).data());
338
339 DiscBaseOutputModule::attachVectorDofData_(baseWriter, velocity_[phaseIdx], name);
340 }
341 }
342
343 if (potentialGradientOutput_()) {
344 size_t numDof = this->simulator_.model().numGridDof();
345
346 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
347 // first, divide the velocity field by the
348 // respective finite volume's surface area
349 for (unsigned i = 0; i < numDof; ++i)
350 potentialGradient_[phaseIdx][i] /= potentialWeight_[phaseIdx][i];
351 // commit the phase velocity
352 char name[512];
353 snprintf(name, 512, "gradP_%s", FluidSystem::phaseName(phaseIdx).data());
354
355 DiscBaseOutputModule::attachVectorDofData_(baseWriter,
356 potentialGradient_[phaseIdx],
357 name);
358 }
359 }
360 }
361
370 virtual bool needExtensiveQuantities() const final
371 {
372 return velocityOutput_() || potentialGradientOutput_();
373 }
374
375private:
376 static bool extrusionFactorOutput_()
377 {
378 static bool val = Parameters::Get<Parameters::VtkWriteExtrusionFactor>();
379 return val;
380 }
381
382 static bool pressureOutput_()
383 {
384 static bool val = Parameters::Get<Parameters::VtkWritePressures>();
385 return val;
386 }
387
388 static bool densityOutput_()
389 {
390 static bool val = Parameters::Get<Parameters::VtkWriteDensities>();
391 return val;
392 }
393
394 static bool saturationOutput_()
395 {
396 static bool val = Parameters::Get<Parameters::VtkWriteSaturations>();
397 return val;
398 }
399
400 static bool mobilityOutput_()
401 {
402 static bool val = Parameters::Get<Parameters::VtkWriteMobilities>();
403 return val;
404 }
405
406 static bool relativePermeabilityOutput_()
407 {
408 static bool val = Parameters::Get<Parameters::VtkWriteRelativePermeabilities>();
409 return val;
410 }
411
412 static bool viscosityOutput_()
413 {
414 static bool val = Parameters::Get<Parameters::VtkWriteViscosities>();
415 return val;
416 }
417
418 static bool averageMolarMassOutput_()
419 {
420 static bool val = Parameters::Get<Parameters::VtkWriteAverageMolarMasses>();
421 return val;
422 }
423
424 static bool porosityOutput_()
425 {
426 static bool val = Parameters::Get<Parameters::VtkWritePorosity>();
427 return val;
428 }
429
430 static bool intrinsicPermeabilityOutput_()
431 {
432 static bool val = Parameters::Get<Parameters::VtkWriteIntrinsicPermeabilities>();
433 return val;
434 }
435
436 static bool velocityOutput_()
437 {
438 static bool val = Parameters::Get<Parameters::VtkWriteFilterVelocities>();
439 return val;
440 }
441
442 static bool potentialGradientOutput_()
443 {
444 static bool val = Parameters::Get<Parameters::VtkWritePotentialGradients>();
445 return val;
446 }
447
448 ScalarBuffer extrusionFactor_;
449 PhaseBuffer pressure_;
450 PhaseBuffer density_;
451 PhaseBuffer saturation_;
452 PhaseBuffer mobility_;
453 PhaseBuffer relativePermeability_;
454 PhaseBuffer viscosity_;
455 PhaseBuffer averageMolarMass_;
456
457 ScalarBuffer porosity_;
458 TensorBuffer intrinsicPermeability_;
459
460 PhaseVectorBuffer velocity_;
461 PhaseBuffer velocityWeight_;
462
463 PhaseVectorBuffer potentialGradient_;
464 PhaseBuffer potentialWeight_;
465};
466
467} // namespace Opm
468
469#endif
The base class for writer modules.
Definition: baseoutputmodule.hh:67
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:85
BaseOutputWriter::TensorBuffer TensorBuffer
Definition: baseoutputmodule.hh:87
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
void commitTensorBuffer_(BaseOutputWriter &baseWriter, const char *name, TensorBuffer &buffer, BufferType bufferType=DofBuffer)
Add a buffer containing tensorial quantities to the result file.
Definition: baseoutputmodule.hh:288
std::array< ScalarBuffer, numPhases > PhaseBuffer
Definition: baseoutputmodule.hh:90
BaseOutputWriter::VectorBuffer VectorBuffer
Definition: baseoutputmodule.hh:86
void resizeTensorBuffer_(TensorBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a tensorial quantity.
Definition: baseoutputmodule.hh:166
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 multiple fluid phase...
Definition: vtkmultiphasemodule.hh:85
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities seen on an element.
Definition: vtkmultiphasemodule.hh:195
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkmultiphasemodule.hh:119
VtkMultiPhaseModule(const Simulator &simulator)
Definition: vtkmultiphasemodule.hh:112
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkmultiphasemodule.hh:151
virtual bool needExtensiveQuantities() const final
Returns true iff the module needs to access the extensive quantities of a context to do its job.
Definition: vtkmultiphasemodule.hh:370
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkmultiphasemodule.hh:299
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: vtkmultiphasemodule.hh:55
static constexpr bool value
Definition: vtkmultiphasemodule.hh:55
Definition: vtkmultiphasemodule.hh:50
static constexpr bool value
Definition: vtkmultiphasemodule.hh:50
Definition: vtkmultiphasemodule.hh:48
static constexpr bool value
Definition: vtkmultiphasemodule.hh:48
Definition: vtkmultiphasemodule.hh:59
static constexpr bool value
Definition: vtkmultiphasemodule.hh:59
Definition: vtkmultiphasemodule.hh:57
static constexpr bool value
Definition: vtkmultiphasemodule.hh:57
Definition: vtkmultiphasemodule.hh:52
static constexpr bool value
Definition: vtkmultiphasemodule.hh:52
Definition: vtkmultiphasemodule.hh:56
static constexpr bool value
Definition: vtkmultiphasemodule.hh:56
Definition: vtkmultiphasemodule.hh:58
static constexpr bool value
Definition: vtkmultiphasemodule.hh:58
Definition: vtkmultiphasemodule.hh:49
static constexpr bool value
Definition: vtkmultiphasemodule.hh:49
Definition: vtkmultiphasemodule.hh:53
static constexpr bool value
Definition: vtkmultiphasemodule.hh:53
Definition: vtkmultiphasemodule.hh:51
static constexpr bool value
Definition: vtkmultiphasemodule.hh:51
Definition: vtkmultiphasemodule.hh:54
static constexpr bool value
Definition: vtkmultiphasemodule.hh:54