StandardWellEquations.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
23#ifndef OPM_STANDARDWELL_EQUATIONS_HEADER_INCLUDED
24#define OPM_STANDARDWELL_EQUATIONS_HEADER_INCLUDED
25
28#include <opm/common/TimingMacros.hpp>
29#include <dune/common/dynmatrix.hh>
30#include <dune/common/dynvector.hh>
31#include <dune/istl/bcrsmatrix.hh>
32#include <dune/istl/bvector.hh>
33
34namespace Opm
35{
36
37class ParallelWellInfo;
38template<class Scalar, int numEq> class StandardWellEquationAccess;
39#if COMPILE_BDA_BRIDGE
41#endif
43template<class Scalar> class WellState;
44
45template<class Scalar, 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::DynamicVector<Scalar>;
55 using BVectorWell = Dune::BlockVector<VectorBlockWellType>;
56
57 // the matrix type for the diagonal matrix D
58 using DiagMatrixBlockWellType = Dune::DynamicMatrix<Scalar>;
59 using DiagMatWell = Dune::BCRSMatrix<DiagMatrixBlockWellType>;
60
61 // the matrix type for the non-diagonal matrix B and C^T
62 using OffDiagMatrixBlockWellType = Dune::DynamicMatrix<Scalar>;
63 using OffDiagMatWell = Dune::BCRSMatrix<OffDiagMatrixBlockWellType>;
64
65 // block vector type
66 using BVector = Dune::BlockVector<Dune::FieldVector<Scalar,numEq>>;
67
68 StandardWellEquations(const ParallelWellInfo& parallel_well_info);
69
75 void init(const int num_cells,
76 const int numWellEq,
77 const int numPerfs,
78 const std::vector<int>& cells);
79
81 void clear();
82
84 void apply(const BVector& x, BVector& Ax) const;
85
87 void apply(BVector& r) const;
88
90 void solve(BVectorWell& dx_well) const;
91
93 void solve(const BVectorWell& rhs_well, BVectorWell& x_well) const;
94
96 void invert();
97
100 void recoverSolutionWell(const BVector& x, BVectorWell& xw) const;
101
102#if COMPILE_BDA_BRIDGE
104 void extract(const int numStaticWellEq,
105 WellContributions& wellContribs) const;
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 bhp_var_index,
120 const WellState<Scalar>& well_state) const;
121
123 unsigned int getNumBlocks() const;
124
127
129 const BVectorWell& residual() const
130 {
131 return resWell_;
132 }
133
134private:
135 friend class StandardWellEquationAccess<Scalar,numEq>;
136
137 // two off-diagonal matrices
138 OffDiagMatWell duneB_;
139 OffDiagMatWell duneC_;
140 // diagonal matrix for the well
141 DiagMatWell invDuneD_;
142 DiagMatWell duneD_;
143
144 // Wrapper for the parallel application of B for distributed wells
146
147 // residuals of the well equations
148 BVectorWell resWell_;
149
150 // several vector used in the matrix calculation
151 mutable BVectorWell Bx_;
152 mutable BVectorWell invDrw_;
153};
154
155}
156
157#endif // OPM_STANDARDWELL_EQUATIONS_HEADER_INCLUDED
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:184
Definition: StandardWellEquations.hpp:38
Definition: StandardWellEquations.hpp:47
Dune::DynamicVector< Scalar > VectorBlockWellType
Definition: StandardWellEquations.hpp:54
Dune::DynamicMatrix< Scalar > DiagMatrixBlockWellType
Definition: StandardWellEquations.hpp:58
Dune::BCRSMatrix< OffDiagMatrixBlockWellType > OffDiagMatWell
Definition: StandardWellEquations.hpp:63
const BVectorWell & residual() const
Returns a const reference to the residual.
Definition: StandardWellEquations.hpp:129
void extract(SparseMatrixAdapter &jacobian) const
Add the matrices of this well to the sparse matrix adapter.
void solve(const BVectorWell &rhs_well, BVectorWell &x_well) const
Apply inverted D matrix to rhs and store in vector.
void extractCPRPressureMatrix(PressureMatrix &jacobian, const BVector &weights, const int pressureVarIndex, const bool use_well_weights, const WellInterfaceGeneric &well, const int bhp_var_index, const WellState< Scalar > &well_state) const
Extract CPR pressure matrix.
void invert()
Invert D matrix.
void apply(const BVector &x, BVector &Ax) const
Apply linear operator to vector.
Dune::BlockVector< VectorBlockWellType > BVectorWell
Definition: StandardWellEquations.hpp:55
Dune::BCRSMatrix< DiagMatrixBlockWellType > DiagMatWell
Definition: StandardWellEquations.hpp:59
void apply(BVector &r) const
Apply linear operator to vector.
void init(const int num_cells, const int numWellEq, const int numPerfs, const std::vector< int > &cells)
Setup sparsity pattern for the matrices.
unsigned int getNumBlocks() const
Get the number of blocks of the C and B matrices.
void sumDistributed(Parallel::Communication comm)
Sum with off-process contribution.
void solve(BVectorWell &dx_well) const
Apply inverted D matrix to residual and store in vector.
Dune::BlockVector< Dune::FieldVector< Scalar, numEq > > BVector
Definition: StandardWellEquations.hpp:66
void recoverSolutionWell(const BVector &x, BVectorWell &xw) const
Recover well solution.
Dune::DynamicMatrix< Scalar > OffDiagMatrixBlockWellType
Definition: StandardWellEquations.hpp:62
StandardWellEquations(const ParallelWellInfo &parallel_well_info)
void clear()
Set all coefficients to 0.
Definition: WellContributions.hpp:52
Definition: WellInterfaceGeneric.hpp:50
Definition: WellState.hpp:62
A wrapper around the B matrix for distributed wells.
Definition: WellHelpers.hpp:53
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: BlackoilPhases.hpp:27