Opm Namespace Reference

Namespaces

namespace  detail
 
namespace  parameter
 
namespace  polymer_reorder
 

Classes

class  BlackoilPolymerModel
 
struct  BlackoilPolymerSolutionState
 
class  CompressibleTpfaPolymer
 
class  FullyImplicitCompressiblePolymerSolver
 
class  GravityColumnSolverPolymer
 
class  IncompPropertiesDefaultPolymer
 
class  IncompTpfaPolymer
 
struct  ModelTraits< BlackoilPolymerModel< Grid > >
 Providing types by template specialisation of ModelTraits for BlackoilPolymerModel. More...
 
class  PolymerBlackoilState
 
class  PolymerInflowBasic
 Basic polymer injection behaviour class. This class gives all injectors the same polymer concentration, during a single time interval. Amount and interval can be specified. More...
 
class  PolymerInflowFromDeck
 Polymer injection behaviour class using deck WPOLYMER. This class reads the accumulated WPOLYMER lines from the deck, and applies the last row given for each well. More...
 
class  PolymerInflowInterface
 Interface for classes encapsulating polymer inflow information. More...
 
class  PolymerProperties
 
class  PolymerPropsAd
 
class  PolymerState
 Simulator state for a two-phase simulator with polymer. More...
 
class  SimulatorCompressiblePolymer
 Class collecting all necessary components for a two-phase simulation. More...
 
class  SimulatorFullyImplicitBlackoilPolymer
 
class  SimulatorFullyImplicitCompressiblePolymer
 Class collecting all necessary components for a two-phase simulation. More...
 
class  SimulatorPolymer
 Class collecting all necessary components for a two-phase simulation. More...
 
struct  SimulatorTraits< SimulatorFullyImplicitBlackoilPolymer< GridT > >
 
struct  SimulatorTraits< SimulatorFullyImplicitCompressiblePolymer< GridT > >
 
class  SinglePointUpwindTwoPhasePolymer
 
class  TransportSolverTwophaseCompressiblePolymer
 
class  TransportSolverTwophasePolymer
 
class  WellStateFullyImplicitBlackoilPolymer
 

Functions

void computeTotalMobility (const Opm::IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< int > &cells, const std::vector< double > &s, const std::vector< double > &c, const std::vector< double > &cmax, std::vector< double > &totmob)
 Computes total mobility for a set of s/c values. More...
 
void computeTotalMobilityOmega (const Opm::IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< int > &cells, const std::vector< double > &s, const std::vector< double > &c, const std::vector< double > &cmax, std::vector< double > &totmob, std::vector< double > &omega)
 Computes total mobility and omega for a set of s/c values. More...
 
void computeFractionalFlow (const Opm::IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< int > &cells, const std::vector< double > &s, const std::vector< double > &c, const std::vector< double > &cmax, std::vector< double > &fractional_flows)
 
void computeFractionalFlow (const Opm::BlackoilPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< int > &cells, const std::vector< double > &p, const std::vector< double > &T, const std::vector< double > &z, const std::vector< double > &s, const std::vector< double > &c, const std::vector< double > &cmax, std::vector< double > &fractional_flows)
 
void computeInjectedProduced (const IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const PolymerState &state, const std::vector< double > &transport_src, const std::vector< double > &inj_c, const double dt, double *injected, double *produced, double &polyinj, double &polyprod)
 Computes injected and produced volumes of all phases, and injected and produced polymer mass. Note 1: assumes that only the first phase is injected. Note 2: assumes that transport has been done with an implicit method, i.e. that the current state gives the mobilities used for the preceding timestep. More...
 
void computeInjectedProduced (const BlackoilPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const PolymerBlackoilState &state, const std::vector< double > &transport_src, const std::vector< double > &inj_c, const double dt, double *injected, double *produced, double &polyinj, double &polyprod)
 Computes injected and produced volumes of all phases, and injected and produced polymer mass - in the compressible case. Note 1: assumes that only the first phase is injected. Note 2: assumes that transport has been done with an implicit method, i.e. that the current state gives the mobilities used for the preceding timestep. More...
 
