// -*- 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 . 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 #include #include #include #include #include #include namespace Opm::Properties { // create new type tag for the VTK tracer output namespace TTag { struct VtkEclTracer {}; } // create the property tags needed for the tracer model template struct VtkWriteEclTracerConcentration { using type = UndefinedProperty; }; // set default values for what quantities to output template struct VtkWriteEclTracerConcentration { static constexpr bool value = false; }; } // namespace Opm::Properties namespace Opm { /*! * \ingroup Vtk * * \brief VTK output module for the tracer model's parameters. */ template class VtkEclTracerModule : public BaseOutputModule { using ParentType = BaseOutputModule; using Simulator = GetPropType; using Scalar = GetPropType; using Evaluation = GetPropType; using ElementContext = GetPropType; using GridView = GetPropType; using FluidSystem = GetPropType; static constexpr int vtkFormat = getPropValue(); using VtkMultiWriter = ::Opm::VtkMultiWriter; using ScalarBuffer = typename ParentType::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; tracerIdxresizeScalarBuffer_(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(&baseWriter); if (!vtkWriter) return; if (eclTracerConcentrationOutput_()){ const auto& tracerModel = this->simulator_.problem().tracerModel(); for(size_t tracerIdx=0; tracerIdxcommitScalarBuffer_(baseWriter,tmp.c_str(), eclTracerConcentration_[tracerIdx]); } } } private: static bool eclTracerConcentrationOutput_() { static bool val = EWOMS_GET_PARAM(TypeTag, bool, VtkWriteEclTracerConcentration); return val; } std::vector eclTracerConcentration_; }; } // namespace Opm #endif