WellHelpers.hpp
Go to the documentation of this file.
1/*
2 Copyright 2016 SINTEF ICT, Applied Mathematics.
3 Copyright 2016 Statoil ASA.
4 Copyright 2020 OPM-OP 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_WELLHELPERS_HEADER_INCLUDED
24#define OPM_WELLHELPERS_HEADER_INCLUDED
25
26#include <dune/istl/bcrsmatrix.hh>
27#include <dune/common/dynmatrix.hh>
28
29#include <array>
30
31namespace Opm {
32
33class ParallelWellInfo;
34struct WellProductionControls;
35struct WellInjectionControls;
36enum class WellProducerCMode;
37enum class WellInjectorCMode;
38
39namespace wellhelpers {
40
51template<typename Scalar>
53{
54public:
55 using Block = Dune::DynamicMatrix<Scalar>;
56 using Matrix = Dune::BCRSMatrix<Block>;
57
58 ParallelStandardWellB(const Matrix& B, const ParallelWellInfo& 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
68private:
69 const Matrix& B_;
70 const ParallelWellInfo& parallel_well_info_;
71};
72
73double computeHydrostaticCorrection(const double well_ref_depth,
74 const double vfp_ref_depth,
75 const double rho, const double gravity);
76
77
79template<typename Scalar, typename Comm>
80void sumDistributedWellEntries(Dune::DynamicMatrix<Scalar>& mat,
81 Dune::DynamicVector<Scalar>& vec,
82 const Comm& comm);
83
84
85// explicit transpose of a dense matrix due to compilation problems
86// used for calculating quasiimpes well weights
87template <class DenseMatrix>
88DenseMatrix transposeDenseDynMatrix(const DenseMatrix& M);
89
91bool rateControlWithZeroProdTarget(const WellProductionControls& controls,
92 WellProducerCMode mode);
93
95bool rateControlWithZeroInjTarget(const WellInjectionControls& controls,
96 WellInjectorCMode mode);
97
98} // namespace wellhelpers
99} // namespace Opm
100
101#endif
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:270
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:184
A wrapper around the B matrix for distributed wells.
Definition: WellHelpers.hpp:53
void mv(const X &x, Y &y) const
y = A x
void mmv(const X &x, Y &y) const
y = A x
Dune::DynamicMatrix< Scalar > Block
Definition: WellHelpers.hpp:55
ParallelStandardWellB(const Matrix &B, const ParallelWellInfo &parallel_well_info)
Dune::BCRSMatrix< Block > Matrix
Definition: WellHelpers.hpp:56
bool rateControlWithZeroInjTarget(const WellInjectionControls &controls, WellInjectorCMode mode)
Helper to check whether the well is under zero injection rate control.
double computeHydrostaticCorrection(const double well_ref_depth, const double vfp_ref_depth, const double rho, const double gravity)
bool rateControlWithZeroProdTarget(const WellProductionControls &controls, WellProducerCMode mode)
Helper to check whether the well is under zero production rate control.
DenseMatrix transposeDenseDynMatrix(const DenseMatrix &M)
void sumDistributedWellEntries(Dune::DynamicMatrix< Scalar > &mat, Dune::DynamicVector< Scalar > &vec, const Comm &comm)
Sums entries of the diagonal Matrix for distributed wells.
Definition: BlackoilPhases.hpp:27