opm-simulators/opm/simulators/aquifers/BlackoilAquiferModel.hpp

118 lines
3.7 KiB
C++
Raw Normal View History

/*
File adapted from BlackoilWellModel.hpp
Copyright 2017 TNO - Heat Transfer & Fluid Dynamics, Modelling & Optimization of the Subsurface
Copyright 2017 Statoil ASA.
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 3 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/>.
*/
#ifndef OPM_BLACKOILAQUIFERMODEL_HEADER_INCLUDED
#define OPM_BLACKOILAQUIFERMODEL_HEADER_INCLUDED
#include <ebos/eclbaseaquifermodel.hh>
#include <opm/parser/eclipse/EclipseState/Aquancon.hpp>
#include <opm/parser/eclipse/EclipseState/AquiferCT.hpp>
#include <opm/parser/eclipse/EclipseState/Aquifetp.hpp>
#include <opm/output/data/Aquifer.hpp>
#include <opm/simulators/aquifers/AquiferCarterTracy.hpp>
#include <opm/simulators/aquifers/AquiferFetkovich.hpp>
2020-12-11 08:33:25 -06:00
#include <opm/grid/CpGrid.hpp>
#include <opm/grid/polyhedralgrid.hh>
#if HAVE_DUNE_ALUGRID
#include <dune/alugrid/grid.hh>
#endif
#include <opm/material/densead/Math.hpp>
#include <vector>
namespace Opm
{
/// Class for handling the blackoil well model.
template <typename TypeTag>
class BlackoilAquiferModel
{
using Simulator = GetPropType<TypeTag, Properties::Simulator>;
using RateVector = GetPropType<TypeTag, Properties::RateVector>;
2020-12-11 08:33:25 -06:00
constexpr bool supportsFaceTag(const Dune::CpGrid&){ return true;}
constexpr bool supportsFaceTag(const Dune::PolyhedralGrid<3, 3>&){ return true;}
#if HAVE_DUNE_ALUGRID
constexpr bool supportsFaceTag(const Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming>&){ return true;}
#endif
template<class G>
constexpr bool supportsFaceTag(const G&){ return false;}
public:
explicit BlackoilAquiferModel(Simulator& simulator);
void initialSolutionApplied();
void initFromRestart(const std::vector<data::AquiferData>& aquiferSoln);
2018-10-30 09:50:36 -05:00
void beginEpisode();
void beginTimeStep();
void beginIteration();
// add the water rate due to aquifers to the source term.
template <class Context>
void addToSource(RateVector& rates, const Context& context, unsigned spaceIdx, unsigned timeIdx) const;
void endIteration();
void endTimeStep();
void endEpisode();
Opm::data::Aquifers aquiferData() const;
template <class Restarter>
void serialize(Restarter& res);
template <class Restarter>
void deserialize(Restarter& res);
protected:
// --------- Types ---------
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
typedef AquiferCarterTracy<TypeTag> AquiferCarterTracy_object;
typedef AquiferFetkovich<TypeTag> AquiferFetkovich_object;
2018-10-30 09:50:36 -05:00
Simulator& simulator_;
// TODO: probably we can use one variable to store both types of aquifers, because
// they share the base class
mutable std::vector<AquiferCarterTracy_object> aquifers_CarterTracy;
mutable std::vector<AquiferFetkovich_object> aquifers_Fetkovich;
// This initialization function is used to connect the parser objects with the ones needed by AquiferCarterTracy
void init();
bool aquiferActive() const;
bool aquiferCarterTracyActive() const;
bool aquiferFetkovichActive() const;
};
} // namespace Opm
#include "BlackoilAquiferModel_impl.hpp"
2018-06-05 02:16:34 -05:00
#endif