double computePolymerMass (const std::vector< double > &pv, const std::vector< double > &s, const std::vector< double > &c, const double dps)
 Computes total (free) polymer mass over all grid cells. More...
 
double computePolymerAdsorbed (const IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< double > &pv, const std::vector< double > &cmax)
 Computes total absorbed polymer mass over all grid cells. More...
 
double computePolymerAdsorbed (const UnstructuredGrid &grid, const BlackoilPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const PolymerBlackoilState &state, const RockCompressibility *rock_comp)
 Computes total absorbed polymer mass over all grid cells. With compressibility. More...
 

Function Documentation

◆ computeFractionalFlow() [1/2]

void Opm::computeFractionalFlow ( const Opm::BlackoilPropertiesInterface &  props,
const Opm::PolymerProperties polyprops,
const std::vector< int > &  cells,
const std::vector< double > &  p,
const std::vector< double > &  T,
const std::vector< double > &  z,
const std::vector< double > &  s,
const std::vector< double > &  c,
const std::vector< double > &  cmax,
std::vector< double > &  fractional_flows 
)

Computes the fractional flow for each cell in the cells argument

Parameters
[in]propsrock and fluid properties
[in]polypropspolymer properties
[in]cellscells with which the saturation values are associated
[in]ppressure (one value per cell)
[in]Ttemperature (one value per cell)
[in]zsurface-volume values (for all P phases)
[in]ssaturation values (for all phases)
[in]cconcentration values
[in]cmaxmax polymer concentration experienced by cell
[out]fractional_flowthe fractional flow for each phase for each cell.

◆ computeFractionalFlow() [2/2]

void Opm::computeFractionalFlow ( const Opm::IncompPropertiesInterface &  props,
const Opm::PolymerProperties polyprops,
const std::vector< int > &  cells,
const std::vector< double > &  s,
const std::vector< double > &  c,
const std::vector< double > &  cmax,
std::vector< double > &  fractional_flows 
)

Computes the fractional flow for each cell in the cells argument

Parameters
[in]propsrock and fluid properties
[in]polypropspolymer properties
[in]cellscells with which the saturation values are associated
[in]ssaturation values (for all phases)
[in]cconcentration values
[in]cmaxmax polymer concentration experienced by cell
[out]fractional_flowthe fractional flow for each phase for each cell.

◆ computeInjectedProduced() [1/2]

void Opm::computeInjectedProduced ( const BlackoilPropertiesInterface &  props,
const Opm::PolymerProperties polyprops,
const PolymerBlackoilState state,
const std::vector< double > &  transport_src,
const std::vector< double > &  inj_c,
const double  dt,
double *  injected,
double *  produced,
double &  polyinj,
double &  polyprod 
)

Computes injected and produced volumes of all phases, and injected and produced polymer mass - in the compressible case. Note 1: assumes that only the first phase is injected. Note 2: assumes that transport has been done with an implicit method, i.e. that the current state gives the mobilities used for the preceding timestep.

Parameters
[in]propsfluid and rock properties.
[in]polypropspolymer properties
[in]statestate variables (pressure, fluxes etc.)
[in]srcif < 0: total reservoir volume outflow, if > 0: first phase surface volume inflow.
[in]inj_cinjected concentration by cell
[in]dttimestep used
[out]injectedmust point to a valid array with P elements, where P = s.size()/src.size().
[out]producedmust also point to a valid array with P elements.
[out]polyinjinjected mass of polymer
[out]polyprodproduced mass of polymer

◆ computeInjectedProduced() [2/2]

void Opm::computeInjectedProduced ( const IncompPropertiesInterface &  props,
const Opm::PolymerProperties polyprops,
const PolymerState state,
const std::vector< double > &  transport_src,
const std::vector< double > &  inj_c,
const double  dt,
double *  injected,
double *  produced,
double &  polyinj,
double &  polyprod 
)

