opm-upscaling
Loading...
Searching...
No Matches
EulerUpstreamImplicit.hpp
1//===========================================================================
2//
3// File: EulerUpstreamImplicit.hpp
4//
5// Created: Tue Jun 16 14:07:53 2009
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8// B�rd Skaflestad <bard.skaflestad@sintef.no>
9//
10// $Date$
11//
12// $Revision$
13//
14//===========================================================================
15
16/*
17Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
18Copyright 2009, 2010 Statoil ASA.
19
20This file is part of The Open Reservoir Simulator Project (OpenRS).
21
22OpenRS is free software: you can redistribute it and/or modify
23it under the terms of the GNU General Public License as published by
24the Free Software Foundation, either version 3 of the License, or
25(at your option) any later version.
26
27OpenRS is distributed in the hope that it will be useful,
28but WITHOUT ANY WARRANTY; without even the implied warranty of
29MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30GNU General Public License for more details.
31
32You should have received a copy of the GNU General Public License
33along with OpenRS. If not, see <http://www.gnu.org/licenses/>.
34*/
35
36#ifndef OPENRS_EULERUPSTREAMIMPLICIT_HEADER
37#define OPENRS_EULERUPSTREAMIMPLICIT_HEADER
38
39#include <opm/common/utility/numeric/SparseVector.hpp>
40#include <opm/common/utility/parameters/ParameterGroup.hpp>
41
42#include <opm/porsol/common/ImplicitTransportDefs.hpp>
43#include <opm/porsol/common/BCRSMatrixBlockAssembler.hpp>
44
45#include <opm/grid/common/GridAdapter.hpp>
46
47#include <opm/core/transport/implicit/ImplicitAssembly.hpp>
48#include <opm/core/transport/implicit/ImplicitTransport.hpp>
49#include <opm/core/transport/implicit/JacobianSystem.hpp>
51
52
53namespace Opm {
54
57 template <class GridInterface, class ReservoirProperties, class BoundaryConditions>
59 {
60 public:
68 const ReservoirProperties& resprop,
69 const BoundaryConditions& boundary);
73 void init(const Opm::ParameterGroup& param);
77 void init(const Opm::ParameterGroup& param,
78 const GridInterface& grid,
79 const ReservoirProperties& resprop,
80 const BoundaryConditions& boundary);
84 void initObj(const GridInterface& grid,
85 const ReservoirProperties& resprop,
86 const BoundaryConditions& boundary);
90 void display();
91
98 template <class PressureSolution>
99 bool transportSolve(std::vector<double>& saturation,
100 const double time,
101 const typename GridInterface::Vector& gravity,
102 const PressureSolution& pressure_sol,
103 const Opm::SparseVector<double>& injection_rates) const;
104
105 protected:
106 // typedef typename GridInterface::CellIterator CIt;
107 // typedef typename CIt::FaceIterator FIt;
108 // typedef typename FIt::Vector Vector;
109
110 template <class PressureSolution>
111 void smallTimeStep(std::vector<double>& saturation,
112 const double time,
113 const typename GridInterface::Vector& gravity,
114 const PressureSolution& pressure_sol,
115 const Opm::SparseVector<double>& injection_rates) const;
116
117 void checkAndPossiblyClampSat(std::vector<double>& s) const;
118
119
120
121 // Interface to implicit solver
122 typedef Opm::TwophaseFluidWrapper TwophaseFluid;
124 typedef Dune::FieldVector<double, 1> ScalarVectorBlockType;
125 typedef Dune::FieldMatrix<double, 1, 1> ScalarMatrixBlockType;
126
127 typedef Dune::BlockVector<ScalarVectorBlockType> ScalarBlockVector;
128 typedef Dune::BCRSMatrix <ScalarMatrixBlockType> ScalarBCRSMatrix;
129
131 typedef Opm::ImplicitTransportDefault::JacobianSystem < ScalarBCRSMatrix, NVecColl > JacSys;
132
133 typedef Opm::LinearSolverBICGSTAB LinearSolver;
134 typedef Opm::ImplicitTransport<TransportModel ,
135 JacSys ,
141
142 // should be initialized by param
143 GridAdapter mygrid_;
144 ReservoirProperties myrp_;
145
146 bool check_sat_;
147 bool clamp_sat_;
148 int max_repeats_;
149 std::vector<double> porevol_;
150
151 std::vector<int> periodic_cells_;
152 std::vector<int> periodic_faces_;
153 std::vector<int> periodic_nbfaces_;
154 std::vector<int> periodic_hfaces_;
155 std::vector<int> direclet_cells_;
156 std::vector<double> direclet_sat_;
157 std::vector<double> direclet_hfaces_;
158 std::vector<double> htrans_;
160 };
161
162} // namespace Opm
163
164#include "EulerUpstreamImplicit_impl.hpp"
165
166#endif // OPENRS_EULERUPSTREAM_HEADER
Numerical model and support classes needed to model transport of two incompressible fluid phases.
EulerUpstreamImplicit()
Definition EulerUpstreamImplicit_impl.hpp:60
void init(const Opm::ParameterGroup &param, const GridInterface &grid, const ReservoirProperties &resprop, const BoundaryConditions &boundary)
void initObj(const GridInterface &grid, const ReservoirProperties &resprop, const BoundaryConditions &boundary)
Definition EulerUpstreamImplicit_impl.hpp:98
bool transportSolve(std::vector< double > &saturation, const double time, const typename GridInterface::Vector &gravity, const PressureSolution &pressure_sol, const Opm::SparseVector< double > &injection_rates) const
Solve transport equation.
void display()
Definition EulerUpstreamImplicit_impl.hpp:217
void init(const Opm::ParameterGroup &param)
Definition EulerUpstreamImplicit_impl.hpp:74
EulerUpstreamImplicit(const GridInterface &grid, const ReservoirProperties &resprop, const BoundaryConditions &boundary)
Definition JacobianSystem.hpp:102
Definition JacobianSystem.hpp:84
Definition JacobianSystem.hpp:62
Definition JacobianSystem.hpp:73
Definition ImplicitTransport.hpp:105
Definition ImplicitTransportDefs.hpp:162
Definition ImplicitTransportDefs.hpp:43
Definition SinglePointUpwindTwoPhase.hpp:278
Definition ImplicitTransportDefs.hpp:115
Inverting small matrices.
Definition ImplicitAssembly.hpp:43
Definition ImplicitTransport.hpp:45