Opm::FlowDiagnostics::Toolbox Class Reference

Toolbox for running flow diagnostics. More...

#include <Toolbox.hpp>

Classes

struct  Forward
 
struct  Reverse
 

Public Member Functions

 Toolbox (const ConnectivityGraph &connectivity)
 Construct from known neighbourship relation. More...
 
 ~Toolbox ()
 Destructor. More...
 
 Toolbox (Toolbox &&rhs)
 Move constructor. More...
 
Toolboxoperator= (Toolbox &&rhs)
 Move assignment. More...
 
void assignPoreVolume (const std::vector< double > &pv)
 Assign pore volumes associated with each active cell. More...
 
void assignConnectionFlux (const ConnectionValues &flux)
 Assign fluxes associated with each connection. More...
 
void assignInflowFlux (const std::map< CellSetID, CellSetValues > &inflow_flux)
 
Forward computeInjectionDiagnostics (const std::vector< CellSet > &start_sets)
 
Reverse computeProductionDiagnostics (const std::vector< CellSet > &start_sets)
 

Detailed Description

Toolbox for running flow diagnostics.

Constructor & Destructor Documentation

◆ Toolbox() [1/2]

Opm::FlowDiagnostics::Toolbox::Toolbox ( const ConnectivityGraph connectivity)
explicit

Construct from known neighbourship relation.

◆ ~Toolbox()

Opm::FlowDiagnostics::Toolbox::~Toolbox ( )

Destructor.

◆ Toolbox() [2/2]

Opm::FlowDiagnostics::Toolbox::Toolbox ( Toolbox &&  rhs)

Move constructor.

Member Function Documentation

◆ assignConnectionFlux()

void Opm::FlowDiagnostics::Toolbox::assignConnectionFlux ( const ConnectionValues flux)

Assign fluxes associated with each connection.

Referenced by example::Setup::selectReportStep().

◆ assignInflowFlux()

void Opm::FlowDiagnostics::Toolbox::assignInflowFlux ( const std::map< CellSetID, CellSetValues > &  inflow_flux)

Assign inflow fluxes, typically from wells.

Inflow fluxes (injection) should be positive, outflow fluxes (production) should be negative, both should be given in the inflow_flux argument passed to this method. Values from a single well should typically be associated with a single CellSetID and be a single CellSetValues object.

Referenced by example::Setup::selectReportStep().

◆ assignPoreVolume()

void Opm::FlowDiagnostics::Toolbox::assignPoreVolume ( const std::vector< double > &  pv)

Assign pore volumes associated with each active cell.

Referenced by example::initToolbox().

◆ computeInjectionDiagnostics()

Forward Opm::FlowDiagnostics::Toolbox::computeInjectionDiagnostics ( const std::vector< CellSet > &  start_sets)

Compute forward time-of-flight and tracer solutions.

An element of

start_sets

provides a set of starting locations for a single tracer.

Forward time-of-flight is the time needed for a neutral fluid particle to flow from the nearest fluid source to an arbitrary point in the model. The tracer solutions identify cells that are flooded by injectors.

You must have called assignPoreVolume() and assignConnectionFlux() before calling this method.

The IDs of the

start_sets

must be unique.

◆ computeProductionDiagnostics()

Reverse Opm::FlowDiagnostics::Toolbox::computeProductionDiagnostics ( const std::vector< CellSet > &  start_sets)

Compute reverse time-of-flight and tracer solutions.

An element of

start_sets

provides a set of starting locations for a single tracer.

Reverse time-of-flight is the time needed for a neutral fluid particle to flow from an arbitrary point to the nearest fluid sink in the model. The tracer solutions identify cells that are drained by producers.

You must have called assignPoreVolume() and assignConnectionFlux() before calling this method.

The IDs of the

start_sets

must be unique.

◆ operator=()

Toolbox & Opm::FlowDiagnostics::Toolbox::operator= ( Toolbox &&  rhs)

Move assignment.


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