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  Copyright (C) 2011-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef EWOMS_BASE_OUTPUT_MODULE_HH
26 #define EWOMS_BASE_OUTPUT_MODULE_HH
27 
28 #include "baseoutputwriter.hh"
29 
31 
33 #include <opm/common/ErrorMacros.hpp>
34 
35 #include <dune/istl/bvector.hh>
36 #include <dune/common/fvector.hh>
37 
38 #include <vector>
39 #include <sstream>
40 #include <string>
41 #include <array>
42 
43 #include <cstdio>
44 
45 
46 namespace Ewoms {
47 namespace Properties {
48 // forward definition of property tags
49 NEW_PROP_TAG(NumPhases);
50 NEW_PROP_TAG(NumComponents);
51 NEW_PROP_TAG(NumEq);
52 
53 NEW_PROP_TAG(Model);
54 NEW_PROP_TAG(Simulator);
55 NEW_PROP_TAG(Scalar);
56 NEW_PROP_TAG(Evaluation);
57 NEW_PROP_TAG(GridView);
58 NEW_PROP_TAG(ElementContext);
59 NEW_PROP_TAG(FluidSystem);
60 NEW_PROP_TAG(DiscBaseOutputModule);
61 } // namespace Properties
62 
70 template<class TypeTag>
72 {
73  typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
74  typedef typename GET_PROP_TYPE(TypeTag, Model) Model;
75  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
76  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
77  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
78  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
79  typedef typename GET_PROP_TYPE(TypeTag, DiscBaseOutputModule) DiscBaseOutputModule;
80 
81  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
82  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
83  enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
84  enum { dim = GridView::dimension };
85  enum { dimWorld = GridView::dimensionworld };
86 
87  typedef BaseOutputWriter::Tensor Tensor;
88 
89 public:
93 
94  typedef std::array<ScalarBuffer, numEq> EqBuffer;
95  typedef std::array<ScalarBuffer, numPhases> PhaseBuffer;
96  typedef std::array<ScalarBuffer, numComponents> ComponentBuffer;
97  typedef std::array<std::array<ScalarBuffer, numComponents>, numPhases> PhaseComponentBuffer;
98 
99  typedef std::array<VectorBuffer, numPhases> PhaseVectorBuffer;
100 
101  BaseOutputModule(const Simulator &simulator)
102  : simulator_(simulator)
103  {}
104 
106  {}
107 
116  virtual void allocBuffers() = 0;
117 
126  virtual void processElement(const ElementContext &elemCtx) = 0;
127 
131  virtual void commitBuffers(BaseOutputWriter &writer) = 0;
132 
133 protected:
134  enum BufferType {
137 
140 
143  };
144 
148  void resizeScalarBuffer_(ScalarBuffer &buffer,
149  BufferType bufferType = DofBuffer)
150  {
151  Scalar n;
152  if (bufferType == VertexBuffer)
153  n = simulator_.gridView().size(dim);
154  else if (bufferType == ElementBuffer)
155  n = simulator_.gridView().size(0);
156  else if (bufferType == DofBuffer)
157  n = simulator_.model().numGridDof();
158  else
159  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
160 
161  buffer.resize(n);
162  std::fill(buffer.begin(), buffer.end(), 0.0);
163  }
164 
168  void resizeTensorBuffer_(TensorBuffer &buffer,
169  BufferType bufferType = DofBuffer)
170  {
171  Scalar n;
172  if (bufferType == VertexBuffer)
173  n = simulator_.gridView().size(dim);
174  else if (bufferType == ElementBuffer)
175  n = simulator_.gridView().size(0);
176  else if (bufferType == DofBuffer)
177  n = simulator_.model().numGridDof();
178  else
179  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
180 
181  buffer.resize(n);
182  Tensor nullMatrix(dimWorld, dimWorld, 0.0);
183  std::fill(buffer.begin(), buffer.end(), nullMatrix);
184  }
185 
190  void resizeEqBuffer_(EqBuffer &buffer,
191  BufferType bufferType = DofBuffer)
192  {
193  Scalar n;
194  if (bufferType == VertexBuffer)
195  n = simulator_.gridView().size(dim);
196  else if (bufferType == ElementBuffer)
197  n = simulator_.gridView().size(0);
198  else if (bufferType == DofBuffer)
199  n = simulator_.model().numGridDof();
200  else
201  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
202 
203  for (int i = 0; i < numEq; ++i) {
204  buffer[i].resize(n);
205  std::fill(buffer[i].begin(), buffer[i].end(), 0.0);
206  }
207  }
208 
213  void resizePhaseBuffer_(PhaseBuffer &buffer,
214  BufferType bufferType = DofBuffer)
215  {
216  Scalar n;
217  if (bufferType == VertexBuffer)
218  n = simulator_.gridView().size(dim);
219  else if (bufferType == ElementBuffer)
220  n = simulator_.gridView().size(0);
221  else if (bufferType == DofBuffer)
222  n = simulator_.model().numGridDof();
223  else
224  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
225 
226  for (int i = 0; i < numPhases; ++i) {
227  buffer[i].resize(n);
228  std::fill(buffer[i].begin(), buffer[i].end(), 0.0);
229  }
230  }
231 
236  void resizeComponentBuffer_(ComponentBuffer &buffer,
237  BufferType bufferType = DofBuffer)
238  {
239  Scalar n;
240  if (bufferType == VertexBuffer)
241  n = simulator_.gridView().size(dim);
242  else if (bufferType == ElementBuffer)
243  n = simulator_.gridView().size(0);
244  else if (bufferType == DofBuffer)
245  n = simulator_.model().numGridDof();
246  else
247  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
248 
249  for (int i = 0; i < numComponents; ++i) {
250  buffer[i].resize(n);
251  std::fill(buffer[i].begin(), buffer[i].end(), 0.0);
252  }
253  }
254 
259  void resizePhaseComponentBuffer_(PhaseComponentBuffer &buffer,
260  BufferType bufferType = DofBuffer)
261  {
262  Scalar n;
263  if (bufferType == VertexBuffer)
264  n = simulator_.gridView().size(dim);
265  else if (bufferType == ElementBuffer)
266  n = simulator_.gridView().size(0);
267  else if (bufferType == DofBuffer)
268  n = simulator_.model().numGridDof();
269  else
270  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
271 
272  for (int i = 0; i < numPhases; ++i) {
273  for (int j = 0; j < numComponents; ++j) {
274  buffer[i][j].resize(n);
275  std::fill(buffer[i][j].begin(), buffer[i][j].end(), 0.0);
276  }
277  }
278  }
279 
284  const char *name,
285  ScalarBuffer &buffer,
286  BufferType bufferType = DofBuffer)
287  {
288  if (bufferType == DofBuffer)
289  DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer, name);
290  else if (bufferType == VertexBuffer)
291  attachScalarVertexData_(baseWriter, buffer, name);
292  else if (bufferType == ElementBuffer)
293  attachScalarElementData_(baseWriter, buffer, name);
294  else
295  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
296  }
297 
302  const char *name,
303  VectorBuffer &buffer,
304  BufferType bufferType = DofBuffer)
305  {
306  if (bufferType == DofBuffer)
307  DiscBaseOutputModule::attachVectorDofData_(baseWriter, buffer, name);
308  else if (bufferType == VertexBuffer)
309  attachVectorVertexData_(baseWriter, buffer, name);
310  else if (bufferType == ElementBuffer)
311  attachVectorElementData_(baseWriter, buffer, name);
312  else
313  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
314  }
315 
320  const char *name,
321  TensorBuffer &buffer,
322  BufferType bufferType = DofBuffer)
323  {
324  if (bufferType == DofBuffer)
325  DiscBaseOutputModule::attachTensorDofData_(baseWriter, buffer, name);
326  else if (bufferType == VertexBuffer)
327  attachTensorVertexData_(baseWriter, buffer, name);
328  else if (bufferType == ElementBuffer)
329  attachTensorElementData_(baseWriter, buffer, name);
330  else
331  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
332  }
333 
338  const char *pattern,
339  EqBuffer &buffer,
340  BufferType bufferType = DofBuffer)
341  {
342  char name[512];
343  for (int i = 0; i < numEq; ++i) {
344  std::string eqName = simulator_.model().primaryVarName(i);
345  snprintf(name, 512, pattern, eqName.c_str());
346 
347  if (bufferType == DofBuffer)
348  DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer[i], name);
349  else if (bufferType == VertexBuffer)
350  attachScalarVertexData_(baseWriter, buffer[i], name);
351  else if (bufferType == ElementBuffer)
352  attachScalarElementData_(baseWriter, buffer[i], name);
353  else
354  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
355  }
356  }
357 
362  const char *pattern,
363  EqBuffer &buffer,
364  BufferType bufferType = DofBuffer)
365  {
366  char name[512];
367  for (int i = 0; i < numEq; ++i) {
368  std::ostringstream oss;
369  oss << i;
370  snprintf(name, 512, pattern, oss.str().c_str());
371 
372  if (bufferType == DofBuffer)
373  DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer[i], name);
374  else if (bufferType == VertexBuffer)
375  attachScalarVertexData_(baseWriter, buffer[i], name);
376  else if (bufferType == ElementBuffer)
377  attachScalarElementData_(baseWriter, buffer[i], name);
378  else
379  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
380  }
381  }
382 
387  const char *pattern,
388  PhaseBuffer &buffer,
389  BufferType bufferType = DofBuffer)
390  {
391  char name[512];
392  for (int i = 0; i < numPhases; ++i) {
393  snprintf(name, 512, pattern, FluidSystem::phaseName(i));
394 
395  if (bufferType == DofBuffer)
396  DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer[i], name);
397  else if (bufferType == VertexBuffer)
398  attachScalarVertexData_(baseWriter, buffer[i], name);
399  else if (bufferType == ElementBuffer)
400  attachScalarElementData_(baseWriter, buffer[i], name);
401  else
402  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
403  }
404  }
405 
410  const char *pattern,
411  ComponentBuffer &buffer,
412  BufferType bufferType = DofBuffer)
413  {
414  char name[512];
415  for (int i = 0; i < numComponents; ++i) {
416  snprintf(name, 512, pattern, FluidSystem::componentName(i));
417 
418  if (bufferType == DofBuffer)
419  DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer[i], name);
420  else if (bufferType == VertexBuffer)
421  attachScalarVertexData_(baseWriter, buffer[i], name);
422  else if (bufferType == ElementBuffer)
423  attachScalarElementData_(baseWriter, buffer[i], name);
424  else
425  OPM_THROW(std::logic_error, "bufferType must be one of Dof, Vertex or Element");
426  }
427  }
428 
433  const char *pattern,
434  PhaseComponentBuffer &buffer,
435  BufferType bufferType = DofBuffer)
436  {
437  char name[512];
438  for (int i= 0; i < numPhases; ++i) {
439  for (int j = 0; j < numComponents; ++j) {
440  snprintf(name, 512, pattern,
441  FluidSystem::phaseName(i),
442  FluidSystem::componentName(j));
443 
444  if (bufferType == DofBuffer)
445  DiscBaseOutputModule::attachScalarDofData_(baseWriter, buffer[i][j], name);
446  else if (bufferType == VertexBuffer)
447  attachScalarVertexData_(baseWriter, buffer[i][j], name);
448  else if (bufferType == ElementBuffer)
449  attachScalarElementData_(baseWriter, buffer[i][j], name);
450  else
451  OPM_THROW(std::logic_error,
452  "bufferType must be one of Dof, Vertex or Element");
453  }
454  }
455  }
456 
458  ScalarBuffer &buffer,
459  const char *name)
460  { baseWriter.attachScalarElementData(buffer, name); }
461 
463  ScalarBuffer &buffer,
464  const char *name)
465  { baseWriter.attachScalarVertexData(buffer, name); }
466 
468  VectorBuffer &buffer,
469  const char *name)
470  { baseWriter.attachVectorElementData(buffer, name); }
471 
473  VectorBuffer &buffer,
474  const char *name)
475  { baseWriter.attachVectorVertexData(buffer, name); }
476 
478  TensorBuffer &buffer,
479  const char *name)
480  { baseWriter.attachTensorElementData(buffer, name); }
481 
483  TensorBuffer &buffer,
484  const char *name)
485  { baseWriter.attachTensorVertexData(buffer, name); }
486 
487  const Simulator &simulator_;
488 };
489 
490 } // namespace Ewoms
491 
492 #endif
Buffer contains data associated with the grid's elements.
Definition: baseoutputmodule.hh:142
virtual void attachVectorElementData(VectorBuffer &buf, std::string name)=0
Add a vectorial element centered quantity to the output.
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:319
void attachVectorVertexData_(BaseOutputWriter &baseWriter, VectorBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:472
std::array< ScalarBuffer, numComponents > ComponentBuffer
Definition: baseoutputmodule.hh:96
void resizeTensorBuffer_(TensorBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a tensorial quantity.
Definition: baseoutputmodule.hh:168
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a scalar quantity.
Definition: baseoutputmodule.hh:148
void resizePhaseBuffer_(PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase-specific quantity.
Definition: baseoutputmodule.hh:213
std::vector< Vector > VectorBuffer
Definition: baseoutputwriter.hh:48
void resizePhaseComponentBuffer_(PhaseComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase and component specific buffer.
Definition: baseoutputmodule.hh:259
void resizeComponentBuffer_(ComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a component specific quantity.
Definition: baseoutputmodule.hh:236
std::vector< Scalar > ScalarBuffer
Definition: baseoutputwriter.hh:47
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:361
virtual void attachVectorVertexData(VectorBuffer &buf, std::string name)=0
Add a vectorial vertex centered vector field to the output.
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
void attachScalarVertexData_(BaseOutputWriter &baseWriter, ScalarBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:462
std::vector< Tensor > TensorBuffer
Definition: baseoutputwriter.hh:49
virtual void attachTensorVertexData(TensorBuffer &buf, std::string name)=0
Add a tensorial vertex centered tensor field to the output.
BaseOutputWriter::TensorBuffer TensorBuffer
Definition: baseoutputmodule.hh:92
void attachVectorElementData_(BaseOutputWriter &baseWriter, VectorBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:467
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:337
std::array< std::array< ScalarBuffer, numComponents >, numPhases > PhaseComponentBuffer
Definition: baseoutputmodule.hh:97
std::array< VectorBuffer, numPhases > PhaseVectorBuffer
Definition: baseoutputmodule.hh:99
BufferType
Definition: baseoutputmodule.hh:134
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:432
NEW_PROP_TAG(Grid)
The type of the DUNE grid.
This file provides the infrastructure to retrieve run-time parameters.
void attachTensorElementData_(BaseOutputWriter &baseWriter, TensorBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:477
void commitComponentBuffer_(BaseOutputWriter &baseWriter, const char *pattern, ComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Add a component-specific buffer to the result file.
Definition: baseoutputmodule.hh:409
The base class for all output writers.
Definition: baseoutputwriter.hh:41
virtual void processElement(const ElementContext &elemCtx)=0
Modify the internal buffers according to the intensive quanties relevant for an element.
BaseOutputWriter::VectorBuffer VectorBuffer
Definition: baseoutputmodule.hh:91
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:73
Definition: baseauxiliarymodule.hh:35
void resizeEqBuffer_(EqBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a equation specific quantity.
Definition: baseoutputmodule.hh:190
BaseOutputModule(const Simulator &simulator)
Definition: baseoutputmodule.hh:101
The base class for all output writers.
Model & model()
Return the physical model used in the simulation.
Definition: simulator.hh:176
virtual ~BaseOutputModule()
Definition: baseoutputmodule.hh:105
std::array< ScalarBuffer, numPhases > PhaseBuffer
Definition: baseoutputmodule.hh:95
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:283
virtual void attachTensorElementData(TensorBuffer &buf, std::string name)=0
Add a tensorial element centered quantity to the output.
virtual void allocBuffers()=0
Allocate memory for the scalar fields we would like to write to disk.
GridView & gridView()
Return the grid view for which the simulation is done.
Definition: simulator.hh:164
void attachScalarElementData_(BaseOutputWriter &baseWriter, ScalarBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:457
virtual void attachScalarVertexData(ScalarBuffer &buf, std::string name)=0
Add a scalar vertex centered vector field to the output.
virtual void attachScalarElementData(ScalarBuffer &buf, std::string name)=0
Add a scalar element centered quantity to the output.
Provides the magic behind the eWoms property system.
Buffer contains data associated with the grid's vertices.
Definition: baseoutputmodule.hh:139
virtual void commitBuffers(BaseOutputWriter &writer)=0
Add all buffers to the VTK output writer.
void commitPhaseBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Add a phase-specific buffer to the result file.
Definition: baseoutputmodule.hh:386
Dune::DynamicMatrix< double > Tensor
Definition: baseoutputwriter.hh:46
void attachTensorVertexData_(BaseOutputWriter &baseWriter, TensorBuffer &buffer, const char *name)
Definition: baseoutputmodule.hh:482
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:301
The base class for writer modules.
Definition: baseoutputmodule.hh:71
const Simulator & simulator_
Definition: baseoutputmodule.hh:487
std::array< ScalarBuffer, numEq > EqBuffer
Definition: baseoutputmodule.hh:94
Buffer contains data associated with the degrees of freedom.
Definition: baseoutputmodule.hh:136
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:90