EulerUpstreamImplicit.hpp
Go to the documentation of this file.
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
44
45#include <opm/grid/common/GridAdapter.hpp>
46
51
52
53namespace Opm {
54
57 template <class GridInterface, class ReservoirProperties, class BoundaryConditions>
59 {
60 public:
67 EulerUpstreamImplicit (const GridInterface& grid,
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
96 template <class PressureSolution>
97 bool transportSolve(std::vector<double>& saturation,
98 const double time,
99 const typename GridInterface::Vector& gravity,
100 const PressureSolution& pressure_sol,
101 const Opm::SparseVector<double>& injection_rates) const;
102
103 protected:
104 // typedef typename GridInterface::CellIterator CIt;
105 // typedef typename CIt::FaceIterator FIt;
106 // typedef typename FIt::Vector Vector;
107
108 template <class PressureSolution>
109 void smallTimeStep(std::vector<double>& saturation,
110 const double time,
111 const typename GridInterface::Vector& gravity,
112 const PressureSolution& pressure_sol,
113 const Opm::SparseVector<double>& injection_rates) const;
114
115 void checkAndPossiblyClampSat(std::vector<double>& s) const;
116
117
118
119 // Interface to implicit solver
122 typedef Dune::FieldVector<double, 1> ScalarVectorBlockType;
123 typedef Dune::FieldMatrix<double, 1, 1> ScalarMatrixBlockType;
124
125 typedef Dune::BlockVector<ScalarVectorBlockType> ScalarBlockVector;
126 typedef Dune::BCRSMatrix <ScalarMatrixBlockType> ScalarBCRSMatrix;
127
130
133 JacSys ,
139
140 // should be initialized by param
141 GridAdapter mygrid_;
142 ReservoirProperties myrp_;
143
147 std::vector<double> porevol_;
148
149 std::vector<int> periodic_cells_;
150 std::vector<int> periodic_faces_;
151 std::vector<int> periodic_nbfaces_;
152 std::vector<int> periodic_hfaces_;
153 std::vector<int> direclet_cells_;
154 std::vector<double> direclet_sat_;
155 std::vector<double> direclet_hfaces_;
156 std::vector<double> htrans_;
158 };
159
160} // namespace Opm
161
163
164#endif // OPENRS_EULERUPSTREAM_HEADER
Definition: EulerUpstreamImplicit.hpp:59
Dune::FieldMatrix< double, 1, 1 > ScalarMatrixBlockType
Definition: EulerUpstreamImplicit.hpp:123
bool clamp_sat_
Definition: EulerUpstreamImplicit.hpp:145
std::vector< double > direclet_sat_
Definition: EulerUpstreamImplicit.hpp:154
Opm::LinearSolverBICGSTAB LinearSolver
Definition: EulerUpstreamImplicit.hpp:131
std::vector< int > periodic_faces_
Definition: EulerUpstreamImplicit.hpp:150
EulerUpstreamImplicit()
Definition: EulerUpstreamImplicit_impl.hpp:60
Opm::TwophaseFluidWrapper TwophaseFluid
Definition: EulerUpstreamImplicit.hpp:120
std::vector< double > direclet_hfaces_
Definition: EulerUpstreamImplicit.hpp:155
std::vector< int > direclet_cells_
Definition: EulerUpstreamImplicit.hpp:153
std::vector< double > htrans_
Definition: EulerUpstreamImplicit.hpp:156
void checkAndPossiblyClampSat(std::vector< double > &s) const
Definition: EulerUpstreamImplicit_impl.hpp:340
Dune::BCRSMatrix< ScalarMatrixBlockType > ScalarBCRSMatrix
Definition: EulerUpstreamImplicit.hpp:126
void smallTimeStep(std::vector< double > &saturation, const double time, const typename GridInterface::Vector &gravity, const PressureSolution &pressure_sol, const Opm::SparseVector< double > &injection_rates) const
Dune::BlockVector< ScalarVectorBlockType > ScalarBlockVector
Definition: EulerUpstreamImplicit.hpp:125
void init(const Opm::ParameterGroup &param, const GridInterface &grid, const ReservoirProperties &resprop, const BoundaryConditions &boundary)
std::vector< int > periodic_cells_
Definition: EulerUpstreamImplicit.hpp:149
Opm::SinglePointUpwindTwoPhase< TwophaseFluid > TransportModel
Definition: EulerUpstreamImplicit.hpp:121
int max_repeats_
Definition: EulerUpstreamImplicit.hpp:146
Opm::ImplicitTransport< TransportModel, JacSys, Opm::MaxNormDune, Opm::ImplicitTransportDefault::VectorNegater, Opm::ImplicitTransportDefault::VectorZero, Opm::ImplicitTransportDefault::MatrixZero, Opm::ImplicitTransportDefault::VectorAssign > TransportSolver
Definition: EulerUpstreamImplicit.hpp:138
void initObj(const GridInterface &grid, const ReservoirProperties &resprop, const BoundaryConditions &boundary)
Definition: EulerUpstreamImplicit_impl.hpp:98
Dune::FieldVector< double, 1 > ScalarVectorBlockType
Definition: EulerUpstreamImplicit.hpp:122
Opm::ImplicitTransportDefault::NewtonVectorCollection< ScalarBlockVector > NVecColl
Definition: EulerUpstreamImplicit.hpp:128
bool check_sat_
Definition: EulerUpstreamImplicit.hpp:144
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, evolving.
std::vector< int > periodic_nbfaces_
Definition: EulerUpstreamImplicit.hpp:151
void display()
Definition: EulerUpstreamImplicit_impl.hpp:217
std::vector< double > porevol_
Definition: EulerUpstreamImplicit.hpp:147
void init(const Opm::ParameterGroup &param)
Definition: EulerUpstreamImplicit_impl.hpp:74
Opm::ImplicitTransportDetails::NRControl ctrl_
Definition: EulerUpstreamImplicit.hpp:157
Opm::ImplicitTransportDefault::JacobianSystem< ScalarBCRSMatrix, NVecColl > JacSys
Definition: EulerUpstreamImplicit.hpp:129
std::vector< int > periodic_hfaces_
Definition: EulerUpstreamImplicit.hpp:152
GridAdapter mygrid_
Definition: EulerUpstreamImplicit.hpp:141
EulerUpstreamImplicit(const GridInterface &grid, const ReservoirProperties &resprop, const BoundaryConditions &boundary)
ReservoirProperties myrp_
Definition: EulerUpstreamImplicit.hpp:142
Definition: JacobianSystem.hpp:222
Definition: JacobianSystem.hpp:108
Definition: JacobianSystem.hpp:89
Definition: JacobianSystem.hpp:65
Definition: JacobianSystem.hpp:78
Definition: ImplicitTransport.hpp:105
Definition: ImplicitTransportDefs.hpp:162
Definition: ImplicitTransportDefs.hpp:43
Definition: SinglePointUpwindTwoPhase.hpp:278
Definition: ImplicitTransportDefs.hpp:115
Dune::BlockVector< Dune::FieldVector< double, 1 > > Vector
A vector holding our RHS.
Definition: matrixops.hpp:33
Definition: ImplicitAssembly.hpp:43
Definition: ImplicitTransport.hpp:45