MSWellHelpers.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 Statoil ASA.
4 Copyright 2020 Equinor ASA.
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_MSWELLHELPERS_HEADER_INCLUDED
24#define OPM_MSWELLHELPERS_HEADER_INCLUDED
25
26#include <dune/istl/matrix.hh>
27
28namespace Dune {
29template<class Matrix> class UMFPack;
30}
31
32namespace Opm {
33
34template<class Scalar> class ParallelWellInfo;
35class DeferredLogger;
36class SICD;
37
38namespace mswellhelpers
39{
40
52 template<class MatrixType>
54 {
55 public:
56 using Scalar = typename MatrixType::field_type;
57 ParallellMSWellB(const MatrixType& B,
58 const ParallelWellInfo<Scalar>& parallel_well_info);
59
61 template<class X, class Y>
62 void mv (const X& x, Y& y) const;
63
65 template<class X, class Y>
66 void mmv (const X& x, Y& y) const;
67
68 private:
69 const MatrixType& B_;
70 const ParallelWellInfo<Scalar>& parallel_well_info_;
71 };
72
74 template <typename MatrixType, typename VectorType>
75 VectorType
77 VectorType x);
78
79
80
82 template <typename VectorType, typename MatrixType>
83 Dune::Matrix<typename MatrixType::block_type>
84 invertWithUMFPack(const int size,
85 const int bsize,
86 Dune::UMFPack<MatrixType>& linsolver);
87
88
89
90 // obtain y = D^-1 * x with a BICSSTAB iterative solver
91 template <typename MatrixType, typename VectorType>
92 VectorType
93 invDX(const MatrixType& D, VectorType x, DeferredLogger& deferred_logger);
94
95 // calculating the friction pressure loss
96 // l is the segment length
97 // area is the segment cross area
98 // diameter is the segment inner diameter
99 // w is mass flow rate through the segment
100 // density is density
101 // roughness is the absolute roughness
102 // mu is the average phase viscosity
103 template <typename ValueType, typename Scalar>
104 ValueType frictionPressureLoss(const Scalar l, const Scalar diameter,
105 const Scalar area, const Scalar roughness,
106 const ValueType& density,
107 const ValueType& w, const ValueType& mu);
108
109
110 template <typename ValueType, typename Scalar>
111 ValueType valveContrictionPressureLoss(const ValueType& mass_rate,
112 const ValueType& density,
113 const Scalar area_con, const Scalar cv);
114
115
116 template <typename ValueType, typename Scalar>
117 ValueType velocityHead(const Scalar area, const ValueType& mass_rate,
118 const ValueType& density);
119
120
121 // calculating the viscosity of oil-water emulsion at local conditons
122 template <typename ValueType, typename Scalar>
123 ValueType emulsionViscosity(const ValueType& water_fraction,
124 const ValueType& water_viscosity,
125 const ValueType& oil_fraction,
126 const ValueType& oil_viscosity,
127 const SICD& sicd);
128
129} // namespace mswellhelpers
130} // namespace Opm
131
132#endif
Definition: MSWellHelpers.hpp:29
Definition: DeferredLogger.hpp:57
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:195
A wrapper around the B matrix for distributed MS wells.
Definition: MSWellHelpers.hpp:54
typename MatrixType::field_type Scalar
Definition: MSWellHelpers.hpp:56
void mv(const X &x, Y &y) const
y = A x
ParallellMSWellB(const MatrixType &B, const ParallelWellInfo< Scalar > &parallel_well_info)
void mmv(const X &x, Y &y) const
y = A x
Definition: fvbaseprimaryvariables.hh:141
ValueType velocityHead(const Scalar area, const ValueType &mass_rate, const ValueType &density)
ValueType valveContrictionPressureLoss(const ValueType &mass_rate, const ValueType &density, const Scalar area_con, const Scalar cv)
VectorType applyUMFPack(Dune::UMFPack< MatrixType > &linsolver, VectorType x)
Applies umfpack and checks for singularity.
ValueType frictionPressureLoss(const Scalar l, const Scalar diameter, const Scalar area, const Scalar roughness, const ValueType &density, const ValueType &w, const ValueType &mu)
Dune::Matrix< typename MatrixType::block_type > invertWithUMFPack(const int size, const int bsize, Dune::UMFPack< MatrixType > &linsolver)
Applies umfpack and checks for singularity.
ValueType emulsionViscosity(const ValueType &water_fraction, const ValueType &water_viscosity, const ValueType &oil_fraction, const ValueType &oil_viscosity, const SICD &sicd)
VectorType invDX(const MatrixType &D, VectorType x, DeferredLogger &deferred_logger)
Definition: blackoilboundaryratevector.hh:39