mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
vtkprimaryvarsparams: move parameters to dedicated struct with a translation unit
This commit is contained in:
parent
be4026ecce
commit
023328bed2
@ -73,6 +73,7 @@ list (APPEND MAIN_SOURCE_FILES
|
||||
opm/models/io/vtkenergyparams.cpp
|
||||
opm/models/io/vtkmultiphaseparams.cpp
|
||||
opm/models/io/vtkphasepresenceparams.cpp
|
||||
opm/models/io/vtkprimaryvarsparams.cpp
|
||||
opm/models/io/restart.cpp
|
||||
opm/models/parallel/mpiutil.cpp
|
||||
opm/models/parallel/tasklets.cpp
|
||||
@ -684,6 +685,7 @@ list (APPEND PUBLIC_HEADER_FILES
|
||||
opm/models/io/vtkphasepresencemodule.hpp
|
||||
opm/models/io/vtkphasepresenceparams.hpp
|
||||
opm/models/io/vtkprimaryvarsmodule.hpp
|
||||
opm/models/io/vtkprimaryvarsparams.hpp
|
||||
opm/models/io/vtkptflashmodule.hh
|
||||
opm/models/io/vtkscalarfunction.hh
|
||||
opm/models/io/vtktemperaturemodule.hh
|
||||
|
@ -31,18 +31,11 @@
|
||||
|
||||
#include <opm/models/io/baseoutputmodule.hh>
|
||||
#include <opm/models/io/vtkmultiwriter.hh>
|
||||
#include <opm/models/io/vtkprimaryvarsparams.hpp>
|
||||
|
||||
#include <opm/models/utils/parametersystem.hpp>
|
||||
#include <opm/models/utils/propertysystem.hh>
|
||||
|
||||
namespace Opm::Parameters {
|
||||
|
||||
struct VtkWritePrimaryVars { static constexpr bool value = false; };
|
||||
struct VtkWriteProcessRank { static constexpr bool value = false; };
|
||||
struct VtkWriteDofIndex { static constexpr bool value = false; };
|
||||
|
||||
} // namespace Opm::Properties
|
||||
|
||||
namespace Opm {
|
||||
|
||||
/*!
|
||||
@ -70,19 +63,16 @@ class VtkPrimaryVarsModule : public BaseOutputModule<TypeTag>
|
||||
public:
|
||||
VtkPrimaryVarsModule(const Simulator& simulator)
|
||||
: ParentType(simulator)
|
||||
{ }
|
||||
{
|
||||
params_.read();
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief Register all run-time parameters for the Vtk output module.
|
||||
*/
|
||||
static void registerParameters()
|
||||
{
|
||||
Parameters::Register<Parameters::VtkWritePrimaryVars>
|
||||
("Include the primary variables into the VTK output files");
|
||||
Parameters::Register<Parameters::VtkWriteProcessRank>
|
||||
("Include the MPI process rank into the VTK output files");
|
||||
Parameters::Register<Parameters::VtkWriteDofIndex>
|
||||
("Include the index of the degrees of freedom into the VTK output files");
|
||||
VtkPrimaryVarsParams::registerParameters();
|
||||
}
|
||||
|
||||
/*!
|
||||
@ -91,13 +81,16 @@ public:
|
||||
*/
|
||||
void allocBuffers()
|
||||
{
|
||||
if (primaryVarsOutput_())
|
||||
if (params_.primaryVarsOutput_) {
|
||||
this->resizeEqBuffer_(primaryVars_);
|
||||
if (processRankOutput_())
|
||||
}
|
||||
if (params_.processRankOutput_) {
|
||||
this->resizeScalarBuffer_(processRank_,
|
||||
/*bufferType=*/ParentType::ElementBuffer);
|
||||
if (dofIndexOutput_())
|
||||
}
|
||||
if (params_.dofIndexOutput_) {
|
||||
this->resizeScalarBuffer_(dofIndex_);
|
||||
}
|
||||
}
|
||||
|
||||
/*!
|
||||
@ -112,19 +105,22 @@ public:
|
||||
|
||||
const auto& elementMapper = elemCtx.model().elementMapper();
|
||||
unsigned elemIdx = static_cast<unsigned>(elementMapper.index(elemCtx.element()));
|
||||
if (processRankOutput_() && !processRank_.empty())
|
||||
if (params_.processRankOutput_ && !processRank_.empty()) {
|
||||
processRank_[elemIdx] = static_cast<unsigned>(this->simulator_.gridView().comm().rank());
|
||||
}
|
||||
|
||||
for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
|
||||
unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
|
||||
const auto& priVars = elemCtx.primaryVars(i, /*timeIdx=*/0);
|
||||
|
||||
if (dofIndexOutput_())
|
||||
if (params_.dofIndexOutput_) {
|
||||
dofIndex_[I] = I;
|
||||
}
|
||||
|
||||
for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
|
||||
if (primaryVarsOutput_() && !primaryVars_[eqIdx].empty())
|
||||
if (params_.primaryVarsOutput_ && !primaryVars_[eqIdx].empty()) {
|
||||
primaryVars_[eqIdx][I] = priVars[eqIdx];
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
@ -134,42 +130,30 @@ public:
|
||||
*/
|
||||
void commitBuffers(BaseOutputWriter& baseWriter)
|
||||
{
|
||||
VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
|
||||
VtkMultiWriter* vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
|
||||
if (!vtkWriter) {
|
||||
return;
|
||||
}
|
||||
|
||||
if (primaryVarsOutput_())
|
||||
if (params_.primaryVarsOutput_) {
|
||||
this->commitPriVarsBuffer_(baseWriter, "PV_%s", primaryVars_);
|
||||
if (processRankOutput_())
|
||||
}
|
||||
if (params_.processRankOutput_) {
|
||||
this->commitScalarBuffer_(baseWriter,
|
||||
"process rank",
|
||||
processRank_,
|
||||
/*bufferType=*/ParentType::ElementBuffer);
|
||||
if (dofIndexOutput_())
|
||||
}
|
||||
if (params_.dofIndexOutput_) {
|
||||
this->commitScalarBuffer_(baseWriter, "DOF index", dofIndex_);
|
||||
}
|
||||
}
|
||||
|
||||
private:
|
||||
static bool primaryVarsOutput_()
|
||||
{
|
||||
static bool val = Parameters::Get<Parameters::VtkWritePrimaryVars>();
|
||||
return val;
|
||||
}
|
||||
static bool processRankOutput_()
|
||||
{
|
||||
static bool val = Parameters::Get<Parameters::VtkWriteProcessRank>();
|
||||
return val;
|
||||
}
|
||||
static bool dofIndexOutput_()
|
||||
{
|
||||
static bool val = Parameters::Get<Parameters::VtkWriteDofIndex>();
|
||||
return val;
|
||||
}
|
||||
|
||||
EqBuffer primaryVars_;
|
||||
ScalarBuffer processRank_;
|
||||
ScalarBuffer dofIndex_;
|
||||
VtkPrimaryVarsParams params_{};
|
||||
EqBuffer primaryVars_{};
|
||||
ScalarBuffer processRank_{};
|
||||
ScalarBuffer dofIndex_{};
|
||||
};
|
||||
|
||||
} // namespace Opm
|
||||
|
48
opm/models/io/vtkprimaryvarsparams.cpp
Normal file
48
opm/models/io/vtkprimaryvarsparams.cpp
Normal file
@ -0,0 +1,48 @@
|
||||
// -*- 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/vtkprimaryvarsparams.hpp>
|
||||
|
||||
#include <opm/models/utils/parametersystem.hpp>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
void VtkPrimaryVarsParams::registerParameters()
|
||||
{
|
||||
Parameters::Register<Parameters::VtkWritePrimaryVars>
|
||||
("Include the primary variables into the VTK output files");
|
||||
Parameters::Register<Parameters::VtkWriteProcessRank>
|
||||
("Include the MPI process rank into the VTK output files");
|
||||
Parameters::Register<Parameters::VtkWriteDofIndex>
|
||||
("Include the index of the degrees of freedom into the VTK output files");
|
||||
}
|
||||
|
||||
void VtkPrimaryVarsParams::read()
|
||||
{
|
||||
primaryVarsOutput_ = Parameters::Get<Parameters::VtkWritePrimaryVars>();
|
||||
processRankOutput_ = Parameters::Get<Parameters::VtkWriteProcessRank>();
|
||||
dofIndexOutput_ = Parameters::Get<Parameters::VtkWriteDofIndex>();
|
||||
}
|
||||
|
||||
} // namespace Opm
|
58
opm/models/io/vtkprimaryvarsparams.hpp
Normal file
58
opm/models/io/vtkprimaryvarsparams.hpp
Normal file
@ -0,0 +1,58 @@
|
||||
// -*- 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::VtkPrimaryVarsModule
|
||||
*/
|
||||
#ifndef OPM_VTK_PRIMARY_VARS_PARAMS_HPP
|
||||
#define OPM_VTK_PRIMARY_VARS_PARAMS_HPP
|
||||
|
||||
namespace Opm::Parameters {
|
||||
|
||||
struct VtkWritePrimaryVars { static constexpr bool value = false; };
|
||||
struct VtkWriteProcessRank { static constexpr bool value = false; };
|
||||
struct VtkWriteDofIndex { static constexpr bool value = false; };
|
||||
|
||||
} // namespace Opm::Parameters
|
||||
|
||||
namespace Opm {
|
||||
|
||||
/*!
|
||||
* \brief Struct holding the parameters for VtkPrimaryPhaseModule.
|
||||
*/
|
||||
struct VtkPrimaryVarsParams
|
||||
{
|
||||
//! \brief Registers the parameters in parameter system.
|
||||
static void registerParameters();
|
||||
|
||||
//! \brief Reads the parameter values from the parameter system.
|
||||
void read();
|
||||
|
||||
bool primaryVarsOutput_;
|
||||
bool processRankOutput_;
|
||||
bool dofIndexOutput_;
|
||||
};
|
||||
|
||||
} // namespace Opm
|
||||
|
||||
#endif // OPM_VTK_PRIMARY_VARS_PARAMS_HPP
|
Loading…
Reference in New Issue
Block a user