opm-simulators/opm/autodiff/SimulatorFullyImplicitBlackoilEbos.hpp
Robert Kloefkorn 3db63b0a22 add flow_ebos, an ebos based simulator
it uses ebos for linearization of the mass balance equations and the
current flow code from opm-simulators for all the rest. currently, the
results match the ones from plain `flow` for SPE1, SPE9 and Norne, but
performance is not optimal: on SPE9, converting from and to the legacy
data structures takes about a third of the time to do the actual mass
balance assembly. nevertheless `flow_ebos` is almost as fast as plain
`flow` for SPE9. (for Norne `flow_ebos` is about 15% slower, even
though the results match quite closely. the reason for this is that it
requires more iterations for some reason.)
2016-08-09 18:38:23 +02:00

74 lines
2.8 KiB
C++

/*
Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
Copyright 2015 Andreas Lauser
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_SIMULATORFULLYIMPLICITBLACKOIL_HEADER_INCLUDED
#define OPM_SIMULATORFULLYIMPLICITBLACKOIL_HEADER_INCLUDED
#include <opm/autodiff/SimulatorBase.hpp>
#include <opm/autodiff/NonlinearSolver.hpp>
#include <opm/autodiff/BlackoilModelEbos.hpp>
namespace Opm {
template <class GridT>
class SimulatorFullyImplicitBlackoilEbos;
class StandardWells;
template <class GridT>
struct SimulatorTraits<SimulatorFullyImplicitBlackoilEbos<GridT> >
{
typedef WellStateFullyImplicitBlackoil WellState;
typedef BlackoilState ReservoirState;
typedef BlackoilOutputWriter OutputWriter;
typedef GridT Grid;
typedef BlackoilModelEbos<Grid> Model;
typedef NonlinearSolver<Model> Solver;
typedef StandardWells WellModel;
};
/// a simulator for the blackoil model
template <class GridT>
class SimulatorFullyImplicitBlackoilEbos
: public SimulatorBase<SimulatorFullyImplicitBlackoilEbos<GridT> >
{
typedef SimulatorBase<SimulatorFullyImplicitBlackoilEbos<GridT> > Base;
public:
// forward the constructor to the base class
SimulatorFullyImplicitBlackoilEbos(const parameter::ParameterGroup& param,
const typename Base::Grid& grid,
DerivedGeology& geo,
BlackoilPropsAdInterface& props,
const RockCompressibility* rock_comp_props,
NewtonIterationBlackoilInterface& linsolver,
const double* gravity,
const bool disgas,
const bool vapoil,
std::shared_ptr<EclipseState> eclipse_state,
BlackoilOutputWriter& output_writer,
const std::vector<double>& threshold_pressures_by_face)
: Base(param, grid, geo, props, rock_comp_props, linsolver, gravity, disgas, vapoil,
eclipse_state, output_writer, threshold_pressures_by_face)
{}
};
} // namespace Opm
#endif // OPM_SIMULATORFULLYIMPLICITBLACKOIL_HEADER_INCLUDED