opm-simulators/ebos/vtkecltracermodule.hh
Bernd Flemisch 21df1cbe31 [properties] adapt to changes in the property system
`NEW_PROP_TAG` is now a definition and not just a declaration.
Eliminate superfluous declarations, include headers with definitions.
Make one necessary forward declaration explicit.
2020-05-18 15:54:26 +02:00

169 lines
5.5 KiB
C++

// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 2 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
Consult the COPYING file in the top-level source directory of this
module for the precise wording of the license and the list of
copyright holders.
*/
/*!
* \file
* \copydoc Opm::VtkEclTracerModule
*/
#ifndef EWOMS_VTK_ECL_TRACER_MODULE_HH
#define EWOMS_VTK_ECL_TRACER_MODULE_HH
#include <opm/models/io/vtkmultiwriter.hh>
#include <opm/models/io/baseoutputmodule.hh>
#include <opm/models/utils/propertysystem.hh>
#include <opm/models/utils/parametersystem.hh>
#include <opm/models/blackoil/blackoilproperties.hh>
#include <dune/common/fvector.hh>
#include <cstdio>
BEGIN_PROPERTIES
// create new type tag for the VTK tracer output
NEW_TYPE_TAG(VtkEclTracer);
// create the property tags needed for the tracer model
NEW_PROP_TAG(VtkWriteEclTracerConcentration);
// set default values for what quantities to output
SET_BOOL_PROP(VtkEclTracer, VtkWriteEclTracerConcentration, false);
END_PROPERTIES
namespace Opm {
/*!
* \ingroup Vtk
*
* \brief VTK output module for the tracer model's parameters.
*/
template <class TypeTag>
class VtkEclTracerModule : public BaseOutputModule<TypeTag>
{
typedef BaseOutputModule<TypeTag> ParentType;
typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
static const int vtkFormat = GET_PROP_VALUE(TypeTag, VtkOutputFormat);
typedef Opm::VtkMultiWriter<GridView, vtkFormat> VtkMultiWriter;
typedef typename ParentType::ScalarBuffer ScalarBuffer;
public:
VtkEclTracerModule(const Simulator& simulator)
: ParentType(simulator)
{ }
/*!
* \brief Register all run-time parameters for the tracer VTK output
* module.
*/
static void registerParameters()
{
EWOMS_REGISTER_PARAM(TypeTag, bool, VtkWriteEclTracerConcentration,
"Include the tracer concentration "
"in the VTK output files");
}
/*!
* \brief Allocate memory for the scalar fields we would like to
* write to the VTK file.
*/
void allocBuffers()
{
if (eclTracerConcentrationOutput_()){
const auto& tracerModel = this->simulator_.problem().tracerModel();
eclTracerConcentration_.resize(tracerModel.numTracers());
for(size_t tracerIdx=0; tracerIdx<eclTracerConcentration_.size();++tracerIdx){
this->resizeScalarBuffer_(eclTracerConcentration_[tracerIdx]);
}
}
}
/*!
* \brief Modify the internal buffers according to the intensive quantities relevant for
* an element
*/
void processElement(const ElementContext& elemCtx)
{
if (!EWOMS_GET_PARAM(TypeTag, bool, EnableVtkOutput))
return;
const auto& tracerModel = elemCtx.problem().tracerModel();
for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
if (eclTracerConcentrationOutput_()){
for(size_t tracerIdx=0; tracerIdx<eclTracerConcentration_.size();++tracerIdx){
eclTracerConcentration_[tracerIdx][globalDofIdx] = tracerModel.tracerConcentration(tracerIdx, globalDofIdx);
}
}
}
}
/*!
* \brief Add all buffers to the VTK output writer.
*/
void commitBuffers(BaseOutputWriter& baseWriter)
{
VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
if (!vtkWriter)
return;
if (eclTracerConcentrationOutput_()){
const auto& tracerModel = this->simulator_.problem().tracerModel();
for(size_t tracerIdx=0; tracerIdx<eclTracerConcentration_.size();++tracerIdx){
const std::string tmp = "tracerConcentration_" + tracerModel.tracerName(tracerIdx);
this->commitScalarBuffer_(baseWriter,tmp.c_str(), eclTracerConcentration_[tracerIdx]);
}
}
}
private:
static bool eclTracerConcentrationOutput_()
{
static bool val = EWOMS_GET_PARAM(TypeTag, bool, VtkWriteEclTracerConcentration);
return val;
}
std::vector<ScalarBuffer> eclTracerConcentration_;
};
} // namespace Opm
#endif