move ebos/eclfluxmodule.hh to opm/simulators/flow

This commit is contained in:
Arne Morten Kvarving 2024-02-02 10:46:44 +01:00
parent 2d604e12a7
commit bc2dd1110e
4 changed files with 19 additions and 20 deletions

View File

@ -412,7 +412,6 @@ list (APPEND PUBLIC_HEADER_FILES
ebos/ebos.hh
ebos/eclbasevanguard.hh
ebos/eclcpgridvanguard.hh
ebos/eclfluxmodule.hh
ebos/eclgenericcpgridvanguard.hh
ebos/eclgenericproblem.hh
ebos/eclgenericproblem_impl.hh
@ -456,6 +455,7 @@ list (APPEND PUBLIC_HEADER_FILES
opm/simulators/flow/LogOutputHelper.hpp
opm/simulators/flow/Main.hpp
opm/simulators/flow/MixingRateControls.hpp
opm/simulators/flow/NewTranFluxModule.hpp
opm/simulators/flow/NonlinearSolver.hpp
opm/simulators/flow/OutputBlackoilModule.hpp
opm/simulators/flow/partitionCells.hpp

View File

@ -35,7 +35,6 @@
#include <dune/common/fmatrix.hh>
#include <ebos/eclcpgridvanguard.hh>
#include <ebos/eclfluxmodule.hh>
#include <ebos/eclgenericproblem.hh>
#include <ebos/eclnewtonmethod.hh>
#include <ebos/eclproblem_properties.hh>
@ -80,6 +79,7 @@
#include <opm/simulators/flow/EquilInitializer.hpp>
#include <opm/simulators/flow/FIBlackoilModel.hpp>
#include <opm/simulators/flow/FlowThresholdPressure.hpp>
#include <opm/simulators/flow/NewTranFluxModule.hpp>
#include <opm/simulators/flow/OutputBlackoilModule.hpp>
#include <opm/simulators/flow/TracerModel.hpp>
#include <opm/simulators/flow/VtkTracerModule.hpp>

View File

@ -29,7 +29,6 @@
#define ECL_PROBLEM_PROPERTIES_HH
#include <ebos/eclcpgridvanguard.hh>
#include <ebos/eclfluxmodule.hh>
#include <ebos/eclnewtonmethod.hh>
#if HAVE_DAMARIS
@ -48,6 +47,7 @@
#include <opm/simulators/flow/DummyGradientCalculator.hpp>
#include <opm/simulators/flow/EclWriter.hpp>
#include <opm/simulators/flow/FIBlackoilModel.hpp>
#include <opm/simulators/flow/NewTranFluxModule.hpp>
#include <opm/simulators/flow/OutputBlackoilModule.hpp>
#include <opm/simulators/flow/VtkTracerModule.hpp>
@ -495,7 +495,7 @@ struct EnableStorageCache<TypeTag, TTag::EclBaseProblem> {
// Use the "velocity module" which uses the Eclipse "NEWTRAN" transmissibilities
template<class TypeTag>
struct FluxModule<TypeTag, TTag::EclBaseProblem> {
using type = EclTransFluxModule<TypeTag>;
using type = NewTranFluxModule<TypeTag>;
};
// Use the dummy gradient calculator in order not to do unnecessary work.

View File

@ -23,13 +23,13 @@
/*!
* \file
*
* \brief This file contains the flux module which is used for ECL problems
* \brief This file contains the flux module which is used for flow problems
*
* This approach to fluxes is very specific to two-point flux approximation and applies
* what the Eclipse Technical Description calls the "NEWTRAN" transmissibility approach.
*/
#ifndef EWOMS_ECL_FLUX_MODULE_HH
#define EWOMS_ECL_FLUX_MODULE_HH
#ifndef OPM_NEWTRAN_FLUX_MODULE_HPP
#define OPM_NEWTRAN_FLUX_MODULE_HPP
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
@ -50,24 +50,24 @@
namespace Opm {
template <class TypeTag>
class EclTransIntensiveQuantities;
class NewTranIntensiveQuantities;
template <class TypeTag>
class EclTransExtensiveQuantities;
class NewTranExtensiveQuantities;
template <class TypeTag>
class EclTransBaseProblem;
class NewTranBaseProblem;
/*!
* \ingroup EclBlackOilSimulator
* \brief Specifies a flux module which uses ECL transmissibilities.
*/
template <class TypeTag>
struct EclTransFluxModule
struct NewTranFluxModule
{
using FluxIntensiveQuantities = EclTransIntensiveQuantities<TypeTag>;
using FluxExtensiveQuantities = EclTransExtensiveQuantities<TypeTag>;
using FluxBaseProblem = EclTransBaseProblem<TypeTag>;
using FluxIntensiveQuantities = NewTranIntensiveQuantities<TypeTag>;
using FluxExtensiveQuantities = NewTranExtensiveQuantities<TypeTag>;
using FluxBaseProblem = NewTranBaseProblem<TypeTag>;
/*!
* \brief Register all run-time parameters for the flux module.
@ -82,7 +82,7 @@ struct EclTransFluxModule
* transmissibility based volume flux calculation.
*/
template <class TypeTag>
class EclTransBaseProblem
class NewTranBaseProblem
{ };
/*!
@ -90,7 +90,7 @@ class EclTransBaseProblem
* \brief Provides the intensive quantities for the ECL flux module
*/
template <class TypeTag>
class EclTransIntensiveQuantities
class NewTranIntensiveQuantities
{
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
protected:
@ -103,7 +103,7 @@ protected:
* \brief Provides the ECL flux module
*/
template <class TypeTag>
class EclTransExtensiveQuantities
class NewTranExtensiveQuantities
{
using Implementation = GetPropType<TypeTag, Properties::ExtensiveQuantities>;
@ -532,8 +532,7 @@ public:
if (upIdx[phaseIdx] == interiorDofIdx) {
// this is slightly hacky because in the automatic differentiation case, it
// only works for the element centered finite volume method. for ebos this
// does not matter, though.
// only works for the element centered finite volume method.
const auto& up = intQuantsIn;
// deal with water induced rock compaction
@ -589,4 +588,4 @@ private:
} // namespace Opm
#endif
#endif // OPM_NEWTRAN_FLUX_MODULE_HPP