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;
41#if COMPILE_BDA_BRIDGE
42class WellContributions;
43#endif
44class WellInterfaceGeneric;
45class 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
94 void createSolver();
95
97 BVectorWell solve() const;
98
101 void recoverSolutionWell(const BVector& x, BVectorWell& xw) const;
102
103#if COMPILE_BDA_BRIDGE
106#endif
107
109 template<class SparseMatrixAdapter>
110 void extract(SparseMatrixAdapter& jacobian) const;
111
113 template<class PressureMatrix>
114 void extractCPRPressureMatrix(PressureMatrix& jacobian,
115 const BVector& weights,
116 const int pressureVarIndex,
117 const bool /*use_well_weights*/,
118 const WellInterfaceGeneric& well,
119 const int seg_pressure_var_ind,
120 const WellState& well_state) const;
121
123 const BVectorWell& residual() const
124 {
125 return resWell_;
126 }
127
128 private:
129 friend class MultisegmentWellEquationAccess<Scalar,numWellEq,numEq>;
130 // two off-diagonal matrices
131 OffDiagMatWell duneB_;
132 OffDiagMatWell duneC_;
133 // "diagonal" matrix for the well. It has offdiagonal entries for inlets and outlets.
134 DiagMatWell duneD_;
135
139 mutable std::shared_ptr<Dune::UMFPack<DiagMatWell>> duneDSolver_;
140
141 // residuals of the well equations
142 BVectorWell resWell_;
143
145};
146
147}
148
149#endif // OPM_MULTISEGMENTWELLWELL_EQUATIONS_HEADER_INCLUDED
Definition AquiferInterface.hpp:35
Class administering assembler access to equation system.
Definition MultisegmentWellAssemble.cpp:45
Definition MultisegmentWellEquations.hpp:49
void clear()
Set all coefficients to 0.
Definition MultisegmentWellEquations.cpp:128
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:303
void apply(const BVector &x, BVector &Ax) const
Apply linear operator to vector.
Definition MultisegmentWellEquations.cpp:139
void extract(SparseMatrixAdapter &jacobian) const
Add the matrices of this well to the sparse matrix adapter.
Definition MultisegmentWellEquations.cpp:267
void recoverSolutionWell(const BVector &x, BVectorWell &xw) const
Recover well solution.
Definition MultisegmentWellEquations.cpp:186
void createSolver()
Compute the LU-decomposition of D matrix.
Definition MultisegmentWellEquations.cpp:163
BVectorWell solve() const
Apply inverted D matrix to residual and return result.
Definition MultisegmentWellEquations.cpp:179
const BVectorWell & residual() const
Returns a const reference to the residual.
Definition MultisegmentWellEquations.hpp:123
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:56
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