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
32#include <memory>
33
34namespace Opm {
35
36class DeferredLogger;
37class SICD;
38
39namespace mswellhelpers
40{
41
43 template <typename MatrixType, typename VectorType>
44 VectorType
46 VectorType x);
47
48
49
51 template <typename VectorType, typename MatrixType>
52 Dune::Matrix<typename MatrixType::block_type>
53 invertWithUMFPack(const int size,
54 const int bsize,
55 Dune::UMFPack<MatrixType>& linsolver);
56
57
58
59 // obtain y = D^-1 * x with a BICSSTAB iterative solver
60 template <typename MatrixType, typename VectorType>
61 VectorType
62 invDX(const MatrixType& D, VectorType x, DeferredLogger& deferred_logger);
63
64 // calculating the friction pressure loss
65 // l is the segment length
66 // area is the segment cross area
67 // diameter is the segment inner diameter
68 // w is mass flow rate through the segment
69 // density is density
70 // roughness is the absolute roughness
71 // mu is the average phase viscosity
72 template <typename ValueType>
73 ValueType frictionPressureLoss(const double l, const double diameter,
74 const double area, const double roughness,
75 const ValueType& density,
76 const ValueType& w, const ValueType& mu);
77
78
79 template <typename ValueType>
80 ValueType valveContrictionPressureLoss(const ValueType& mass_rate,
81 const ValueType& density,
82 const double area_con, const double cv);
83
84
85 template <typename ValueType>
86 ValueType velocityHead(const double area, const ValueType& mass_rate,
87 const ValueType& density);
88
89
90 // calculating the viscosity of oil-water emulsion at local conditons
91 template <typename ValueType>
92 ValueType emulsionViscosity(const ValueType& water_fraction,
93 const ValueType& water_viscosity,
94 const ValueType& oil_fraction,
95 const ValueType& oil_viscosity,
96 const SICD& sicd);
97
98} // namespace mswellhelpers
99} // namespace Opm
100
101#endif
Definition: MSWellHelpers.hpp:29
Definition: DeferredLogger.hpp:57
Definition: SupportsFaceTag.hpp:27
VectorType applyUMFPack(Dune::UMFPack< MatrixType > &linsolver, VectorType x)
Applies umfpack and checks for singularity.
Dune::Matrix< typename MatrixType::block_type > invertWithUMFPack(const int size, const int bsize, Dune::UMFPack< MatrixType > &linsolver)
Applies umfpack and checks for singularity.
ValueType valveContrictionPressureLoss(const ValueType &mass_rate, const ValueType &density, const double area_con, const double cv)
ValueType velocityHead(const double area, const ValueType &mass_rate, const ValueType &density)
ValueType frictionPressureLoss(const double l, const double diameter, const double area, const double roughness, const ValueType &density, const ValueType &w, const ValueType &mu)
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: BlackoilPhases.hpp:27