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