WellInterfaceIndices.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 Statoil ASA.
4 Copyright 2017 IRIS
5 Copyright 2019 Norce
6
7 This file is part of the Open Porous Media project (OPM).
8
9 OPM is free software: you can redistribute it and/or modify
10 it under the terms of the GNU General Public License as published by
11 the Free Software Foundation, either version 3 of the License, or
12 (at your option) any later version.
13
14 OPM is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18
19 You should have received a copy of the GNU General Public License
20 along with OPM. If not, see <http://www.gnu.org/licenses/>.
21*/
22
23#ifndef OPM_WELLINTERFACE_INDICES_HEADER_INCLUDED
24#define OPM_WELLINTERFACE_INDICES_HEADER_INCLUDED
25
26#include <opm/material/densead/Evaluation.hpp>
27
29
30namespace Opm {
31
32template<class FluidSystem, class Indices>
34{
35public:
36 using WellInterfaceFluidSystem<FluidSystem>::Gas;
37 using WellInterfaceFluidSystem<FluidSystem>::Oil;
38 using WellInterfaceFluidSystem<FluidSystem>::Water;
39 using Scalar = typename FluidSystem::Scalar;
40 using Eval = DenseAd::Evaluation<Scalar, /*size=*/Indices::numEq>;
42
43 int flowPhaseToModelCompIdx(const int phaseIdx) const;
44 int modelCompIdxToFlowCompIdx(const int compIdx) const;
45 int flowPhaseToModelPhaseIdx(const int phaseIdx) const;
46 Scalar scalingFactor(const int phaseIdx) const;
47
48 template <class EvalWell>
49 Eval restrictEval(const EvalWell& in) const
50 {
51 Eval out = 0.0;
52 out.setValue(in.value());
53 for (int eqIdx = 0; eqIdx < Indices::numEq; ++eqIdx) {
54 out.setDerivative(eqIdx, in.derivative(eqIdx));
55 }
56 return out;
57 }
58
59protected:
60 WellInterfaceIndices(const Well& well,
61 const ParallelWellInfo<Scalar>& parallel_well_info,
62 const int time_step,
63 const ModelParameters& param,
65 const int pvtRegionIdx,
66 const int num_components,
67 const int num_phases,
68 const int index_of_well,
69 const std::vector<PerforationData<Scalar>>& perf_data);
70};
71
72}
73
74#endif // OPM_WELLINTERFACE_INDICES_HEADER_INCLUDED
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:195
Definition: RateConverter.hpp:71
Definition: WellInterfaceFluidSystem.hpp:51
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:63
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:62
typename FluidSystem::Scalar Scalar
Definition: WellInterfaceFluidSystem.hpp:59
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:64
typename WellInterfaceGeneric< Scalar >::ModelParameters ModelParameters
Definition: WellInterfaceFluidSystem.hpp:60
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:124
Definition: WellInterfaceIndices.hpp:34
int flowPhaseToModelCompIdx(const int phaseIdx) const
typename FluidSystem::Scalar Scalar
Definition: WellInterfaceIndices.hpp:39
Eval restrictEval(const EvalWell &in) const
Definition: WellInterfaceIndices.hpp:49
int flowPhaseToModelPhaseIdx(const int phaseIdx) const
int modelCompIdxToFlowCompIdx(const int compIdx) const
Scalar scalingFactor(const int phaseIdx) const
DenseAd::Evaluation< Scalar, Indices::numEq > Eval
Definition: WellInterfaceIndices.hpp:40
WellInterfaceIndices(const Well &well, const ParallelWellInfo< Scalar > &parallel_well_info, const int time_step, const ModelParameters &param, const typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType &rate_converter, const int pvtRegionIdx, const int num_components, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Definition: blackoilboundaryratevector.hh:39
Static data associated with a well perforation.
Definition: PerforationData.hpp:30