#include <BlackoilPolymerModel.hpp>

Inheritance diagram for Opm::BlackoilPolymerModel< Grid >:
Inheritance graph

Public Types

typedef BlackoilModelBase
< Grid, BlackoilPolymerModel
< Grid > > 
Base
 
typedef Base::ReservoirState ReservoirState
 
typedef Base::WellState WellState
 

Public Member Functions

 BlackoilPolymerModel (const typename Base::ModelParameters &param, const Grid &grid, const BlackoilPropsAdInterface &fluid, const DerivedGeology &geo, const RockCompressibility *rock_comp_props, const PolymerPropsAd &polymer_props_ad, const Wells *wells, const NewtonIterationBlackoilInterface &linsolver, EclipseStateConstPtr eclipse_state, const bool has_disgas, const bool has_vapoil, const bool has_polymer, const bool has_plyshlog, const bool has_shrate, const std::vector< double > &wells_rep_radius, const std::vector< double > &wells_perf_length, const std::vector< double > &wells_bore_diameter, const bool terminal_output)
 
void prepareStep (const double dt, ReservoirState &reservoir_state, WellState &well_state)
 
void afterStep (const double dt, ReservoirState &reservoir_state, WellState &well_state)
 
void updateState (const V &dx, ReservoirState &reservoir_state, WellState &well_state)
 
void assemble (const ReservoirState &reservoir_state, WellState &well_state, const bool initial_assembly)
 

Protected Types

enum  { Concentration = CanonicalVariablePositions::Next }
 
typedef Base::SolutionState SolutionState
 
typedef Base::DataBlock DataBlock
 

Protected Member Functions

void makeConstantState (SolutionState &state) const
 
std::vector< V > variableStateInitials (const ReservoirState &x, const WellState &xw) const
 
std::vector< int > variableStateIndices () const
 
SolutionState variableStateExtractVars (const ReservoirState &x, const std::vector< int > &indices, std::vector< ADB > &vars) const
 
void computeAccum (const SolutionState &state, const int aix)
 
void assembleMassBalanceEq (const SolutionState &state)
 
void addWellContributionToMassBalanceEq (const std::vector< ADB > &cq_s, const SolutionState &state, WellState &xw)
 
void computeMassFlux (const int actph, const V &transi, const ADB &kr, const ADB &p, const SolutionState &state)
 
void computeCmax (ReservoirState &state)
 
ADB computeMc (const SolutionState &state) const
 
const std::vector< PhasePresence > phaseCondition () const
 
void computeWaterShearVelocityFaces (const V &transi, const std::vector< ADB > &kr, const std::vector< ADB > &phasePressure, const SolutionState &state, std::vector< double > &water_vel, std::vector< double > &visc_mult)
 
void computeWaterShearVelocityWells (const SolutionState &state, WellState &xw, const ADB &cq_sw, std::vector< double > &water_vel_wells, std::vector< double > &visc_mult_wells)
 

Protected Attributes

const PolymerPropsAdpolymer_props_ad_
 
const bool has_polymer_
 
const bool has_plyshlog_
 
const bool has_shrate_
 
const int poly_pos_
 
cmax_
 
std::vector< double > wells_rep_radius_
 
std::vector< double > wells_perf_length_
 
std::vector< double > wells_bore_diameter_
 
std::vector< double > shear_mult_faces_
 
std::vector< double > shear_mult_wells_
 

Friends

class BlackoilModelBase< Grid, BlackoilPolymerModel< Grid > >
 

Detailed Description

template<class Grid>
class Opm::BlackoilPolymerModel< Grid >

A model implementation for three-phase black oil with polymer.

The simulator is capable of handling three-phase problems where gas can be dissolved in oil and vice versa, with polymer in the water phase. It uses an industry-standard TPFA discretization with per-phase upwind weighting of mobilities.

It uses automatic differentiation via the class AutoDiffBlock to simplify assembly of the jacobian matrix.

Member Typedef Documentation

template<class Grid>
typedef BlackoilModelBase<Grid, BlackoilPolymerModel<Grid> > Opm::BlackoilPolymerModel< Grid >::Base
template<class Grid>
typedef Base::DataBlock Opm::BlackoilPolymerModel< Grid >::DataBlock
protected
template<class Grid>
typedef Base::ReservoirState Opm::BlackoilPolymerModel< Grid >::ReservoirState
template<class Grid>
typedef Base::SolutionState Opm::BlackoilPolymerModel< Grid >::SolutionState
protected
template<class Grid>
typedef Base::WellState Opm::BlackoilPolymerModel< Grid >::WellState

Member Enumeration Documentation

template<class Grid>
anonymous enum
protected
Enumerator
Concentration 

Constructor & Destructor Documentation

template<class Grid>
Opm::BlackoilPolymerModel< Grid >::BlackoilPolymerModel ( const typename Base::ModelParameters &  param,
const Grid &  grid,
const BlackoilPropsAdInterface &  fluid,
const DerivedGeology &  geo,
const RockCompressibility *  rock_comp_props,
const PolymerPropsAd polymer_props_ad,
const Wells *  wells,
const NewtonIterationBlackoilInterface &  linsolver,
EclipseStateConstPtr  eclipse_state,
const bool  has_disgas,
const bool  has_vapoil,
const bool  has_polymer,
const bool  has_plyshlog,
const bool  has_shrate,
const std::vector< double > &  wells_rep_radius,
const std::vector< double > &  wells_perf_length,
const std::vector< double > &  wells_bore_diameter,
const bool  terminal_output 
)

Construct the model. It will retain references to the arguments of this functions, and they are expected to remain in scope for the lifetime of the solver.

