vtkmultiphasemodule.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_MULTI_PHASE_MODULE_HPP
28#define OPM_VTK_MULTI_PHASE_MODULE_HPP
29
30#include <dune/common/fvector.hh>
31
32#include <opm/material/common/MathToolbox.hpp>
33#include <opm/material/common/Valgrind.hpp>
34
36
40
43
44#include <cstdio>
45
46namespace Opm {
47
66template<class TypeTag>
68{
70
74
78
79 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
81
82 enum { dimWorld = GridView::dimensionworld };
83 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
84
86 using VectorBuffer = typename ParentType::VectorBuffer;
88 using PhaseBuffer = typename ParentType::PhaseBuffer;
89
90 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
91
92 using PhaseVectorBuffer = std::array<VectorBuffer, numPhases>;
93
94public:
95 VtkMultiPhaseModule(const Simulator& simulator)
96 : ParentType(simulator)
97 {
98 params_.read();
99 }
100
104 static void registerParameters()
105 {
107 }
108
114 {
115 if (params_.extrusionFactorOutput_) {
116 this->resizeScalarBuffer_(extrusionFactor_);
117 }
118 if (params_.pressureOutput_) {
119 this->resizePhaseBuffer_(pressure_);
120 }
121 if (params_.densityOutput_) {
122 this->resizePhaseBuffer_(density_);
123 }
124 if (params_.saturationOutput_) {
125 this->resizePhaseBuffer_(saturation_);
126 }
127 if (params_.mobilityOutput_) {
128 this->resizePhaseBuffer_(mobility_);
129 }
130 if (params_.relativePermeabilityOutput_) {
131 this->resizePhaseBuffer_(relativePermeability_);
132 }
133 if (params_.viscosityOutput_) {
134 this->resizePhaseBuffer_(viscosity_);
135 }
136 if (params_.averageMolarMassOutput_) {
137 this->resizePhaseBuffer_(averageMolarMass_);
138 }
139
140 if (params_.porosityOutput_) {
141 this->resizeScalarBuffer_(porosity_);
142 }
143 if (params_.intrinsicPermeabilityOutput_) {
144 this->resizeTensorBuffer_(intrinsicPermeability_);
145 }
146
147 if (params_.velocityOutput_) {
148 size_t nDof = this->simulator_.model().numGridDof();
149 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
150 velocity_[phaseIdx].resize(nDof);
151 for (unsigned dofIdx = 0; dofIdx < nDof; ++ dofIdx) {
152 velocity_[phaseIdx][dofIdx].resize(dimWorld);
153 velocity_[phaseIdx][dofIdx] = 0.0;
154 }
155 }
156 this->resizePhaseBuffer_(velocityWeight_);
157 }
158
159 if (params_.potentialGradientOutput_) {
160 size_t nDof = this->simulator_.model().numGridDof();
161 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
162 potentialGradient_[phaseIdx].resize(nDof);
163 for (unsigned dofIdx = 0; dofIdx < nDof; ++ dofIdx) {
164 potentialGradient_[phaseIdx][dofIdx].resize(dimWorld);
165 potentialGradient_[phaseIdx][dofIdx] = 0.0;
166 }
167 }
168
169 this->resizePhaseBuffer_(potentialWeight_);
170 }
171 }
172
177 void processElement(const ElementContext& elemCtx)
178 {
179 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
180 return;
181 }
182
183 const auto& problem = elemCtx.problem();
184 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
185 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
186 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
187 const auto& fs = intQuants.fluidState();
188
189 if (params_.extrusionFactorOutput_) {
190 extrusionFactor_[I] = intQuants.extrusionFactor();
191 }
192 if (params_.porosityOutput_) {
193 porosity_[I] = getValue(intQuants.porosity());
194 }
195
196 if (params_.intrinsicPermeabilityOutput_) {
197 const auto& K = problem.intrinsicPermeability(elemCtx, i, /*timeIdx=*/0);
198 for (unsigned rowIdx = 0; rowIdx < K.rows; ++rowIdx) {
199 for (unsigned colIdx = 0; colIdx < K.cols; ++colIdx) {
200 intrinsicPermeability_[I][rowIdx][colIdx] = K[rowIdx][colIdx];
201 }
202 }
203 }
204
205 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
206 if (!FluidSystem::phaseIsActive(phaseIdx)) {
207 continue;
208 }
209 if (params_.pressureOutput_) {
210 pressure_[phaseIdx][I] = getValue(fs.pressure(phaseIdx));
211 }
212 if (params_.densityOutput_) {
213 density_[phaseIdx][I] = getValue(fs.density(phaseIdx));
214 }
215 if (params_.saturationOutput_) {
216 saturation_[phaseIdx][I] = getValue(fs.saturation(phaseIdx));
217 }
218 if (params_.mobilityOutput_) {
219 mobility_[phaseIdx][I] = getValue(intQuants.mobility(phaseIdx));
220 }
221 if (params_.relativePermeabilityOutput_) {
222 relativePermeability_[phaseIdx][I] = getValue(intQuants.relativePermeability(phaseIdx));
223 }
224 if (params_.viscosityOutput_) {
225 viscosity_[phaseIdx][I] = getValue(fs.viscosity(phaseIdx));
226 }
227 if (params_.averageMolarMassOutput_) {
228 averageMolarMass_[phaseIdx][I] = getValue(fs.averageMolarMass(phaseIdx));
229 }
230 }
231 }
232
233 if (params_.potentialGradientOutput_) {
234 // calculate velocities if requested
235 for (unsigned faceIdx = 0; faceIdx < elemCtx.numInteriorFaces(/*timeIdx=*/0); ++ faceIdx) {
236 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, /*timeIdx=*/0);
237
238 unsigned i = extQuants.interiorIndex();
239 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
240
241 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
242 Scalar weight = extQuants.extrusionFactor();
243
244 potentialWeight_[phaseIdx][I] += weight;
245
246 const auto& inputPGrad = extQuants.potentialGrad(phaseIdx);
247 DimVector pGrad;
248 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
249 pGrad[dimIdx] = getValue(inputPGrad[dimIdx])*weight;
250 }
251 potentialGradient_[phaseIdx][I] += pGrad;
252 } // end for all phases
253 } // end for all faces
254 }
255
256 if (params_.velocityOutput_) {
257 // calculate velocities if requested
258 for (unsigned faceIdx = 0; faceIdx < elemCtx.numInteriorFaces(/*timeIdx=*/0); ++ faceIdx) {
259 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, /*timeIdx=*/0);
260
261 unsigned i = extQuants.interiorIndex();
262 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
263
264 unsigned j = extQuants.exteriorIndex();
265 unsigned J = elemCtx.globalSpaceIndex(j, /*timeIdx=*/0);
266
267 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
268 Scalar weight = std::max<Scalar>(1e-16,
269 std::abs(getValue(extQuants.volumeFlux(phaseIdx))));
270 Valgrind::CheckDefined(extQuants.extrusionFactor());
271 assert(extQuants.extrusionFactor() > 0);
272 weight *= extQuants.extrusionFactor();
273
274 const auto& inputV = extQuants.filterVelocity(phaseIdx);
275 DimVector v;
276 for (unsigned k = 0; k < dimWorld; ++k) {
277 v[k] = getValue(inputV[k]);
278 }
279 if (v.two_norm() > 1e-20) {
280 weight /= v.two_norm();
281 }
282 v *= weight;
283
284 velocity_[phaseIdx][I] += v;
285 velocity_[phaseIdx][J] += v;
286
287 velocityWeight_[phaseIdx][I] += weight;
288 velocityWeight_[phaseIdx][J] += weight;
289 } // end for all phases
290 } // end for all faces
291 }
292 }
293
298 {
299 VtkMultiWriter* vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
300 if (!vtkWriter) {
301 return;
302 }
303
304 if (params_.extrusionFactorOutput_) {
305 this->commitScalarBuffer_(baseWriter, "extrusionFactor", extrusionFactor_);
306 }
307 if (params_.pressureOutput_) {
308 this->commitPhaseBuffer_(baseWriter, "pressure_%s", pressure_);
309 }
310 if (params_.densityOutput_) {
311 this->commitPhaseBuffer_(baseWriter, "density_%s", density_);
312 }
313 if (params_.saturationOutput_) {
314 this->commitPhaseBuffer_(baseWriter, "saturation_%s", saturation_);
315 }
316 if (params_.mobilityOutput_) {
317 this->commitPhaseBuffer_(baseWriter, "mobility_%s", mobility_);
318 }
319 if (params_.relativePermeabilityOutput_) {
320 this->commitPhaseBuffer_(baseWriter, "relativePerm_%s", relativePermeability_);
321 }
322 if (params_.viscosityOutput_) {
323 this->commitPhaseBuffer_(baseWriter, "viscosity_%s", viscosity_);
324 }
325 if (params_.averageMolarMassOutput_) {
326 this->commitPhaseBuffer_(baseWriter, "averageMolarMass_%s", averageMolarMass_);
327 }
328
329 if (params_.porosityOutput_) {
330 this->commitScalarBuffer_(baseWriter, "porosity", porosity_);
331 }
332 if (params_.intrinsicPermeabilityOutput_) {
333 this->commitTensorBuffer_(baseWriter, "intrinsicPerm", intrinsicPermeability_);
334 }
335
336 if (params_.velocityOutput_) {
337 size_t numDof = this->simulator_.model().numGridDof();
338
339 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
340 // first, divide the velocity field by the
341 // respective finite volume's surface area
342 for (unsigned i = 0; i < numDof; ++i) {
343 velocity_[phaseIdx][i] /= velocityWeight_[phaseIdx][i];
344 }
345 // commit the phase velocity
346 char name[512];
347 snprintf(name, 512, "filterVelocity_%s", FluidSystem::phaseName(phaseIdx).data());
348
349 DiscBaseOutputModule::attachVectorDofData_(baseWriter, velocity_[phaseIdx], name);
350 }
351 }
352
353 if (params_.potentialGradientOutput_) {
354 size_t numDof = this->simulator_.model().numGridDof();
355
356 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
357 // first, divide the velocity field by the
358 // respective finite volume's surface area
359 for (unsigned i = 0; i < numDof; ++i) {
360 potentialGradient_[phaseIdx][i] /= potentialWeight_[phaseIdx][i];
361 }
362 // commit the phase velocity
363 char name[512];
364 snprintf(name, 512, "gradP_%s", FluidSystem::phaseName(phaseIdx).data());
365
366 DiscBaseOutputModule::attachVectorDofData_(baseWriter,
367 potentialGradient_[phaseIdx],
368 name);
369 }
370 }
371 }
372
381 bool needExtensiveQuantities() const final
382 {
383 return params_.velocityOutput_ || params_.potentialGradientOutput_;
384 }
385
386private:
387 VtkMultiPhaseParams params_{};
388 ScalarBuffer extrusionFactor_{};
389 PhaseBuffer pressure_{};
390 PhaseBuffer density_{};
391 PhaseBuffer saturation_{};
392 PhaseBuffer mobility_{};
393 PhaseBuffer relativePermeability_{};
394 PhaseBuffer viscosity_{};
395 PhaseBuffer averageMolarMass_{};
396
397 ScalarBuffer porosity_{};
398 TensorBuffer intrinsicPermeability_{};
399
400 PhaseVectorBuffer velocity_{};
401 PhaseBuffer velocityWeight_{};
402
403 PhaseVectorBuffer potentialGradient_{};
404 PhaseBuffer potentialWeight_{};
405};
406
407} // namespace Opm
408
409#endif // OPM_VTK_MULTI_PHASE_MODULE_HPP
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.hpp:68
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities seen on an element.
Definition: vtkmultiphasemodule.hpp:177
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkmultiphasemodule.hpp:104
VtkMultiPhaseModule(const Simulator &simulator)
Definition: vtkmultiphasemodule.hpp:95
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkmultiphasemodule.hpp:113
bool needExtensiveQuantities() const final
Returns true iff the module needs to access the extensive quantities of a context to do its job.
Definition: vtkmultiphasemodule.hpp:381
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkmultiphasemodule.hpp:297
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 VtkMultiPhaseModule.
Definition: vtkmultiphaseparams.hpp:54
bool averageMolarMassOutput_
Definition: vtkmultiphaseparams.hpp:68
bool pressureOutput_
Definition: vtkmultiphaseparams.hpp:62
bool extrusionFactorOutput_
Definition: vtkmultiphaseparams.hpp:61
bool porosityOutput_
Definition: vtkmultiphaseparams.hpp:69
bool velocityOutput_
Definition: vtkmultiphaseparams.hpp:71
bool densityOutput_
Definition: vtkmultiphaseparams.hpp:63
static void registerParameters()
Registers the parameters in parameter system.
bool saturationOutput_
Definition: vtkmultiphaseparams.hpp:64
bool potentialGradientOutput_
Definition: vtkmultiphaseparams.hpp:72
bool intrinsicPermeabilityOutput_
Definition: vtkmultiphaseparams.hpp:70
void read()
Reads the parameter values from the parameter system.
bool viscosityOutput_
Definition: vtkmultiphaseparams.hpp:67
bool relativePermeabilityOutput_
Definition: vtkmultiphaseparams.hpp:66
bool mobilityOutput_
Definition: vtkmultiphaseparams.hpp:65