baseoutputmodule.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_BASE_OUTPUT_MODULE_HH
28#define EWOMS_BASE_OUTPUT_MODULE_HH
29
30#include "baseoutputwriter.hh"
31
37
38#include <dune/istl/bvector.hh>
39#include <dune/common/fvector.hh>
40
41#include <array>
42#include <sstream>
43#include <string>
44#include <string_view>
45#include <vector>
46
47#include <cstdio>
48
49namespace Opm::Properties {
50
51template <class TypeTag, class MyTypeTag>
52struct FluidSystem;
53
54} // namespace Opm::Properties
55
56namespace Opm {
57
58#if __GNUC__ || __clang__
59#pragma GCC diagnostic push
60#pragma GCC diagnostic ignored "-Wpragmas"
61#pragma GCC diagnostic ignored "-Wformat-nonliteral"
62#endif
63
71template<class TypeTag>
73{
81
82 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
83 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
84 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
85 enum { dim = GridView::dimension };
86 enum { dimWorld = GridView::dimensionworld };
87 using Vector = BaseOutputWriter::Vector;
88 using Tensor = BaseOutputWriter::Tensor;
89
90public:
94
95 using EqBuffer = std::array<ScalarBuffer, numEq>;
96 using PhaseBuffer = std::array<ScalarBuffer, numPhases>;
97 using ComponentBuffer = std::array<ScalarBuffer, numComponents>;
98 using PhaseComponentBuffer = std::array<std::array<ScalarBuffer, numComponents>, numPhases>;
99
100 using PhaseVectorBuffer = std::array<VectorBuffer, numPhases>;
101
102 BaseOutputModule(const Simulator& simulator)
103 : simulator_(simulator)
104 {}
105
107 {}
108
117 virtual void allocBuffers() = 0;
118
127 virtual void processElement(const ElementContext& elemCtx) = 0;
128
132 virtual void commitBuffers(BaseOutputWriter& writer) = 0;
133
144 virtual bool needExtensiveQuantities() const
145 { return false; }
146
147protected:
151
154
157 };
158
163 BufferType bufferType = DofBuffer)
164 {
165 size_t n;
166 if (bufferType == VertexBuffer)
167 n = static_cast<size_t>(simulator_.gridView().size(dim));
168 else if (bufferType == ElementBuffer)
169 n = static_cast<size_t>(simulator_.gridView().size(0));
170 else if (bufferType == DofBuffer)
171 n = simulator_.model().numGridDof();
172 else
173 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
174
175 buffer.resize(n);
176 std::fill(buffer.begin(), buffer.end(), 0.0);
177 }
178
183 BufferType bufferType = DofBuffer)
184 {
185 size_t n;
186 if (bufferType == VertexBuffer)
187 n = static_cast<size_t>(simulator_.gridView().size(dim));
188 else if (bufferType == ElementBuffer)
189 n = static_cast<size_t>(simulator_.gridView().size(0));
190 else if (bufferType == DofBuffer)
191 n = simulator_.model().numGridDof();
192 else
193 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
194
195 buffer.resize(n);
196 Tensor nullMatrix(dimWorld, dimWorld, 0.0);
197 std::fill(buffer.begin(), buffer.end(), nullMatrix);
198 }
199
201 BufferType bufferType = DofBuffer)
202 {
203 size_t n;
204 if (bufferType == VertexBuffer)
205 n = static_cast<size_t>(simulator_.gridView().size(dim));
206 else if (bufferType == ElementBuffer)
207 n = static_cast<size_t>(simulator_.gridView().size(0));
208 else if (bufferType == DofBuffer)
209 n = simulator_.model().numGridDof();
210 else
211 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
212
213 buffer.resize(n);
214 Vector zerovector(dimWorld,0.0);
215 zerovector = 0.0;
216 std::fill(buffer.begin(), buffer.end(), zerovector);
217 }
218
224 BufferType bufferType = DofBuffer)
225 {
226 size_t n;
227 if (bufferType == VertexBuffer)
228 n = static_cast<size_t>(simulator_.gridView().size(dim));
229 else if (bufferType == ElementBuffer)
230 n = static_cast<size_t>(simulator_.gridView().size(0));
231 else if (bufferType == DofBuffer)
232 n = simulator_.model().numGridDof();
233 else
234 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
235
236 for (unsigned i = 0; i < numEq; ++i) {
237 buffer[i].resize(n);
238 std::fill(buffer[i].begin(), buffer[i].end(), 0.0);
239 }
240 }
241
247 BufferType bufferType = DofBuffer)
248 {
249 size_t n;
250 if (bufferType == VertexBuffer)
251 n = static_cast<size_t>(simulator_.gridView().size(dim));
252 else if (bufferType == ElementBuffer)
253 n = static_cast<size_t>(simulator_.gridView().size(0));
254 else if (bufferType == DofBuffer)
255 n = simulator_.model().numGridDof();
256 else
257 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
258
259 for (unsigned i = 0; i < numPhases; ++i) {
260 buffer[i].resize(n);
261 std::fill(buffer[i].begin(), buffer[i].end(), 0.0);
262 }
263 }
264
270 BufferType bufferType = DofBuffer)
271 {
272 size_t n;
273 if (bufferType == VertexBuffer)
274 n = static_cast<size_t>(simulator_.gridView().size(dim));
275 else if (bufferType == ElementBuffer)
276 n = static_cast<size_t>(simulator_.gridView().size(0));
277 else if (bufferType == DofBuffer)
278 n = simulator_.model().numGridDof();
279 else
280 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
281
282 for (unsigned i = 0; i < numComponents; ++i) {
283 buffer[i].resize(n);
284 std::fill(buffer[i].begin(), buffer[i].end(), 0.0);
285 }
286 }
287
293 BufferType bufferType = DofBuffer)
294 {
295 size_t n;
296 if (bufferType == VertexBuffer)
297 n = static_cast<size_t>(simulator_.gridView().size(dim));
298 else if (bufferType == ElementBuffer)
299 n = static_cast<size_t>(simulator_.gridView().size(0));
300 else if (bufferType == DofBuffer)
301 n = simulator_.model().numGridDof();
302 else
303 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
304
305 for (unsigned i = 0; i < numPhases; ++i) {
306 for (unsigned j = 0; j < numComponents; ++j) {
307 buffer[i][j].resize(n);
308 std::fill(buffer[i][j].begin(), buffer[i][j].end(), 0.0);
309 }
310 }
311 }
312
317 const char *name,
318 ScalarBuffer& buffer,
319 BufferType bufferType = DofBuffer)
320 {
321 if (bufferType == DofBuffer)
322 DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer, name);
323 else if (bufferType == VertexBuffer)
324 attachScalarVertexData_(baseWriter, buffer, name);
325 else if (bufferType == ElementBuffer)
326 attachScalarElementData_(baseWriter, buffer, name);
327 else
328 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
329 }
330
335 const char *name,
336 VectorBuffer& buffer,
337 BufferType bufferType = DofBuffer)
338 {
339 if (bufferType == DofBuffer)
340 DiscBaseOutputModule::attachVectorDofData_(baseWriter, buffer, name);
341 else if (bufferType == VertexBuffer)
342 attachVectorVertexData_(baseWriter, buffer, name);
343 else if (bufferType == ElementBuffer)
344 attachVectorElementData_(baseWriter, buffer, name);
345 else
346 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
347 }
348
353 const char *name,
354 TensorBuffer& buffer,
355 BufferType bufferType = DofBuffer)
356 {
357 if (bufferType == DofBuffer)
358 DiscBaseOutputModule::attachTensorDofData_(baseWriter, buffer, name);
359 else if (bufferType == VertexBuffer)
360 attachTensorVertexData_(baseWriter, buffer, name);
361 else if (bufferType == ElementBuffer)
362 attachTensorElementData_(baseWriter, buffer, name);
363 else
364 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
365 }
366
371 const char *pattern,
372 EqBuffer& buffer,
373 BufferType bufferType = DofBuffer)
374 {
375 char name[512];
376 for (unsigned i = 0; i < numEq; ++i) {
377 std::string eqName = simulator_.model().primaryVarName(i);
378 snprintf(name, 512, pattern, eqName.c_str());
379
380 if (bufferType == DofBuffer)
381 DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer[i], name);
382 else if (bufferType == VertexBuffer)
383 attachScalarVertexData_(baseWriter, buffer[i], name);
384 else if (bufferType == ElementBuffer)
385 attachScalarElementData_(baseWriter, buffer[i], name);
386 else
387 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
388 }
389 }
390
395 const char *pattern,
396 EqBuffer& buffer,
397 BufferType bufferType = DofBuffer)
398 {
399 char name[512];
400 for (unsigned i = 0; i < numEq; ++i) {
401 std::ostringstream oss;
402 oss << i;
403 snprintf(name, 512, pattern, oss.str().c_str());
404
405 if (bufferType == DofBuffer)
406 DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer[i], name);
407 else if (bufferType == VertexBuffer)
408 attachScalarVertexData_(baseWriter, buffer[i], name);
409 else if (bufferType == ElementBuffer)
410 attachScalarElementData_(baseWriter, buffer[i], name);
411 else
412 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
413 }
414 }
415
420 const char *pattern,
421 PhaseBuffer& buffer,
422 BufferType bufferType = DofBuffer)
423 {
424 char name[512];
425 for (unsigned i = 0; i < numPhases; ++i) {
426 snprintf(name, 512, pattern, FluidSystem::phaseName(i).data());
427
428 if (bufferType == DofBuffer)
429 DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer[i], name);
430 else if (bufferType == VertexBuffer)
431 attachScalarVertexData_(baseWriter, buffer[i], name);
432 else if (bufferType == ElementBuffer)
433 attachScalarElementData_(baseWriter, buffer[i], name);
434 else
435 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
436 }
437 }
438
443 const char *pattern,
444 ComponentBuffer& buffer,
445 BufferType bufferType = DofBuffer)
446 {
447 char name[512];
448 for (unsigned i = 0; i < numComponents; ++i) {
449 snprintf(name, 512, pattern, FluidSystem::componentName(i).data());
450
451 if (bufferType == DofBuffer)
452 DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer[i], name);
453 else if (bufferType == VertexBuffer)
454 attachScalarVertexData_(baseWriter, buffer[i], name);
455 else if (bufferType == ElementBuffer)
456 attachScalarElementData_(baseWriter, buffer[i], name);
457 else
458 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
459 }
460 }
461
466 const char *pattern,
467 PhaseComponentBuffer& buffer,
468 BufferType bufferType = DofBuffer)
469 {
470 char name[512];
471 for (unsigned i= 0; i < numPhases; ++i) {
472 for (unsigned j = 0; j < numComponents; ++j) {
473 snprintf(name, 512, pattern,
474 FluidSystem::phaseName(i).data(),
475 FluidSystem::componentName(j).data());
476
477 if (bufferType == DofBuffer)
478 DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer[i][j], name);
479 else if (bufferType == VertexBuffer)
480 attachScalarVertexData_(baseWriter, buffer[i][j], name);
481 else if (bufferType == ElementBuffer)
482 attachScalarElementData_(baseWriter, buffer[i][j], name);
483 else
484 throw std::logic_error("bufferType must be one of Dof, Vertex or Element");
485 }
486 }
487 }
488
490 ScalarBuffer& buffer,
491 const char *name)
492 { baseWriter.attachScalarElementData(buffer, name); }
493
495 ScalarBuffer& buffer,
496 const char *name)
497 { baseWriter.attachScalarVertexData(buffer, name); }
498
500 VectorBuffer& buffer,
501 const char *name)
502 { baseWriter.attachVectorElementData(buffer, name); }
503
505 VectorBuffer& buffer,
506 const char *name)
507 { baseWriter.attachVectorVertexData(buffer, name); }
508
510 TensorBuffer& buffer,
511 const char *name)
512 { baseWriter.attachTensorElementData(buffer, name); }
513
515 TensorBuffer& buffer,
516 const char *name)
517 { baseWriter.attachTensorVertexData(buffer, name); }
518
519 const Simulator& simulator_;
520};
521
522#if __GNUC__ || __clang__
523#pragma GCC diagnostic pop
524#endif
525
526} // namespace Opm
527
528#endif
Defines a type tags and some fundamental properties all models.
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
void commitComponentBuffer_(BaseOutputWriter &baseWriter, const char *pattern, ComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Add a component-specific buffer to the result file.
Definition: baseoutputmodule.hh:442
void attachVectorElementData_(BaseOutputWriter &baseWriter, VectorBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:499
virtual bool needExtensiveQuantities() const
Returns true iff the module needs to access the extensive quantities of a context to do its job.
Definition: baseoutputmodule.hh:144
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, numComponents > ComponentBuffer
Definition: baseoutputmodule.hh:97
std::array< ScalarBuffer, numPhases > PhaseBuffer
Definition: baseoutputmodule.hh:96
void resizeVectorBuffer_(VectorBuffer &buffer, BufferType bufferType=DofBuffer)
Definition: baseoutputmodule.hh:200
void resizeComponentBuffer_(ComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a component specific quantity.
Definition: baseoutputmodule.hh:269
void attachScalarElementData_(BaseOutputWriter &baseWriter, ScalarBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:489
void attachVectorVertexData_(BaseOutputWriter &baseWriter, VectorBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:504
void commitEqBuffer_(BaseOutputWriter &baseWriter, const char *pattern, EqBuffer &buffer, BufferType bufferType=DofBuffer)
Add a buffer with as many variables as PDEs to the result file.
Definition: baseoutputmodule.hh:394
void resizePhaseComponentBuffer_(PhaseComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase and component specific buffer.
Definition: baseoutputmodule.hh:292
BufferType
Definition: baseoutputmodule.hh:148
@ DofBuffer
Buffer contains data associated with the degrees of freedom.
Definition: baseoutputmodule.hh:150
@ VertexBuffer
Buffer contains data associated with the grid's vertices.
Definition: baseoutputmodule.hh:153
@ ElementBuffer
Buffer contains data associated with the grid's elements.
Definition: baseoutputmodule.hh:156
void commitPriVarsBuffer_(BaseOutputWriter &baseWriter, const char *pattern, EqBuffer &buffer, BufferType bufferType=DofBuffer)
Add a buffer with as many variables as PDEs to the result file.
Definition: baseoutputmodule.hh:370
void resizeEqBuffer_(EqBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a equation specific quantity.
Definition: baseoutputmodule.hh:223
void commitPhaseComponentBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Add a phase and component specific quantities to the output.
Definition: baseoutputmodule.hh:465
std::array< std::array< ScalarBuffer, numComponents >, numPhases > PhaseComponentBuffer
Definition: baseoutputmodule.hh:98
BaseOutputWriter::VectorBuffer VectorBuffer
Definition: baseoutputmodule.hh:92
virtual void processElement(const ElementContext &elemCtx)=0
Modify the internal buffers according to the intensive quanties relevant for an element.
BaseOutputModule(const Simulator &simulator)
Definition: baseoutputmodule.hh:102
void resizeTensorBuffer_(TensorBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a tensorial quantity.
Definition: baseoutputmodule.hh:182
void attachScalarVertexData_(BaseOutputWriter &baseWriter, ScalarBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:494
virtual void commitBuffers(BaseOutputWriter &writer)=0
Add all buffers to the VTK output writer.
void attachTensorElementData_(BaseOutputWriter &baseWriter, TensorBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:509
void commitVectorBuffer_(BaseOutputWriter &baseWriter, const char *name, VectorBuffer &buffer, BufferType bufferType=DofBuffer)
Add a buffer containing vectorial quantities to the result file.
Definition: baseoutputmodule.hh:334
virtual ~BaseOutputModule()
Definition: baseoutputmodule.hh:106
virtual void allocBuffers()=0
Allocate memory for the scalar fields we would like to write to disk.
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 attachTensorVertexData_(BaseOutputWriter &baseWriter, TensorBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:514
std::array< ScalarBuffer, numEq > EqBuffer
Definition: baseoutputmodule.hh:95
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
virtual void attachScalarVertexData(ScalarBuffer &buf, std::string name)=0
Add a scalar vertex centered vector field to the output.
virtual void attachVectorVertexData(VectorBuffer &buf, std::string name)=0
Add a vectorial vertex centered vector field to the output.
Dune::DynamicVector< double > Vector
Definition: baseoutputwriter.hh:47
std::vector< Tensor > TensorBuffer
Definition: baseoutputwriter.hh:51
virtual void attachTensorElementData(TensorBuffer &buf, std::string name)=0
Add a tensorial element centered quantity to the output.
Dune::DynamicMatrix< double > Tensor
Definition: baseoutputwriter.hh:48
virtual void attachTensorVertexData(TensorBuffer &buf, std::string name)=0
Add a tensorial vertex centered tensor field to the output.
std::vector< Vector > VectorBuffer
Definition: baseoutputwriter.hh:50
virtual void attachScalarElementData(ScalarBuffer &buf, std::string name)=0
Add a scalar element centered quantity to the output.
std::vector< Scalar > ScalarBuffer
Definition: baseoutputwriter.hh:49
virtual void attachVectorElementData(VectorBuffer &buf, std::string name)=0
Add a vectorial element centered quantity to the output.
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
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.