Opm::WellsManager Class Reference

#include <WellsManager.hpp>

Public Member Functions

 WellsManager ()
 Default constructor – no wells. More...
 
 WellsManager (struct Wells *W)
 
template<class F2C , class FC >
 WellsManager (const Opm::EclipseState &eclipseState, const Opm::Schedule &schedule, const size_t timeStep, int num_cells, const int *global_cell, const int *cart_dims, int dimensions, const F2C &f2c, FC begin_face_centroids, const DynamicListEconLimited &list_econ_limited, bool is_parallel_run=false, const std::unordered_set< std::string > &deactivated_wells=std::unordered_set< std::string >())
 
 WellsManager (const Opm::EclipseState &eclipseState, const Opm::Schedule &schedule, const size_t timeStep, const UnstructuredGrid &grid)
 
 ~WellsManager ()
 Destructor. More...
 
bool empty () const
 Does the "deck" define any wells? More...
 
const Wellsc_wells () const
 
const WellCollectionwellCollection () const
 Access the well group hierarchy. More...
 
WellCollectionwellCollection ()
 
bool conditionsMet (const std::vector< double > &well_bhp, const std::vector< double > &well_reservoirrates_phase, const std::vector< double > &well_surfacerates_phase)
 
void applyExplicitReinjectionControls (const std::vector< double > &well_reservoirrates_phase, const std::vector< double > &well_surfacerates_phase)
 
template<class C2F , class FC >
 WellsManager (const Opm::EclipseState &eclipseState, const Opm::Schedule &schedule, const size_t timeStep, int number_of_cells, const int *global_cell, const int *cart_dims, int dimensions, const C2F &cell_to_faces, FC begin_face_centroids, const DynamicListEconLimited &list_econ_limited, bool is_parallel_run, const std::unordered_set< std::string > &deactivated_wells)
 

Detailed Description

This class manages a Wells struct in the sense that it encapsulates creation and destruction of the wells data structure. The resulting Wells is available through the c_wells() method.

Constructor & Destructor Documentation

◆ WellsManager() [1/5]

Opm::WellsManager::WellsManager ( )

Default constructor – no wells.

◆ WellsManager() [2/5]

Opm::WellsManager::WellsManager ( struct Wells W)
explicit

Construct from existing wells object. WellsManager is not properly initialised in the sense that the logic to manage control switching does not exist.

Parameters
[in]WExisting wells object.

◆ WellsManager() [3/5]

template<class F2C , class FC >
Opm::WellsManager::WellsManager ( const Opm::EclipseState &  eclipseState,
const Opm::Schedule &  schedule,
const size_t  timeStep,
int  num_cells,
const int *  global_cell,
const int *  cart_dims,
int  dimensions,
const F2C &  f2c,
FC  begin_face_centroids,
const DynamicListEconLimited list_econ_limited,
bool  is_parallel_run = false,
const std::unordered_set< std::string > &  deactivated_wells = std::unordered_set< std::string >() 
)

Construct from input deck and grid. The permeability argument may be zero if the input contain well productivity indices, otherwise it must be given in order to approximate these by the Peaceman formula.

Parameters
deactivated_wellsA set of wells that should be treated like shut wells. E.g. in a a parallel run these would be the wells handeled by another process. Defaults to empty set.

◆ WellsManager() [4/5]

Opm::WellsManager::WellsManager ( const Opm::EclipseState &  eclipseState,
const Opm::Schedule &  schedule,
const size_t  timeStep,
const UnstructuredGrid &  grid 
)

◆ ~WellsManager()

Opm::WellsManager::~WellsManager ( )

Destructor.

◆ WellsManager() [5/5]

template<class C2F , class FC >
Opm::WellsManager::WellsManager ( const Opm::EclipseState &  eclipseState,
const Opm::Schedule &  schedule,
const size_t  timeStep,
int  number_of_cells,
const int *  global_cell,
const int *  cart_dims,
int  dimensions,
const C2F &  cell_to_faces,
FC  begin_face_centroids,
const DynamicListEconLimited list_econ_limited,
bool  is_parallel_run,
const std::unordered_set< std::string > &  deactivated_wells 
)

Member Function Documentation

◆ applyExplicitReinjectionControls()

void Opm::WellsManager::applyExplicitReinjectionControls ( const std::vector< double > &  well_reservoirrates_phase,
const std::vector< double > &  well_surfacerates_phase 
)

Applies explicit reinjection controls. This must be called at each timestep to be correct.

Parameters
[in]well_reservoirrates_phaseA vector containing reservoir rates by phase for each well. Is assumed to be ordered the same way as the related Wells-struct, with all phase rates of a single well adjacent in the array.
[in]well_surfacerates_phaseA vector containing surface rates by phase for each well. Is assumed to be ordered the same way as the related Wells-struct, with all phase rates of a single well adjacent in the array.

◆ c_wells()

const Wells * Opm::WellsManager::c_wells ( ) const

Access the managed Wells. The method is named similarly to c_str() in std::string, to make it clear that we are returning a C-compatible struct.

◆ conditionsMet()

bool Opm::WellsManager::conditionsMet ( const std::vector< double > &  well_bhp,
const std::vector< double > &  well_reservoirrates_phase,
const std::vector< double > &  well_surfacerates_phase 
)

Checks if each condition is met, applies well controls where needed (that is, it either changes the active control of violating wells, or shuts down wells). Only one change is applied per invocation. Typical use will be

solve_pressure();
while(!wells.conditionsMet(...)) {
solve_pressure();
}
Parameters
[in]well_bhpA vector containing the bhp for each well. Is assumed to be ordered the same way as the related Wells-struct.
[in]well_reservoirrates_phaseA vector containing reservoir rates by phase for each well. Is assumed to be ordered the same way as the related Wells-struct, with all phase rates of a single well adjacent in the array.
[in]well_surfacerates_phaseA vector containing surface rates by phase for each well. Is assumed to be ordered the same way as the related Wells-struct, with all phase rates of a single well adjacent in the array.
Returns
true if no violations were found, false otherwise (false also implies a change).

◆ empty()

bool Opm::WellsManager::empty ( ) const

Does the "deck" define any wells?

◆ wellCollection() [1/2]

WellCollection & Opm::WellsManager::wellCollection ( )

◆ wellCollection() [2/2]

const WellCollection & Opm::WellsManager::wellCollection ( ) const

Access the well group hierarchy.


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