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
28#include <dune/common/fmatrix.hh>
29#include <dune/common/fvector.hh>
30#include <dune/istl/bcrsmatrix.hh>
31#include <dune/istl/bvector.hh>
32
33#include <memory>
34
35namespace Dune {
36template<class M> class UMFPack;
37}
38
39namespace Opm
40{
41
42template<class Scalar, typename IndexTraits, int numWellEq, int numEq> class MultisegmentWellEquationAccess;
43template<class Scalar, typename IndexTraits> class MultisegmentWellGeneric;
44#if COMPILE_GPU_BRIDGE
45template<class Scalar> class WellContributions;
46#endif
47template<typename Scalar, typename IndexTraits> class WellInterfaceGeneric;
48template<typename Scalar, typename IndexTraits> class WellState;
49
50template<class Scalar, typename IndexTraits, int numWellEq, int numEq>
52{
53public:
54 // sparsity pattern for the matrices
55 // [A C^T [x = [ res
56 // B D ] x_well] res_well]
57
58 // the vector type for the res_well and x_well
59 using VectorBlockWellType = Dune::FieldVector<Scalar,numWellEq>;
60 using BVectorWell = Dune::BlockVector<VectorBlockWellType>;
61
62 using VectorBlockType = Dune::FieldVector<Scalar,numEq>;
63 using BVector = Dune::BlockVector<VectorBlockType>;
64
65 // the matrix type for the diagonal matrix D
66 using DiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numWellEq>;
67 using DiagMatWell = Dune::BCRSMatrix<DiagMatrixBlockWellType>;
68
69 // the matrix type for the non-diagonal matrix B and C^T
70 using OffDiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numEq>;
71 using OffDiagMatWell = Dune::BCRSMatrix<OffDiagMatrixBlockWellType>;
72
74
80 void init(const int numPerfs,
81 const std::vector<int>& cells,
82 const std::vector<std::vector<int>>& segment_inlets,
83 const std::vector<std::vector<int>>& segment_perforations);
84
86 void clear();
87
89 void apply(const BVector& x, BVector& Ax) const;
90
92 void apply(BVector& r) const;
93
96
99
101 BVectorWell solve(const BVectorWell& rhs) const;
102
105 void recoverSolutionWell(const BVector& x, BVectorWell& xw) const;
106
107#if COMPILE_GPU_BRIDGE
109 void extract(WellContributions<Scalar>& wellContribs) const;
110#endif
111
113 template<class SparseMatrixAdapter>
114 void extract(SparseMatrixAdapter& jacobian) const;
115
117 template<class PressureMatrix>
118 void extractCPRPressureMatrix(PressureMatrix& jacobian,
119 const BVector& weights,
120 const int pressureVarIndex,
121 const bool /*use_well_weights*/,
123 const int seg_pressure_var_ind,
124 const WellState<Scalar, IndexTraits>& well_state) const;
125
128
130 const BVectorWell& residual() const
131 {
132 return resWell_;
133 }
134
135 private:
136 friend class MultisegmentWellEquationAccess<Scalar,IndexTraits,numWellEq,numEq>;
137 // two off-diagonal matrices
138 OffDiagMatWell duneB_;
139 OffDiagMatWell duneC_;
140 // "diagonal" matrix for the well. It has offdiagonal entries for inlets and outlets.
141 DiagMatWell duneD_;
142
146 mutable std::shared_ptr<Dune::UMFPack<DiagMatWell>> duneDSolver_;
147
148 // residuals of the well equations
149 BVectorWell resWell_;
150
152
153 // Store the global index of well perforated cells
154 std::vector<int> cells_;
155
156 const ParallelWellInfo<Scalar>& pw_info_;
157
158 // Wrapper for the parallel application of B for distributed wells
160};
161
162}
163
164#endif // OPM_MULTISEGMENTWELLWELL_EQUATIONS_HEADER_INCLUDED
Definition: MultisegmentWellEquations.hpp:42
Definition: MultisegmentWellEquations.hpp:52
Dune::FieldVector< Scalar, numWellEq > VectorBlockWellType
Definition: MultisegmentWellEquations.hpp:59
void sumDistributed(Parallel::Communication comm)
Sum with off-process contribution.
BVectorWell solve(const BVectorWell &rhs) const
Apply inverted D matrix to rhs and return result.
Dune::BCRSMatrix< DiagMatrixBlockWellType > DiagMatWell
Definition: MultisegmentWellEquations.hpp:67
MultisegmentWellEquations(const MultisegmentWellGeneric< Scalar, IndexTraits > &well, const ParallelWellInfo< Scalar > &pw_info)
Dune::FieldMatrix< Scalar, numWellEq, numEq > OffDiagMatrixBlockWellType
Definition: MultisegmentWellEquations.hpp:70
void apply(BVector &r) const
Apply linear operator to vector.
const BVectorWell & residual() const
Returns a const reference to the residual.
Definition: MultisegmentWellEquations.hpp:130
Dune::BCRSMatrix< OffDiagMatrixBlockWellType > OffDiagMatWell
Definition: MultisegmentWellEquations.hpp:71
void extractCPRPressureMatrix(PressureMatrix &jacobian, const BVector &weights, const int pressureVarIndex, const bool, const WellInterfaceGeneric< Scalar, IndexTraits > &well, const int seg_pressure_var_ind, const WellState< Scalar, IndexTraits > &well_state) const
Extract CPR pressure matrix.
BVectorWell solve() const
Apply inverted D matrix to residual and return result.
Dune::BlockVector< VectorBlockWellType > BVectorWell
Definition: MultisegmentWellEquations.hpp:60
void extract(SparseMatrixAdapter &jacobian) const
Add the matrices of this well to the sparse matrix adapter.
void clear()
Set all coefficients to 0.
void apply(const BVector &x, BVector &Ax) const
Apply linear operator to vector.
void createSolver()
Compute the LU-decomposition of D matrix.
Dune::FieldMatrix< Scalar, numWellEq, numWellEq > DiagMatrixBlockWellType
Definition: MultisegmentWellEquations.hpp:66
Dune::FieldVector< Scalar, numEq > VectorBlockType
Definition: MultisegmentWellEquations.hpp:62
void init(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::BlockVector< VectorBlockType > BVector
Definition: MultisegmentWellEquations.hpp:63
void recoverSolutionWell(const BVector &x, BVectorWell &xw) const
Recover well solution.
Definition: MultisegmentWellGeneric.hpp:39
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
Definition: WellContributions.hpp:51
Definition: WellInterfaceGeneric.hpp:53
Definition: WellState.hpp:66
Definition: fvbaseprimaryvariables.hh:141
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilboundaryratevector.hh:39