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 <boost/unordered_map.hpp>
40
41#include <opm/core/utility/SparseVector.hpp>
42#include <opm/core/utility/parameters/ParameterGroup.hpp>
43
46
47#include <dune/grid/common/GridAdapter.hpp>
48
49#include <opm/core/transport/implicit/ImplicitAssembly.hpp>
50#include <opm/core/transport/implicit/ImplicitTransport.hpp>
51#include <opm/core/transport/implicit/JacobianSystem.hpp>
52#include <opm/core/transport/implicit/SinglePointUpwindTwoPhase.hpp>
53
54
55namespace Opm {
56
59 template <class GridInterface, class ReservoirProperties, class BoundaryConditions>
61 {
62 public:
69 EulerUpstreamImplicit (const GridInterface& grid,
70 const ReservoirProperties& resprop,
71 const BoundaryConditions& boundary);
75 void init(const Opm::parameter::ParameterGroup& param);
79 void init(const Opm::parameter::ParameterGroup& param,
80 const GridInterface& grid,
81 const ReservoirProperties& resprop,
82 const BoundaryConditions& boundary);
86 void initObj(const GridInterface& grid,
87 const ReservoirProperties& resprop,
88 const BoundaryConditions& boundary);
92 void display();
93
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
123 typedef Opm::SinglePointUpwindTwoPhase<TwophaseFluid> TransportModel;
124# if 0
125 typedef Opm::ImplicitTransportDefault::NewtonVectorCollection< ::std::vector<double> > NVecColl;
126 typedef Opm::ImplicitTransportDefault::JacobianSystem < struct CSRMatrix, NVecColl > JacSys;
127 typedef Opm::ImplicitTransportLinAlgSupport::CSRMatrixUmfpackSolver LinearSolver;
128# else
129 typedef Dune::FieldVector<double, 1> ScalarVectorBlockType;
130 typedef Dune::FieldMatrix<double, 1, 1> ScalarMatrixBlockType;
131
132 typedef Dune::BlockVector<ScalarVectorBlockType> ScalarBlockVector;
133 typedef Dune::BCRSMatrix <ScalarMatrixBlockType> ScalarBCRSMatrix;
134
135 typedef Opm::ImplicitTransportDefault::NewtonVectorCollection< ScalarBlockVector > NVecColl;
136 typedef Opm::ImplicitTransportDefault::JacobianSystem < ScalarBCRSMatrix, NVecColl > JacSys;
137
139# endif
140 typedef Opm::ImplicitTransport<TransportModel ,
141 JacSys ,
143 Opm::ImplicitTransportDefault::VectorNegater,
144 Opm::ImplicitTransportDefault::VectorZero ,
145 Opm::ImplicitTransportDefault::MatrixZero ,
146 Opm::ImplicitTransportDefault::VectorAssign > TransportSolver;
147
148 // should be initialized by param
149
150 //TransportModel model_;
151 //TransportSolver tsolver_;
152 GridAdapter mygrid_;
153 ReservoirProperties myrp_;
154
155 //Opm::SimpleFluid2pWrapper< ReservoirProperties > myfluid_;
156 //TwophaseFluid myfluid_;
157
161 std::vector<double> porevol_;
162
163 std::vector<int> periodic_cells_;
164 std::vector<int> periodic_faces_;
165 std::vector<int> periodic_nbfaces_;
166 std::vector<int> periodic_hfaces_;
167 std::vector<int> direclet_cells_;
168 std::vector<double> direclet_sat_;
169 std::vector<double> direclet_hfaces_;
170 //std::vector< double > trans_;
171 std::vector<double> htrans_;
172 Opm::ImplicitTransportDetails::NRControl ctrl_;
173 // Storing residual so that we won't have to reallocate it for every step.
174 //mutable std::vector<double> residual_;
175 };
176
177} // namespace Opm
178
180
181#endif // OPENRS_EULERUPSTREAM_HEADER
Definition: EulerUpstreamImplicit.hpp:61
Dune::FieldMatrix< double, 1, 1 > ScalarMatrixBlockType
Definition: EulerUpstreamImplicit.hpp:130
bool clamp_sat_
Definition: EulerUpstreamImplicit.hpp:159
std::vector< double > direclet_sat_
Definition: EulerUpstreamImplicit.hpp:168
Opm::LinearSolverBICGSTAB LinearSolver
Definition: EulerUpstreamImplicit.hpp:138
std::vector< int > periodic_faces_
Definition: EulerUpstreamImplicit.hpp:164
EulerUpstreamImplicit()
Definition: EulerUpstreamImplicit_impl.hpp:60
Opm::TwophaseFluidWrapper TwophaseFluid
Definition: EulerUpstreamImplicit.hpp:122
std::vector< double > direclet_hfaces_
Definition: EulerUpstreamImplicit.hpp:169
std::vector< int > direclet_cells_
Definition: EulerUpstreamImplicit.hpp:167
std::vector< double > htrans_
Definition: EulerUpstreamImplicit.hpp:171
void checkAndPossiblyClampSat(std::vector< double > &s) const
Definition: EulerUpstreamImplicit_impl.hpp:342
Dune::BCRSMatrix< ScalarMatrixBlockType > ScalarBCRSMatrix
Definition: EulerUpstreamImplicit.hpp:133
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:132
std::vector< int > periodic_cells_
Definition: EulerUpstreamImplicit.hpp:163
void init(const Opm::parameter::ParameterGroup &param, const GridInterface &grid, const ReservoirProperties &resprop, const BoundaryConditions &boundary)
Opm::SinglePointUpwindTwoPhase< TwophaseFluid > TransportModel
Definition: EulerUpstreamImplicit.hpp:123
int max_repeats_
Definition: EulerUpstreamImplicit.hpp:160
Opm::ImplicitTransport< TransportModel, JacSys, Opm::MaxNormDune, Opm::ImplicitTransportDefault::VectorNegater, Opm::ImplicitTransportDefault::VectorZero, Opm::ImplicitTransportDefault::MatrixZero, Opm::ImplicitTransportDefault::VectorAssign > TransportSolver
Definition: EulerUpstreamImplicit.hpp:146
void initObj(const GridInterface &grid, const ReservoirProperties &resprop, const BoundaryConditions &boundary)
Definition: EulerUpstreamImplicit_impl.hpp:98
Dune::FieldVector< double, 1 > ScalarVectorBlockType
Definition: EulerUpstreamImplicit.hpp:129
Opm::ImplicitTransportDefault::NewtonVectorCollection< ScalarBlockVector > NVecColl
Definition: EulerUpstreamImplicit.hpp:135
bool check_sat_
Definition: EulerUpstreamImplicit.hpp:158
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:165
void display()
Definition: EulerUpstreamImplicit_impl.hpp:219
std::vector< double > porevol_
Definition: EulerUpstreamImplicit.hpp:161
void init(const Opm::parameter::ParameterGroup &param)
Definition: EulerUpstreamImplicit_impl.hpp:74
Opm::ImplicitTransportDetails::NRControl ctrl_
Definition: EulerUpstreamImplicit.hpp:172
Opm::ImplicitTransportDefault::JacobianSystem< ScalarBCRSMatrix, NVecColl > JacSys
Definition: EulerUpstreamImplicit.hpp:136
std::vector< int > periodic_hfaces_
Definition: EulerUpstreamImplicit.hpp:166
GridAdapter mygrid_
Definition: EulerUpstreamImplicit.hpp:152
EulerUpstreamImplicit(const GridInterface &grid, const ReservoirProperties &resprop, const BoundaryConditions &boundary)
ReservoirProperties myrp_
Definition: EulerUpstreamImplicit.hpp:153
Definition: ImplicitTransportDefs.hpp:174
Definition: ImplicitTransportDefs.hpp:48
Definition: ImplicitTransportDefs.hpp:120
Definition: BlackoilFluid.hpp:32