20#ifndef OPM_BLACKOILFLUID_HEADER_INCLUDED 
   21#define OPM_BLACKOILFLUID_HEADER_INCLUDED 
   27#include <dune/common/fvector.hh> 
   36    class BlackoilFluidData;
 
   47        void init(Opm::DeckConstPtr deck)
 
   49            fmi_params_.
init(deck);
 
   52            const auto& densityRecord = deck->getKeyword(
"DENSITY").getRecord(0);
 
   53            surface_densities_[
Oil] = densityRecord.getItem(
"OIL").getSIDouble(0);
 
   54            surface_densities_[
Water] = densityRecord.getItem(
"WATER").getSIDouble(0);
 
   55            surface_densities_[
Gas] = densityRecord.getItem(
"GAS").getSIDouble(0);
 
   65            computeEquilibrium(state);
 
   68            for (
int phase = 0; phase < 
numPhases; ++phase) {
 
   81            return surface_densities_;
 
   88            for (
int phase = 0; phase < 
numPhases; ++phase) {
 
   90                    phase_dens[phase] += A[
numComponents*phase + comp]*surface_densities_[comp];
 
  101        template <
class States>
 
  104            const std::vector<PhaseVec>& p = states.phase_pressure;
 
  105            const std::vector<CompVec>& z = states.surface_volume_density;
 
  106            std::vector<PhaseVec>& B = states.formation_volume_factor;
 
  107            std::vector<PhaseVec>& R = states.solution_factor;
 
  114        template <
class States>
 
  118            const std::vector<PhaseVec>& p = states.phase_pressure;
 
  119            const std::vector<CompVec>& z = states.surface_volume_density;
 
  120            std::vector<PhaseVec>& mu = states.viscosity;
 
  126        template <
class States>
 
  129            const std::vector<PhaseVec>& p = states.phase_pressure;
 
  130            const std::vector<CompVec>& z = states.surface_volume_density;
 
  131            std::vector<PhaseVec>& B = states.formation_volume_factor;
 
  132            std::vector<PhaseVec>& dB = states.formation_volume_factor_deriv;
 
  133            std::vector<PhaseVec>& R = states.solution_factor;
 
  134            std::vector<PhaseVec>& dR = states.solution_factor_deriv;
 
  135            std::vector<PhaseVec>& mu = states.viscosity;
 
  136            pvt_.
dBdp(p, z, B, dB);
 
  137            pvt_.
dRdp(p, z, R, dR);
 
  143        template <
class States>
 
  146            int num = states.formation_volume_factor.size();
 
  147            states.state_matrix.resize(num);
 
  148#pragma omp parallel for 
  149            for (
int i = 0; i < num; ++i) {
 
  150                const PhaseVec& B = states.formation_volume_factor[i];
 
  151                const PhaseVec& R = states.solution_factor[i];
 
  168        template <
class States>
 
  171            int num = states.formation_volume_factor.size();
 
  172            states.state_matrix.resize(num);
 
  173            states.phase_volume_density.resize(num);
 
  174            states.total_phase_volume_density.resize(num);
 
  175            states.saturation.resize(num);
 
  176            states.phase_compressibility.resize(num);
 
  177            states.total_compressibility.resize(num);
 
  178            states.experimental_term.resize(num);
 
  179#pragma omp parallel for 
  180            for (
int i = 0; i < num; ++i) {
 
  181                const CompVec& z = states.surface_volume_density[i];
 
  182                const PhaseVec& B = states.formation_volume_factor[i];
 
  183                const PhaseVec& dB = states.formation_volume_factor_deriv[i];
 
  184                const PhaseVec& R = states.solution_factor[i];
 
  185                const PhaseVec& dR = states.solution_factor_deriv[i];
 
  187                PhaseVec& u = states.phase_volume_density[i];
 
  188                double& tot_phase_vol_dens = states.total_phase_volume_density[i];
 
  190                PhaseVec& cp = states.phase_compressibility[i];
 
  191                double& tot_comp = states.total_compressibility[i];
 
  192                double& exp_term = states.experimental_term[i];
 
  193                computeSingleEquilibrium(B, dB, R, dR, z,
 
  194                                         At, u, tot_phase_vol_dens,
 
  195                                         s, cp, tot_comp, exp_term);
 
  201        template <
class States>
 
  204            int num = states.saturation.size();
 
  205            states.relperm.resize(num);
 
  206            states.mobility.resize(num);
 
  207#pragma omp parallel for 
  208            for (
int i = 0; i < num; ++i) {
 
  209                const CompVec& s = states.saturation[i];
 
  210                const PhaseVec& mu = states.viscosity[i];
 
  212                PhaseVec& lambda = states.mobility[i];
 
  214                for (
int phase = 0; phase < 
numPhases; ++phase) {
 
  215                    lambda[phase] = kr[phase]/mu[phase];
 
  223        template <
class States>
 
  226            int num = states.saturation.size();
 
  227            states.relperm.resize(num);
 
  228            states.relperm_deriv.resize(num);
 
  229            states.mobility.resize(num);
 
  230            states.mobility_deriv.resize(num);
 
  231#pragma omp parallel for 
  232            for (
int i = 0; i < num; ++i) {
 
  233                const CompVec& s = states.saturation[i];
 
  234                const PhaseVec& mu = states.viscosity[i];
 
  237                PhaseVec& lambda = states.mobility[i];
 
  241                for (
int phase = 0; phase < 
numPhases; ++phase) {
 
  242                    lambda[phase] = kr[phase]/mu[phase];
 
  245                        dlambda[phase][p2] = dkr[phase][p2]/mu[phase];
 
  264    template <
class Flu
idState>
 
  265    void computeEquilibrium(
FluidState& fluid_state)
 const 
  268        const PhaseVec& p = fluid_state.phase_pressure_;
 
  269        const CompVec& z = fluid_state.surface_volume_;
 
  270        PhaseVec& B = fluid_state.formation_volume_factor_;
 
  274        PhaseVec& R = fluid_state.solution_factor_; 
 
  289        PhaseVec& u = fluid_state.phase_volume_density_;
 
  290        double& tot_phase_vol_dens = fluid_state.total_phase_volume_density_;
 
  291        PhaseVec& s = fluid_state.saturation_;
 
  292        PhaseVec& cp = fluid_state.phase_compressibility_;
 
  293        double& tot_comp = fluid_state.total_compressibility_;
 
  294        double& exp_term = fluid_state.experimental_term_;
 
  296        computeSingleEquilibrium(B, dB, R, dR, z,
 
  297                                 At, u, tot_phase_vol_dens,
 
  298                                 s, cp, tot_comp, exp_term);
 
  301        PhaseVec& mu = fluid_state.viscosity_;
 
  309    static void computeSingleEquilibrium(
const PhaseVec& B,
 
  316                                         double& tot_phase_vol_dens,
 
  343        for (
int phase = 0; phase < 
numPhases; ++phase) {
 
  344            s[phase] = u[phase]/tot_phase_vol_dens;
 
  358        Dune::FMatrixHelp::invertMatrix(At, Ait);
 
  361        Dune::FMatrixHelp::multMatrix(dAt, Ait, Ct);
 
  416        template <
class Gr
id, 
class Rock>
 
  420                        const typename Grid::Vector gravity,
 
  421                        const std::vector<PhaseVec>& cell_pressure,
 
  422                        const std::vector<PhaseVec>& face_pressure,
 
  423                        const std::vector<CompVec>& cell_z,
 
  427            int num_cells = cell_z.size();
 
  428            assert(num_cells == grid.numCells());
 
  431            static_assert(np == nc, 
"");
 
  432            static_assert(np == 3, 
"");
 
  444#pragma omp parallel for 
  445            for (
int cell = 0; cell < num_cells; ++cell) {
 
  446                double pv = rock.
porosity(cell)*grid.cellVolume(cell);
 
  454                                    cell_pressure, face_pressure,
 
  465        template <
class Gr
id>
 
  468                                     const typename Grid::Vector gravity,
 
  469                                     const std::vector<PhaseVec>& cell_pressure,
 
  470                                     const std::vector<PhaseVec>& face_pressure,
 
  471                                     const std::vector<CompVec>& cell_z,
 
  474            int num_faces = face_pressure.size();
 
  475            assert(num_faces == grid.numFaces());
 
  476            bool nonzero_gravity = gravity.two_norm() > 0.0;
 
  482#pragma omp parallel for 
  483            for (
int face = 0; face < num_faces; ++face) {
 
  485                typedef typename Grid::Vector Vec;
 
  486                Vec fc = grid.faceCentroid(face);
 
  487                int c[2] = { grid.faceCell(face, 0), grid.faceCell(face, 1) };
 
  493                for (
int j = 0; j < 2; ++j) {
 
  496                        phase_p[j] = cell_pressure[c[j]];
 
  498                        if (nonzero_gravity) {
 
  500                            cdiff -= grid.cellCentroid(c[j]);
 
  502                            gravcontrib[j] *= (cdiff*gravity);
 
  504                            gravcontrib[j] = 0.0;
 
  508                        phase_p[j] = face_pressure[face];
 
  510                        gravcontrib[j] = 0.0;
 
  520                for (
int phase = 0; phase < 
numPhases; ++phase) {
 
  523                    pot[1][phase] = phase_p[1][phase];
 
  532                double face_z_factor = 0.5;
 
  535                for (
int j = 0; j < 2; ++j) {
 
  537                        face_z += cell_z[c[j]];
 
  553                face_z *= face_z_factor;
 
  557                for (
int phase = 0; phase < 
numPhases; ++phase) {
 
  558                    if (pot[0][phase] == pot[1][phase]) {
 
  560                        double aver = 0.5*(phase_mob[0][phase] + phase_mob[1][phase]);
 
  564                                + phasemob_deriv[1][phase][p2];
 
  568                        int upwind = pot[0][phase] > pot[1][phase] ? 0 : 1;
 
Definition: BlackoilDefs.hpp:31
 
Dune::FieldVector< Scalar, numComponents > CompVec
Definition: BlackoilDefs.hpp:40
 
Dune::FieldVector< Scalar, numPhases > PhaseVec
Definition: BlackoilDefs.hpp:41
 
Dune::FieldMatrix< Scalar, numComponents, numPhases > PhaseToCompMatrix
Definition: BlackoilDefs.hpp:43
 
@ Aqua
Definition: BlackoilDefs.hpp:37
 
@ Vapour
Definition: BlackoilDefs.hpp:37
 
@ Liquid
Definition: BlackoilDefs.hpp:37
 
@ numComponents
Definition: BlackoilDefs.hpp:33
 
Dune::FieldMatrix< Scalar, numPhases, numPhases > PhaseJacobian
Definition: BlackoilDefs.hpp:44
 
@ numPhases
Definition: BlackoilDefs.hpp:34
 
@ Gas
Definition: BlackoilDefs.hpp:36
 
@ Oil
Definition: BlackoilDefs.hpp:36
 
@ Water
Definition: BlackoilDefs.hpp:36
 
Definition: BlackoilFluid.hpp:42
 
void computePvtNoDerivs(States &states) const
Definition: BlackoilFluid.hpp:115
 
void init(Opm::DeckConstPtr deck)
Definition: BlackoilFluid.hpp:47
 
void computeBAndR(States &states) const
Definition: BlackoilFluid.hpp:102
 
const CompVec & surfaceDensities() const
Definition: BlackoilFluid.hpp:79
 
PhaseVec phaseDensities(const double *A) const
Definition: BlackoilFluid.hpp:85
 
void computePvt(States &states) const
Definition: BlackoilFluid.hpp:127
 
void computeMobilities(States &states) const
Definition: BlackoilFluid.hpp:224
 
void computePvtDepending(States &states) const
Definition: BlackoilFluid.hpp:169
 
FluidState computeState(PhaseVec phase_pressure, CompVec z) const
Definition: BlackoilFluid.hpp:58
 
BlackoilFluidData FluidData
Definition: BlackoilFluid.hpp:45
 
void computeStateMatrix(States &states) const
Definition: BlackoilFluid.hpp:144
 
FluidStateBlackoil FluidState
Definition: BlackoilFluid.hpp:44
 
void computeMobilitiesNoDerivs(States &states) const
Definition: BlackoilFluid.hpp:202
 
Definition: BlackoilPVT.hpp:34
 
double getViscosity(double press, const CompVec &surfvol, PhaseIndex phase) const
 
double R(double press, const CompVec &surfvol, PhaseIndex phase) const
 
void init(Opm::DeckConstPtr deck)
 
double B(double press, const CompVec &surfvol, PhaseIndex phase) const
 
double dBdp(double press, const CompVec &surfvol, PhaseIndex phase) const
 
double dRdp(double press, const CompVec &surfvol, PhaseIndex phase) const
 
void init(Opm::DeckConstPtr deck)
Definition: FluidMatrixInteractionBlackoil.hpp:47
 
static void kr(krContainerT &kr, const Params ¶ms, const SatContainerT &saturations, Scalar)
The relative permeability of all phases.
Definition: FluidMatrixInteractionBlackoil.hpp:153
 
static void dkr(krContainerT &dkr, const Params ¶ms, const SatContainerT &saturations, Scalar)
The saturation derivatives of relative permeability of all phases.
Definition: FluidMatrixInteractionBlackoil.hpp:179
 
A property class for porous media rock.
Definition: ImplicitTransportDefs.hpp:80
 
double porosity(int cell_index) const
Read-access to porosity.
Definition: Rock_impl.hpp:76
 
Definition: BlackoilFluid.hpp:32
 
Definition: BlackoilFluid.hpp:406
 
void computeUpwindProperties(const Grid &grid, const BlackoilFluid &fluid, const typename Grid::Vector gravity, const std::vector< PhaseVec > &cell_pressure, const std::vector< PhaseVec > &face_pressure, const std::vector< CompVec > &cell_z, const CompVec &bdy_z)
Definition: BlackoilFluid.hpp:466
 
FaceFluidData face_data
Definition: BlackoilFluid.hpp:413
 
void computeNew(const Grid &grid, const Rock &rock, const BlackoilFluid &fluid, const typename Grid::Vector gravity, const std::vector< PhaseVec > &cell_pressure, const std::vector< PhaseVec > &face_pressure, const std::vector< CompVec > &cell_z, const CompVec &bdy_z, const double dt)
Definition: BlackoilFluid.hpp:417
 
std::vector< double > voldiscr
Definition: BlackoilFluid.hpp:409
 
AllFluidStates cell_data
Definition: BlackoilFluid.hpp:408
 
std::vector< double > relvoldiscr
Definition: BlackoilFluid.hpp:410
 
Multiple fluid states for a black oil model.
Definition: FluidStateBlackoil.hpp:58
 
std::vector< PhaseJacobian > mobility_deriv
Definition: FluidStateBlackoil.hpp:85
 
std::vector< PhaseVec > phase_pressure
Definition: FluidStateBlackoil.hpp:61
 
std::vector< Scalar > total_phase_volume_density
Definition: FluidStateBlackoil.hpp:75
 
std::vector< PhaseToCompMatrix > state_matrix
Definition: FluidStateBlackoil.hpp:73
 
std::vector< PhaseVec > mobility
Definition: FluidStateBlackoil.hpp:84
 
std::vector< CompVec > surface_volume_density
Definition: FluidStateBlackoil.hpp:60
 
Definition: BlackoilFluid.hpp:382
 
std::vector< PhaseVec > gravity_potential
Definition: BlackoilFluid.hpp:401
 
std::vector< PhaseVec > solution_factor
Definition: BlackoilFluid.hpp:389
 
std::vector< PhaseVec > formation_volume_factor
Definition: BlackoilFluid.hpp:388
 
std::vector< PhaseVec > phase_pressure
Definition: BlackoilFluid.hpp:385
 
std::vector< PhaseVec > mobility
Definition: BlackoilFluid.hpp:397
 
std::vector< CompVec > surface_volume_density
Definition: BlackoilFluid.hpp:384
 
std::vector< PhaseJacobian > mobility_deriv
Definition: BlackoilFluid.hpp:398
 
std::vector< PhaseToCompMatrix > state_matrix
Definition: BlackoilFluid.hpp:394
 
Fluid state for a black oil model.
Definition: FluidStateBlackoil.hpp:33
 
CompVec surface_volume_
Definition: FluidStateBlackoil.hpp:35
 
PhaseVec mobility_
Definition: FluidStateBlackoil.hpp:49
 
PhaseVec saturation_
Definition: FluidStateBlackoil.hpp:42
 
PhaseVec phase_pressure_
Definition: FluidStateBlackoil.hpp:36
 
Scalar temperature_
Definition: FluidStateBlackoil.hpp:34
 
PhaseJacobian dmobility_
Definition: FluidStateBlackoil.hpp:50
 
PhaseJacobian drelperm_
Definition: FluidStateBlackoil.hpp:48
 
PhaseVec viscosity_
Definition: FluidStateBlackoil.hpp:46
 
PhaseVec relperm_
Definition: FluidStateBlackoil.hpp:47