mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
vtkblackoilmicpmodule: move parameters to dedicated struct with a translation unit
This commit is contained in:
parent
28bbe60b6e
commit
ce4e57dab3
@ -63,6 +63,7 @@ list (APPEND MAIN_SOURCE_FILES
|
|||||||
opm/models/blackoil/blackoilpolymerparams.cpp
|
opm/models/blackoil/blackoilpolymerparams.cpp
|
||||||
opm/models/blackoil/blackoilsolventparams.cpp
|
opm/models/blackoil/blackoilsolventparams.cpp
|
||||||
opm/models/io/vtkblackoilenergyparams.cpp
|
opm/models/io/vtkblackoilenergyparams.cpp
|
||||||
|
opm/models/io/vtkblackoilmicpparams.cpp
|
||||||
opm/models/io/vtkblackoilparams.cpp
|
opm/models/io/vtkblackoilparams.cpp
|
||||||
opm/models/io/restart.cpp
|
opm/models/io/restart.cpp
|
||||||
opm/models/parallel/mpiutil.cpp
|
opm/models/parallel/mpiutil.cpp
|
||||||
@ -654,6 +655,7 @@ list (APPEND PUBLIC_HEADER_FILES
|
|||||||
opm/models/io/vtkblackoilenergymodule.hpp
|
opm/models/io/vtkblackoilenergymodule.hpp
|
||||||
opm/models/io/vtkblackoilenergyparams.hpp
|
opm/models/io/vtkblackoilenergyparams.hpp
|
||||||
opm/models/io/vtkblackoilmicpmodule.hpp
|
opm/models/io/vtkblackoilmicpmodule.hpp
|
||||||
|
opm/models/io/vtkblackoilmicpparams.hpp
|
||||||
opm/models/io/vtkblackoilmodule.hpp
|
opm/models/io/vtkblackoilmodule.hpp
|
||||||
opm/models/io/vtkblackoilparams.hpp
|
opm/models/io/vtkblackoilparams.hpp
|
||||||
opm/models/io/vtkblackoilpolymermodule.hh
|
opm/models/io/vtkblackoilpolymermodule.hh
|
||||||
|
@ -36,22 +36,12 @@
|
|||||||
#include <opm/models/discretization/common/fvbaseparameters.hh>
|
#include <opm/models/discretization/common/fvbaseparameters.hh>
|
||||||
|
|
||||||
#include <opm/models/io/baseoutputmodule.hh>
|
#include <opm/models/io/baseoutputmodule.hh>
|
||||||
|
#include <opm/models/io/vtkblackoilmicpparams.hpp>
|
||||||
#include <opm/models/io/vtkmultiwriter.hh>
|
#include <opm/models/io/vtkmultiwriter.hh>
|
||||||
|
|
||||||
#include <opm/models/utils/parametersystem.hpp>
|
#include <opm/models/utils/parametersystem.hpp>
|
||||||
#include <opm/models/utils/propertysystem.hh>
|
#include <opm/models/utils/propertysystem.hh>
|
||||||
|
|
||||||
namespace Opm::Parameters {
|
|
||||||
|
|
||||||
// set default values for what quantities to output
|
|
||||||
struct VtkWriteMicrobialConcentration { static constexpr bool value = true; };
|
|
||||||
struct VtkWriteOxygenConcentration { static constexpr bool value = true; };
|
|
||||||
struct VtkWriteUreaConcentration { static constexpr bool value = true; };
|
|
||||||
struct VtkWriteBiofilmConcentration { static constexpr bool value = true; };
|
|
||||||
struct VtkWriteCalciteConcentration { static constexpr bool value = true; };
|
|
||||||
|
|
||||||
} // namespace Opm::Parameters
|
|
||||||
|
|
||||||
namespace Opm {
|
namespace Opm {
|
||||||
/*!
|
/*!
|
||||||
* \ingroup Vtk
|
* \ingroup Vtk
|
||||||
@ -79,7 +69,11 @@ class VtkBlackOilMICPModule : public BaseOutputModule<TypeTag>
|
|||||||
public:
|
public:
|
||||||
VtkBlackOilMICPModule(const Simulator& simulator)
|
VtkBlackOilMICPModule(const Simulator& simulator)
|
||||||
: ParentType(simulator)
|
: ParentType(simulator)
|
||||||
{ }
|
{
|
||||||
|
if constexpr (enableMICP) {
|
||||||
|
params_.read();
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
/*!
|
/*!
|
||||||
* \brief Register all run-time parameters for the multi-phase VTK output
|
* \brief Register all run-time parameters for the multi-phase VTK output
|
||||||
@ -87,22 +81,9 @@ public:
|
|||||||
*/
|
*/
|
||||||
static void registerParameters()
|
static void registerParameters()
|
||||||
{
|
{
|
||||||
if (!enableMICP)
|
if constexpr (enableMICP) {
|
||||||
return;
|
VtkBlackoilMICPParams::registerParameters();
|
||||||
|
}
|
||||||
Parameters::Register<Parameters::VtkWriteMicrobialConcentration>
|
|
||||||
("Include the concentration of the microbial component in the water phase "
|
|
||||||
"in the VTK output files");
|
|
||||||
Parameters::Register<Parameters::VtkWriteOxygenConcentration>
|
|
||||||
("Include the concentration of the oxygen component in the water phase "
|
|
||||||
"in the VTK output files");
|
|
||||||
Parameters::Register<Parameters::VtkWriteUreaConcentration>
|
|
||||||
("Include the concentration of the urea component in the water phase "
|
|
||||||
"in the VTK output files");
|
|
||||||
Parameters::Register<Parameters::VtkWriteBiofilmConcentration>
|
|
||||||
("Include the biofilm volume fraction in the VTK output files");
|
|
||||||
Parameters::Register<Parameters::VtkWriteCalciteConcentration>
|
|
||||||
("Include the calcite volume fraction in the VTK output files");
|
|
||||||
}
|
}
|
||||||
|
|
||||||
/*!
|
/*!
|
||||||
@ -111,23 +92,27 @@ public:
|
|||||||
*/
|
*/
|
||||||
void allocBuffers()
|
void allocBuffers()
|
||||||
{
|
{
|
||||||
if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
|
if constexpr (enableMICP) {
|
||||||
return;
|
if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
|
||||||
|
if (params_.microbialConcentrationOutput_) {
|
||||||
|
this->resizeScalarBuffer_(microbialConcentration_);
|
||||||
|
}
|
||||||
|
if (params_.oxygenConcentrationOutput_) {
|
||||||
|
this->resizeScalarBuffer_(oxygenConcentration_);
|
||||||
|
}
|
||||||
|
if (params_.ureaConcentrationOutput_) {
|
||||||
|
this->resizeScalarBuffer_(ureaConcentration_);
|
||||||
|
}
|
||||||
|
if (params_.biofilmConcentrationOutput_) {
|
||||||
|
this->resizeScalarBuffer_(biofilmConcentration_);
|
||||||
|
}
|
||||||
|
if (params_.calciteConcentrationOutput_) {
|
||||||
|
this->resizeScalarBuffer_(calciteConcentration_);
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
if (!enableMICP)
|
|
||||||
return;
|
|
||||||
|
|
||||||
if (microbialConcentrationOutput_())
|
|
||||||
this->resizeScalarBuffer_(microbialConcentration_);
|
|
||||||
if (oxygenConcentrationOutput_())
|
|
||||||
this->resizeScalarBuffer_(oxygenConcentration_);
|
|
||||||
if (ureaConcentrationOutput_())
|
|
||||||
this->resizeScalarBuffer_(ureaConcentration_);
|
|
||||||
if (biofilmConcentrationOutput_())
|
|
||||||
this->resizeScalarBuffer_(biofilmConcentration_);
|
|
||||||
if (calciteConcentrationOutput_())
|
|
||||||
this->resizeScalarBuffer_(calciteConcentration_);
|
|
||||||
}
|
}
|
||||||
|
|
||||||
/*!
|
/*!
|
||||||
@ -136,37 +121,40 @@ public:
|
|||||||
*/
|
*/
|
||||||
void processElement(const ElementContext& elemCtx)
|
void processElement(const ElementContext& elemCtx)
|
||||||
{
|
{
|
||||||
if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
|
if constexpr (enableMICP) {
|
||||||
return;
|
if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
|
||||||
}
|
return;
|
||||||
|
}
|
||||||
|
|
||||||
if (!enableMICP)
|
for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
|
||||||
return;
|
const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
|
||||||
|
unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
|
||||||
|
|
||||||
for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
|
if (params_.microbialConcentrationOutput_) {
|
||||||
const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
|
microbialConcentration_[globalDofIdx] =
|
||||||
unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
|
scalarValue(intQuants.microbialConcentration());
|
||||||
|
}
|
||||||
|
|
||||||
if (microbialConcentrationOutput_())
|
if (params_.oxygenConcentrationOutput_) {
|
||||||
microbialConcentration_[globalDofIdx] =
|
oxygenConcentration_[globalDofIdx] =
|
||||||
scalarValue(intQuants.microbialConcentration());
|
scalarValue(intQuants.oxygenConcentration());
|
||||||
|
}
|
||||||
|
|
||||||
if (oxygenConcentrationOutput_())
|
if (params_.ureaConcentrationOutput_) {
|
||||||
oxygenConcentration_[globalDofIdx] =
|
ureaConcentration_[globalDofIdx] =
|
||||||
scalarValue(intQuants.oxygenConcentration());
|
10 * scalarValue(intQuants.ureaConcentration());//Multypliging by scaling factor 10 (see WellInterface_impl.hpp)
|
||||||
|
}
|
||||||
|
|
||||||
if (ureaConcentrationOutput_())
|
if (params_.biofilmConcentrationOutput_) {
|
||||||
ureaConcentration_[globalDofIdx] =
|
biofilmConcentration_[globalDofIdx] =
|
||||||
10 * scalarValue(intQuants.ureaConcentration());//Multypliging by scaling factor 10 (see WellInterface_impl.hpp)
|
scalarValue(intQuants.biofilmConcentration());
|
||||||
|
}
|
||||||
if (biofilmConcentrationOutput_())
|
|
||||||
biofilmConcentration_[globalDofIdx] =
|
|
||||||
scalarValue(intQuants.biofilmConcentration());
|
|
||||||
|
|
||||||
if (calciteConcentrationOutput_())
|
|
||||||
calciteConcentration_[globalDofIdx] =
|
|
||||||
scalarValue(intQuants.calciteConcentration());
|
|
||||||
|
|
||||||
|
if (params_.calciteConcentrationOutput_) {
|
||||||
|
calciteConcentration_[globalDofIdx] =
|
||||||
|
scalarValue(intQuants.calciteConcentration());
|
||||||
|
}
|
||||||
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -175,66 +163,41 @@ public:
|
|||||||
*/
|
*/
|
||||||
void commitBuffers(BaseOutputWriter& baseWriter)
|
void commitBuffers(BaseOutputWriter& baseWriter)
|
||||||
{
|
{
|
||||||
VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
|
if constexpr (enableMICP) {
|
||||||
if (!vtkWriter)
|
VtkMultiWriter* vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
|
||||||
return;
|
if (!vtkWriter) {
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
|
||||||
if (!enableMICP)
|
if (params_.microbialConcentrationOutput_) {
|
||||||
return;
|
this->commitScalarBuffer_(baseWriter, "microbial concentration", microbialConcentration_);
|
||||||
|
}
|
||||||
|
|
||||||
if (microbialConcentrationOutput_())
|
if (params_.oxygenConcentrationOutput_) {
|
||||||
this->commitScalarBuffer_(baseWriter, "microbial concentration", microbialConcentration_);
|
this->commitScalarBuffer_(baseWriter, "oxygen concentration", oxygenConcentration_);
|
||||||
|
}
|
||||||
|
|
||||||
if (oxygenConcentrationOutput_())
|
if (params_.ureaConcentrationOutput_) {
|
||||||
this->commitScalarBuffer_(baseWriter, "oxygen concentration", oxygenConcentration_);
|
this->commitScalarBuffer_(baseWriter, "urea concentration", ureaConcentration_);
|
||||||
|
}
|
||||||
|
|
||||||
if (ureaConcentrationOutput_())
|
if (params_.biofilmConcentrationOutput_) {
|
||||||
this->commitScalarBuffer_(baseWriter, "urea concentration", ureaConcentration_);
|
this->commitScalarBuffer_(baseWriter, "biofilm fraction", biofilmConcentration_);
|
||||||
|
}
|
||||||
if (biofilmConcentrationOutput_())
|
|
||||||
this->commitScalarBuffer_(baseWriter, "biofilm fraction", biofilmConcentration_);
|
|
||||||
|
|
||||||
if (calciteConcentrationOutput_())
|
|
||||||
this->commitScalarBuffer_(baseWriter, "calcite fraction", calciteConcentration_);
|
|
||||||
|
|
||||||
|
if (params_.calciteConcentrationOutput_) {
|
||||||
|
this->commitScalarBuffer_(baseWriter, "calcite fraction", calciteConcentration_);
|
||||||
|
}
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
private:
|
private:
|
||||||
static bool microbialConcentrationOutput_()
|
VtkBlackoilMICPParams params_{};
|
||||||
{
|
ScalarBuffer microbialConcentration_{};
|
||||||
static bool val = Parameters::Get<Parameters::VtkWriteMicrobialConcentration>();
|
ScalarBuffer oxygenConcentration_{};
|
||||||
return val;
|
ScalarBuffer ureaConcentration_{};
|
||||||
}
|
ScalarBuffer biofilmConcentration_{};
|
||||||
|
ScalarBuffer calciteConcentration_{};
|
||||||
static bool oxygenConcentrationOutput_()
|
|
||||||
{
|
|
||||||
static bool val = Parameters::Get<Parameters::VtkWriteOxygenConcentration>();
|
|
||||||
return val;
|
|
||||||
}
|
|
||||||
|
|
||||||
static bool ureaConcentrationOutput_()
|
|
||||||
{
|
|
||||||
static bool val = Parameters::Get<Parameters::VtkWriteUreaConcentration>();
|
|
||||||
return val;
|
|
||||||
}
|
|
||||||
|
|
||||||
static bool biofilmConcentrationOutput_()
|
|
||||||
{
|
|
||||||
static bool val = Parameters::Get<Parameters::VtkWriteBiofilmConcentration>();
|
|
||||||
return val;
|
|
||||||
}
|
|
||||||
|
|
||||||
static bool calciteConcentrationOutput_()
|
|
||||||
{
|
|
||||||
static bool val = Parameters::Get<Parameters::VtkWriteCalciteConcentration>();
|
|
||||||
return val;
|
|
||||||
}
|
|
||||||
|
|
||||||
ScalarBuffer microbialConcentration_;
|
|
||||||
ScalarBuffer oxygenConcentration_;
|
|
||||||
ScalarBuffer ureaConcentration_;
|
|
||||||
ScalarBuffer biofilmConcentration_;
|
|
||||||
ScalarBuffer calciteConcentration_;
|
|
||||||
};
|
};
|
||||||
|
|
||||||
} // namespace Opm
|
} // namespace Opm
|
||||||
|
57
opm/models/io/vtkblackoilmicpparams.cpp
Normal file
57
opm/models/io/vtkblackoilmicpparams.cpp
Normal file
@ -0,0 +1,57 @@
|
|||||||
|
// -*- 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.
|
||||||
|
*/
|
||||||
|
|
||||||
|
#include <config.h>
|
||||||
|
#include <opm/models/io/vtkblackoilmicpparams.hpp>
|
||||||
|
|
||||||
|
#include <opm/models/utils/parametersystem.hpp>
|
||||||
|
|
||||||
|
namespace Opm {
|
||||||
|
|
||||||
|
void VtkBlackoilMICPParams::registerParameters()
|
||||||
|
{
|
||||||
|
Parameters::Register<Parameters::VtkWriteMicrobialConcentration>
|
||||||
|
("Include the concentration of the microbial component in the water phase "
|
||||||
|
"in the VTK output files");
|
||||||
|
Parameters::Register<Parameters::VtkWriteOxygenConcentration>
|
||||||
|
("Include the concentration of the oxygen component in the water phase "
|
||||||
|
"in the VTK output files");
|
||||||
|
Parameters::Register<Parameters::VtkWriteUreaConcentration>
|
||||||
|
("Include the concentration of the urea component in the water phase "
|
||||||
|
"in the VTK output files");
|
||||||
|
Parameters::Register<Parameters::VtkWriteBiofilmConcentration>
|
||||||
|
("Include the biofilm volume fraction in the VTK output files");
|
||||||
|
Parameters::Register<Parameters::VtkWriteCalciteConcentration>
|
||||||
|
("Include the calcite volume fraction in the VTK output files");
|
||||||
|
}
|
||||||
|
|
||||||
|
void VtkBlackoilMICPParams::read()
|
||||||
|
{
|
||||||
|
microbialConcentrationOutput_ = Parameters::Get<Parameters::VtkWriteMicrobialConcentration>();
|
||||||
|
oxygenConcentrationOutput_ = Parameters::Get<Parameters::VtkWriteOxygenConcentration>();
|
||||||
|
ureaConcentrationOutput_ = Parameters::Get<Parameters::VtkWriteUreaConcentration>();
|
||||||
|
biofilmConcentrationOutput_ = Parameters::Get<Parameters::VtkWriteBiofilmConcentration>();
|
||||||
|
calciteConcentrationOutput_ = Parameters::Get<Parameters::VtkWriteCalciteConcentration>();
|
||||||
|
}
|
||||||
|
|
||||||
|
} // namespace Opm
|
63
opm/models/io/vtkblackoilmicpparams.hpp
Normal file
63
opm/models/io/vtkblackoilmicpparams.hpp
Normal file
@ -0,0 +1,63 @@
|
|||||||
|
// -*- 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::VtkBlackOilMICPModule
|
||||||
|
*/
|
||||||
|
#ifndef OPM_VTK_BLACK_OIL_MICP_PARAMS_HPP
|
||||||
|
#define OPM_VTK_BLACK_OIL_MICP_PARAMS_HPP
|
||||||
|
|
||||||
|
namespace Opm::Parameters {
|
||||||
|
|
||||||
|
// set default values for what quantities to output
|
||||||
|
struct VtkWriteMicrobialConcentration { static constexpr bool value = true; };
|
||||||
|
struct VtkWriteOxygenConcentration { static constexpr bool value = true; };
|
||||||
|
struct VtkWriteUreaConcentration { static constexpr bool value = true; };
|
||||||
|
struct VtkWriteBiofilmConcentration { static constexpr bool value = true; };
|
||||||
|
struct VtkWriteCalciteConcentration { static constexpr bool value = true; };
|
||||||
|
|
||||||
|
} // namespace Opm::Parameters
|
||||||
|
|
||||||
|
namespace Opm {
|
||||||
|
|
||||||
|
/*!
|
||||||
|
* \brief Struct holding the parameters for VtkBlackoilMICPModule.
|
||||||
|
*/
|
||||||
|
struct VtkBlackoilMICPParams
|
||||||
|
{
|
||||||
|
//! \brief Registers the parameters in parameter system.
|
||||||
|
static void registerParameters();
|
||||||
|
|
||||||
|
//! \brief Reads the parameter values from the parameter system.
|
||||||
|
void read();
|
||||||
|
|
||||||
|
bool microbialConcentrationOutput_;
|
||||||
|
bool oxygenConcentrationOutput_;
|
||||||
|
bool ureaConcentrationOutput_;
|
||||||
|
bool biofilmConcentrationOutput_;
|
||||||
|
bool calciteConcentrationOutput_;
|
||||||
|
};
|
||||||
|
|
||||||
|
} // namespace Opm
|
||||||
|
|
||||||
|
#endif // OPM_VTK_BLACKOIL_MICP_PARAMS_HPP
|
Loading…
Reference in New Issue
Block a user