MultisegmentWellEquations.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 Statoil ASA.
4 Copyright 2016 - 2017 IRIS AS.
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 3 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20*/
21
22#ifndef OPM_MULTISEGMENTWELL_EQUATIONS_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_EQUATIONS_HEADER_INCLUDED
24
25#include <dune/common/fmatrix.hh>
26#include <dune/common/fvector.hh>
27#include <dune/istl/bcrsmatrix.hh>
28#include <dune/istl/bvector.hh>
29
30#include <memory>
31
32namespace Dune {
33template<class M> class UMFPack;
34}
35
36namespace Opm
37{
38
39template<class Scalar, int numWellEq, int numEq> class MultisegmentWellEquationAccess;
40template<class Scalar> class MultisegmentWellGeneric;
41#if COMPILE_BDA_BRIDGE
43#endif
45template<class Scalar> class WellState;
46
47template<class Scalar, int numWellEq, int numEq>
49{
50public:
51 // sparsity pattern for the matrices
52 // [A C^T [x = [ res
53 // B D ] x_well] res_well]
54
55 // the vector type for the res_well and x_well
56 using VectorBlockWellType = Dune::FieldVector<Scalar,numWellEq>;
57 using BVectorWell = Dune::BlockVector<VectorBlockWellType>;
58
59 using VectorBlockType = Dune::FieldVector<Scalar,numEq>;
60 using BVector = Dune::BlockVector<VectorBlockType>;
61
62 // the matrix type for the diagonal matrix D
63 using DiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numWellEq>;
64 using DiagMatWell = Dune::BCRSMatrix<DiagMatrixBlockWellType>;
65
66 // the matrix type for the non-diagonal matrix B and C^T
67 using OffDiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numEq>;
68 using OffDiagMatWell = Dune::BCRSMatrix<OffDiagMatrixBlockWellType>;
69
71
78 void init(const int num_cells,
79 const int numPerfs,
80 const std::vector<int>& cells,
81 const std::vector<std::vector<int>>& segment_inlets,
82 const std::vector<std::vector<int>>& segment_perforations);
83
85 void clear();
86
88 void apply(const BVector& x, BVector& Ax) const;
89
91 void apply(BVector& r) const;
92
95
98
100 BVectorWell solve(const BVectorWell& rhs) const;
101
104 void recoverSolutionWell(const BVector& x, BVectorWell& xw) const;
105
106#if COMPILE_BDA_BRIDGE
108 void extract(WellContributions& wellContribs) const;
109#endif
110
112 template<class SparseMatrixAdapter>
113 void extract(SparseMatrixAdapter& jacobian) const;
114
116 template<class PressureMatrix>
117 void extractCPRPressureMatrix(PressureMatrix& jacobian,
118 const BVector& weights,
119 const int pressureVarIndex,
120 const bool /*use_well_weights*/,
121 const WellInterfaceGeneric& well,
122 const int seg_pressure_var_ind,
123 const WellState<Scalar>& well_state) const;
124
126 const BVectorWell& residual() const
127 {
128 return resWell_;
129 }
130
131 private:
132 friend class MultisegmentWellEquationAccess<Scalar,numWellEq,numEq>;
133 // two off-diagonal matrices
134 OffDiagMatWell duneB_;
135 OffDiagMatWell duneC_;
136 // "diagonal" matrix for the well. It has offdiagonal entries for inlets and outlets.
137 DiagMatWell duneD_;
138
142 mutable std::shared_ptr<Dune::UMFPack<DiagMatWell>> duneDSolver_;
143
144 // residuals of the well equations
145 BVectorWell resWell_;
146
148};
149
150}
151
152#endif // OPM_MULTISEGMENTWELLWELL_EQUATIONS_HEADER_INCLUDED
Definition: MultisegmentWellEquations.hpp:39
Definition: MultisegmentWellEquations.hpp:49
BVectorWell solve() const
Apply inverted D matrix to residual and return result.
Dune::BCRSMatrix< DiagMatrixBlockWellType > DiagMatWell
Definition: MultisegmentWellEquations.hpp:64
Dune::FieldMatrix< Scalar, numWellEq, numWellEq > DiagMatrixBlockWellType
Definition: MultisegmentWellEquations.hpp:63
Dune::FieldVector< Scalar, numEq > VectorBlockType
Definition: MultisegmentWellEquations.hpp:59
void extractCPRPressureMatrix(PressureMatrix &jacobian, const BVector &weights, const int pressureVarIndex, const bool, const WellInterfaceGeneric &well, const int seg_pressure_var_ind, const WellState< Scalar > &well_state) const
Extract CPR pressure matrix.
MultisegmentWellEquations(const MultisegmentWellGeneric< Scalar > &well)
Dune::BlockVector< VectorBlockType > BVector
Definition: MultisegmentWellEquations.hpp:60
Dune::BCRSMatrix< OffDiagMatrixBlockWellType > OffDiagMatWell
Definition: MultisegmentWellEquations.hpp:68
void clear()
Set all coefficients to 0.
Dune::BlockVector< VectorBlockWellType > BVectorWell
Definition: MultisegmentWellEquations.hpp:57
Dune::FieldVector< Scalar, numWellEq > VectorBlockWellType
Definition: MultisegmentWellEquations.hpp:56
void apply(const BVector &x, BVector &Ax) const
Apply linear operator to vector.
void extract(SparseMatrixAdapter &jacobian) const
Add the matrices of this well to the sparse matrix adapter.
void recoverSolutionWell(const BVector &x, BVectorWell &xw) const
Recover well solution.
void createSolver()
Compute the LU-decomposition of D matrix.
const BVectorWell & residual() const
Returns a const reference to the residual.
Definition: MultisegmentWellEquations.hpp:126
void apply(BVector &r) const
Apply linear operator to vector.
BVectorWell solve(const BVectorWell &rhs) const
Apply inverted D matrix to rhs and return result.
void init(const int num_cells, const int numPerfs, const std::vector< int > &cells, const std::vector< std::vector< int > > &segment_inlets, const std::vector< std::vector< int > > &segment_perforations)
Setup sparsity pattern for the matrices.
Dune::FieldMatrix< Scalar, numWellEq, numEq > OffDiagMatrixBlockWellType
Definition: MultisegmentWellEquations.hpp:67
Definition: MultisegmentWellGeneric.hpp:42
Definition: WellContributions.hpp:52
Definition: WellInterfaceGeneric.hpp:50
Definition: WellState.hpp:62
Definition: SupportsFaceTag.hpp:27
Definition: BlackoilPhases.hpp:27