My Project
Loading...
Searching...
No Matches
MultisegmentWellEquations.hpp
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;
41class WellContributions;
42class WellInterfaceGeneric;
43class WellState;
44
45template<class Scalar, int numWellEq, int numEq>
47{
48public:
49 // sparsity pattern for the matrices
50 // [A C^T [x = [ res
51 // B D ] x_well] res_well]
52
53 // the vector type for the res_well and x_well
54 using VectorBlockWellType = Dune::FieldVector<Scalar,numWellEq>;
55 using BVectorWell = Dune::BlockVector<VectorBlockWellType>;
56
57 using VectorBlockType = Dune::FieldVector<Scalar,numEq>;
58 using BVector = Dune::BlockVector<VectorBlockType>;
59
60 // the matrix type for the diagonal matrix D
61 using DiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numWellEq>;
62 using DiagMatWell = Dune::BCRSMatrix<DiagMatrixBlockWellType>;
63
64 // the matrix type for the non-diagonal matrix B and C^T
65 using OffDiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numEq>;
66 using OffDiagMatWell = Dune::BCRSMatrix<OffDiagMatrixBlockWellType>;
67
69
76 void init(const int num_cells,
77 const int numPerfs,
78 const std::vector<int>& cells,
79 const std::vector<std::vector<int>>& segment_inlets,
80 const std::vector<std::vector<int>>& segment_perforations);
81
83 void clear();
84
86 void apply(const BVector& x, BVector& Ax) const;
87
89 void apply(BVector& r) const;
90
92 void createSolver();
93
95 BVectorWell solve() const;
96
99 void recoverSolutionWell(const BVector& x, BVectorWell& xw) const;
100
103
105 template<class SparseMatrixAdapter>
106 void extract(SparseMatrixAdapter& jacobian) const;
107
109 template<class PressureMatrix>
110 void extractCPRPressureMatrix(PressureMatrix& jacobian,
111 const BVector& weights,
112 const int pressureVarIndex,
113 const bool /*use_well_weights*/,
114 const WellInterfaceGeneric& well,
115 const int seg_pressure_var_ind,
116 const WellState& well_state) const;
117
119 const BVectorWell& residual() const
120 {
121 return resWell_;
122 }
123
124 private:
125 friend class MultisegmentWellEquationAccess<Scalar,numWellEq,numEq>;
126 // two off-diagonal matrices
127 OffDiagMatWell duneB_;
128 OffDiagMatWell duneC_;
129 // "diagonal" matrix for the well. It has offdiagonal entries for inlets and outlets.
130 DiagMatWell duneD_;
131
135 mutable std::shared_ptr<Dune::UMFPack<DiagMatWell>> duneDSolver_;
136
137 // residuals of the well equations
138 BVectorWell resWell_;
139
141};
142
143}
144
145#endif // OPM_MULTISEGMENTWELLWELL_EQUATIONS_HEADER_INCLUDED
Definition AquiferInterface.hpp:35
Class administering assembler access to equation system.
Definition MultisegmentWellAssemble.cpp:45
Definition MultisegmentWellEquations.hpp:47
void extract(WellContributions &wellContribs) const
Add the matrices of this well to the WellContributions object.
Definition MultisegmentWellEquations.cpp:194
void clear()
Set all coefficients to 0.
Definition MultisegmentWellEquations.cpp:125
void extractCPRPressureMatrix(PressureMatrix &jacobian, const BVector &weights, const int pressureVarIndex, const bool, const WellInterfaceGeneric &well, const int seg_pressure_var_ind, const WellState &well_state) const
Extract CPR pressure matrix.
Definition MultisegmentWellEquations.cpp:299
void apply(const BVector &x, BVector &Ax) const
Apply linear operator to vector.
Definition MultisegmentWellEquations.cpp:136
void recoverSolutionWell(const BVector &x, BVectorWell &xw) const
Recover well solution.
Definition MultisegmentWellEquations.cpp:183
void createSolver()
Compute the LU-decomposition of D matrix.
Definition MultisegmentWellEquations.cpp:160
BVectorWell solve() const
Apply inverted D matrix to residual and return result.
Definition MultisegmentWellEquations.cpp:176
const BVectorWell & residual() const
Returns a const reference to the residual.
Definition MultisegmentWellEquations.hpp:119
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.
Definition MultisegmentWellEquations.cpp:53
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:52
Definition WellInterfaceGeneric.hpp:51
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition WellState.hpp:60
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27