Computes injected and produced volumes of all phases, and injected and produced polymer mass. Note 1: assumes that only the first phase is injected. Note 2: assumes that transport has been done with an implicit method, i.e. that the current state gives the mobilities used for the preceding timestep.

Parameters
[in]propsfluid and rock properties.
[in]polypropspolymer properties
[in]statestate variables (pressure, fluxes etc.)
[in]srcif < 0: total reservoir volume outflow, if > 0: first phase reservoir volume inflow.
[in]inj_cinjected concentration by cell
[in]dttimestep used
[out]injectedmust point to a valid array with P elements, where P = s.size()/src.size().
[out]producedmust also point to a valid array with P elements.
[out]polyinjinjected mass of polymer
[out]polyprodproduced mass of polymer

◆ computePolymerAdsorbed() [1/2]

double Opm::computePolymerAdsorbed ( const IncompPropertiesInterface &  props,
const Opm::PolymerProperties polyprops,
const std::vector< double > &  pv,
const std::vector< double > &  cmax 
)

Computes total absorbed polymer mass over all grid cells.

Parameters
[in]propsfluid and rock properties.
[in]polypropspolymer properties
[in]pvthe pore volume by cell.
[in]cmaxmax polymer concentration for cell
Returns
total absorbed polymer mass.

◆ computePolymerAdsorbed() [2/2]

double Opm::computePolymerAdsorbed ( const UnstructuredGrid &  grid,
const BlackoilPropertiesInterface &  props,
const Opm::PolymerProperties polyprops,
const PolymerBlackoilState state,
const RockCompressibility *  rock_comp 
)

Computes total absorbed polymer mass over all grid cells. With compressibility.

Parameters
[in]gridgrid
[in]propsfluid and rock properties.
[in]polypropspolymer properties
[in]stateState variables
[in]rock_compRock compressibility (optional)
Returns
total absorbed polymer mass.

◆ computePolymerMass()

double Opm::computePolymerMass ( const std::vector< double > &  pv,
const std::vector< double > &  s,
const std::vector< double > &  c,
const double  dps 
)

Computes total (free) polymer mass over all grid cells.

Parameters
[in]pvthe pore volume by cell.
[in]ssaturation values (for all P phases)
[in]cpolymer concentration
[in]dpsdead pore space
Returns
total polymer mass in grid.

◆ computeTotalMobility()

void Opm::computeTotalMobility ( const Opm::IncompPropertiesInterface &  props,
const Opm::PolymerProperties polyprops,
const std::vector< int > &  cells,
const std::vector< double > &  s,
const std::vector< double > &  c,
const std::vector< double > &  cmax,
std::vector< double > &  totmob 
)

Computes total mobility for a set of s/c values.

Parameters
[in]propsrock and fluid properties
[in]polypropspolymer properties
[in]cellscells with which the saturation values are associated
[in]ssaturation values (for all phases)
[in]cpolymer concentration
[out]totmobtotal mobilities.

◆ computeTotalMobilityOmega()

void Opm::computeTotalMobilityOmega ( const Opm::IncompPropertiesInterface &  props,
const Opm::PolymerProperties polyprops,
const std::vector< int > &  cells,
const std::vector< double > &  s,
const std::vector< double > &  c,
const std::vector< double > &  cmax,
std::vector< double > &  totmob,
std::vector< double > &  omega 
)

Computes total mobility and omega for a set of s/c values.

Parameters
[in]propsrock and fluid properties
[in]polypropspolymer properties
[in]cellscells with which the saturation values are associated
[in]ssaturation values (for all phases)
[in]cpolymer concentration
[in]cmaxmax polymer concentration experienced by cell
[out]totmobtotal mobility
[out]omegamobility-weighted (or fractional-flow weighted) fluid densities.