Parameters
[in]paramparameters
[in]gridgrid data structure
[in]fluidfluid properties
[in]georock properties
[in]rock_comp_propsif non-null, rock compressibility properties
[in]wellswell structure
[in]linsolverlinear solver
[in]has_disgasturn on dissolved gas
[in]has_vapoilturn on vaporized oil feature
[in]has_polymerturn on polymer feature
[in]has_plyshlogtrue when PLYSHLOG keyword available
[in]has_shratetrue when PLYSHLOG keyword available
[in]wells_rep_radiusrepresentative radius of well perforations during shear effects calculation
[in]wells_perf_lengthperforation length for well perforations
[in]wells_bore_diameterwellbore diameters for well performations
[in]terminal_outputrequest output to cout/cerr

References Opm::BlackoilPolymerModel< Grid >::has_polymer_, and Opm::BlackoilPolymerModel< Grid >::poly_pos_.

Member Function Documentation

template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::addWellContributionToMassBalanceEq ( const std::vector< ADB > &  cq_s,
const SolutionState state,
WellState xw 
)
protected
template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::afterStep ( const double  dt,
ReservoirState reservoir_state,
WellState well_state 
)

Called once after each time step.

Parameters
[in]dttime step size
[in,out]reservoir_statereservoir state variables
[in,out]well_statewell state variables
template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::assemble ( const ReservoirState reservoir_state,
WellState well_state,
const bool  initial_assembly 
)

Assemble the residual and Jacobian of the nonlinear system.

Parameters
[in]reservoir_statereservoir state variables
[in,out]well_statewell state variables
[in]initial_assemblypass true if this is the first call to assemble() in this timestep
template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::assembleMassBalanceEq ( const SolutionState state)
protected
template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::computeAccum ( const SolutionState state,
const int  aix 
)
protected

References cmax, and sat.

template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::computeCmax ( ReservoirState state)
protected
template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::computeMassFlux ( const int  actph,
const V &  transi,
const ADB &  kr,
const ADB &  p,
const SolutionState state 
)
protected

References cmax.

template<class Grid >
ADB Opm::BlackoilPolymerModel< Grid >::computeMc ( const SolutionState state) const
protected
template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::computeWaterShearVelocityFaces ( const V &  transi,
const std::vector< ADB > &  kr,
const std::vector< ADB > &  phasePressure,
const SolutionState state,
std::vector< double > &  water_vel,
std::vector< double > &  visc_mult 
)
protected

Computing the water velocity without shear-thinning for the cell faces. The water velocity will be used for shear-thinning calculation.

References cmax.

template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::computeWaterShearVelocityWells ( const SolutionState state,
WellState xw,
const ADB &  cq_sw,
std::vector< double > &  water_vel_wells,
std::vector< double > &  visc_mult_wells 
)
protected

Computing the water velocity without shear-thinning for the well perforations based on the water flux rate. The water velocity will be used for shear-thinning calculation.

template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::makeConstantState ( SolutionState state) const
protected
template<class Grid>
const std::vector<PhasePresence> Opm::BlackoilPolymerModel< Grid >::phaseCondition ( ) const
inlineprotected
template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::prepareStep ( const double  dt,
ReservoirState reservoir_state,
WellState well_state 
)

Called once before each time step.

Parameters
[in]dttime step size
[in,out]reservoir_statereservoir state variables
[in,out]well_statewell state variables
template<class Grid >
void Opm::BlackoilPolymerModel< Grid >::updateState ( const V &  dx,
ReservoirState reservoir_state,
WellState well_state 
)

Apply an update to the primary variables, chopped if appropriate.

Parameters
[in]dxupdates to apply to primary variables
[in,out]reservoir_statereservoir state variables
[in,out]well_statewell state variables
template<class Grid >
BlackoilPolymerModel< Grid >::SolutionState Opm::BlackoilPolymerModel< Grid >::variableStateExtractVars ( const ReservoirState x,
const std::vector< int > &  indices,
std::vector< ADB > &  vars 
) const
protected
template<class Grid >
std::vector< int > Opm::BlackoilPolymerModel< Grid >::variableStateIndices ( ) const
protected
template<class Grid >
std::vector< V > Opm::BlackoilPolymerModel< Grid >::variableStateInitials ( const ReservoirState x,
const WellState xw 
) const
protected

Friends And Related Function Documentation

template<class Grid>
friend class BlackoilModelBase< Grid, BlackoilPolymerModel< Grid > >
friend

Member Data Documentation

template<class Grid>
V Opm::BlackoilPolymerModel< Grid >::cmax_
protected
template<class Grid>
const bool Opm::BlackoilPolymerModel< Grid >::has_plyshlog_
protected
template<class Grid>
const bool Opm::BlackoilPolymerModel< Grid >::has_polymer_
protected
template<class Grid>
const bool Opm::BlackoilPolymerModel< Grid >::has_shrate_
protected
template<class Grid>
const int Opm::BlackoilPolymerModel< Grid >::poly_pos_
protected
template<class Grid>
const PolymerPropsAd& Opm::BlackoilPolymerModel< Grid >::polymer_props_ad_
protected
template<class Grid>
std::vector<double> Opm::BlackoilPolymerModel< Grid >::shear_mult_faces_
protected
template<class Grid>
std::vector<double> Opm::BlackoilPolymerModel< Grid >::shear_mult_wells_
protected
template<class Grid>
std::vector<double> Opm::BlackoilPolymerModel< Grid >::wells_bore_diameter_
protected
template<class Grid>
std::vector<double> Opm::BlackoilPolymerModel< Grid >::wells_perf_length_
protected
template<class Grid>
std::vector<double> Opm::BlackoilPolymerModel< Grid >::wells_rep_radius_
protected

The documentation for this class was generated from the following files: