My Project
Loading...
Searching...
No Matches
Opm::UDQConfig Class Reference

Collection of all user-defined quantities in the current simulation run. More...

#include <UDQConfig.hpp>

Classes

class  DynamicSelector
 Container of entities from a dynamic context. More...
 

Public Types

using RegionSetMatcherFactory = std::function< std::unique_ptr< RegionSetMatcher >()>
 Factory function for constructing region set matchers.
 
using SegmentMatcherFactory = std::function< std::unique_ptr< SegmentMatcher >()>
 Factory function for constructing segment set matchers.
 

Public Member Functions

 UDQConfig ()=default
 Default constructor.
 
 UDQConfig (const UDQParams &params)
 Constructor.
 
 UDQConfig (const UDQParams &params, const RestartIO::RstState &rst_state)
 Constructor.
 
const std::string & unit (const std::string &key) const
 Retrieve unit string for a particular UDQ.
 
bool has_unit (const std::string &keyword) const
 Query whether or not a particular UDQ has an associated unit string.
 
bool has_keyword (const std::string &keyword) const
 Query whether or not a particular UDQ exists in the collection.
 
void add_record (SegmentMatcherFactory create_segment_matcher, const DeckRecord &record, const KeywordLocation &location, std::size_t report_step, const std::optional< DynamicSelector > &dynamic_selector=std::nullopt)
 Incorporate a single UDQ record into the known collection.
 
void add_unit (const std::string &keyword, const std::string &unit)
 Incorporate a unit string for a UDQ.
 
void add_update (const std::string &keyword, std::size_t report_step, const KeywordLocation &location, const std::vector< std::string > &data)
 Incorporate update status change for a UDQ.
 
void add_assign (const std::string &quantity, SegmentMatcherFactory create_segment_matcher, const std::vector< std::string > &selector, double value, std::size_t report_step, const std::optional< DynamicSelector > &dynamic_selector=std::nullopt)
 Incorporate a UDQ assignment.
 
void add_define (const std::string &quantity, const KeywordLocation &location, const std::vector< std::string > &expression, std::size_t report_step)
 Incorporate a UDQ defining expressions.
 
void add_table (const std::string &name, UDT udt)
 Incorporate a user defined table.
 
bool clear_pending_assignments ()
 Clear all pending assignments.
 
void eval_assign (const WellMatcher &wm, const GroupOrder &go, SegmentMatcherFactory create_segment_matcher, SummaryState &st, UDQState &udq_state) const
 Apply all pending assignments.
 
void eval (std::size_t report_step, const WellMatcher &wm, const GroupOrder &go, SegmentMatcherFactory create_segment_matcher, RegionSetMatcherFactory create_region_matcher, SummaryState &st, UDQState &udq_state) const
 Compute new values for all UDQs.
 
const UDQDefinedefine (const std::string &key) const
 Retrieve defining expression and evaluation object for a single UDQ.
 
const UDQAssignassign (const std::string &key) const
 Retrieve any pending assignment object for a single UDQ.
 
std::vector< UDQDefinedefinitions () const
 Retrieve defining expressions and evaluation objects for all known UDQs.
 
std::vector< UDQDefinedefinitions (UDQVarType var_type) const
 Retrieve defining expressions and evaluation objects for all known UDQs of a particular category.
 
std::vector< UDQInputinput () const
 Retrieve unprocessed input objects for all UDQs.
 
void exportTypeCount (std::array< int, static_cast< std::size_t >(UDQVarType::NumTypes)> &count) const
 Export count of all known UDQ categories in the current run.
 
std::size_t size () const
 Total number of active DEFINE and ASSIGN statements.
 
UDQInput operator[] (const std::string &keyword) const
 Unprocessed input object for named quantity.
 
UDQInput operator[] (std::size_t insert_index) const
 Unprocessed input object for enumerated quantity.
 
std::vector< UDQAssignassignments () const
 Retrieve pending assignment objects for all known UDQs.
 
std::vector< UDQAssignassignments (UDQVarType var_type) const
 Retrieve pending assignment objects for all known UDQs of a particular category.
 
const UDQParamsparams () const
 Retrieve run's active UDQ parameters.
 
const UDQFunctionTablefunction_table () const
 Retrieve run's active UDQ function table.
 
const std::unordered_map< std::string, UDT > & tables () const
 Retrieve run's active user defined tables.
 
bool operator== (const UDQConfig &config) const
 Equality predicate.
 
void required_summary (std::unordered_set< std::string > &summary_keys) const
 Export all summary vectors needed to compute values for the current collection of user defined quantities.
 
template<class Serializer >
void serializeOp (Serializer &serializer)
 Convert between byte array and object representation.
 

Static Public Member Functions

static UDQConfig serializationTestObject ()
 Create a serialisation test object.
 

Detailed Description

Collection of all user-defined quantities in the current simulation run.

Constructor & Destructor Documentation

◆ UDQConfig() [1/2]

Opm::UDQConfig::UDQConfig ( const UDQParams params)
explicit

Constructor.

Main constructor for a base run.

Parameters
[in]paramsUDQ parameters from UDQPARAM keyword.

◆ UDQConfig() [2/2]

Opm::UDQConfig::UDQConfig ( const UDQParams params,
const RestartIO::RstState rst_state 
)

Constructor.

Main constructor for a restarted simulation run.

Parameters
[in]paramsUDQ parameters from UDQPARAM keyword.
[in]rst_stateObject state from restart file information.

Member Function Documentation

◆ add_assign()

void Opm::UDQConfig::add_assign ( const std::string &  quantity,
SegmentMatcherFactory  create_segment_matcher,
const std::vector< std::string > &  selector,
double  value,
std::size_t  report_step,
const std::optional< DynamicSelector > &  dynamic_selector = std::nullopt 
)

Incorporate a UDQ assignment.

Implements the ASSIGN statement.

Parameters
[in]quantityUnqualified UDQ name such as FUNNY, GUITAR, WURST, or SUSHI.
[in]create_segment_matcherFactory function for constructing segment set matchers.
[in]selectorUDQ set element selection to which this assignment applies. Might be a well name pattern if quantity is a well level UDQ.
[in]valueNumeric value from ASSIGN record.
[in]report_stepTime at which this assignment statement is encountered.
[in]dynamic_selectorNamed entities in a dynamic context, such as the wells matching an ACTIONX condition. Nullopt if no such dynamic entities apply to this UDQ record.

◆ add_define()

void Opm::UDQConfig::add_define ( const std::string &  quantity,
const KeywordLocation location,
const std::vector< std::string > &  expression,
std::size_t  report_step 
)

Incorporate a UDQ defining expressions.

Implements the DEFINE statement.

Parameters
[in]quantityUnqualified UDQ name such as FUNNY, GUITAR, WURST, or SUSHI.
[in]locationInput file/line information for the UDQ record. Mostly for diagnostic purposes.
[in]expressionDefining expression for the UDQ. Function add_define() parses this expression into an abstract syntax tree which is used in subsequent evaluation contexts.
[in]report_stepTime at which this assignment statement is encountered.

◆ add_record()

void Opm::UDQConfig::add_record ( SegmentMatcherFactory  create_segment_matcher,
const DeckRecord record,
const KeywordLocation location,
std::size_t  report_step,
const std::optional< DynamicSelector > &  dynamic_selector = std::nullopt 
)

Incorporate a single UDQ record into the known collection.

Parameters
[in]create_segment_matcherFactory function for constructing segment set matchers.
[in]recordUDQ keyword record, such as a DEFINE, ASSIGN, UPDATE, or UNIT statement.
[in]locationInput file/line information for the UDQ record. Mostly for diagnostic purposes.
[in]report_stepTime at which this record is encountered.
[in]dynamic_selectorNamed entities in a dynamic context, such as the wells matching an ACTIONX condition. Nullopt if no such dynamic entities apply to this UDQ record.

◆ add_table()

void Opm::UDQConfig::add_table ( const std::string &  name,
UDT  udt 
)

Incorporate a user defined table.

Implements the UDT keyword.

Parameters
[in]nameName of user defined table.
[in]udtTabulated values.

◆ add_unit()

void Opm::UDQConfig::add_unit ( const std::string &  keyword,
const std::string &  unit 
)

Incorporate a unit string for a UDQ.

Implements the UNIT statement.

Parameters
[in]keywordUnqualified UDQ name such as FUNNY, GUITAR, WURST, or SUSHI.
[in]unitUnit string for UDQ keyword.

◆ add_update()

void Opm::UDQConfig::add_update ( const std::string &  keyword,
std::size_t  report_step,
const KeywordLocation location,
const std::vector< std::string > &  data 
)

Incorporate update status change for a UDQ.

Implements the UPDATE statement

Parameters
[in]keywordUnqualified UDQ name such as FUNNY, GUITAR, WURST, or SUSHI.
[in]report_stepTime at which this record is encountered.
[in]locationInput file/line information for the UDQ record. Mostly for diagnostic purposes.
[in]dataUpdate status. Should be a single element vector containing one of the status strings ON, OFF, or NEXT.

◆ assign()

const UDQAssign & Opm::UDQConfig::assign ( const std::string &  key) const

Retrieve any pending assignment object for a single UDQ.

Parameters
[in]keyUnqualified UDQ name such as FUNNY, GUITAR, WURST, or SUSHI.
Returns
Pending assignment object for key. Throws an exception of type
std::out_of_range
if no such object exists for the UDQ key.

◆ assignments()

std::vector< UDQAssign > Opm::UDQConfig::assignments ( UDQVarType  var_type) const

Retrieve pending assignment objects for all known UDQs of a particular category.

Parameters
[in]var_typeUDQ category.
Returns
Pending assignment objects for all known UDQs of category var_type.

◆ clear_pending_assignments()

bool Opm::UDQConfig::clear_pending_assignments ( )

Clear all pending assignments.

Clears all internal data structures of any assignment records. Typically called at the end of a report step in order to signify that all assignments have been applied.

Returns
Whether or not there were any active assignments in the internal representation. Allows client code to take action, if needed, based on the knowledge that all assignments have been applied and to prepare for the next report step.

◆ define()

const UDQDefine & Opm::UDQConfig::define ( const std::string &  key) const

Retrieve defining expression and evaluation object for a single UDQ.

Parameters
[in]keyUnqualified UDQ name such as FUNNY, GUITAR, WURST, or SUSHI.
Returns
Defining expression and evaluation object for key. Throws an exception of type
std::out_of_range
if no such object exists for the UDQ key.

◆ definitions()

std::vector< UDQDefine > Opm::UDQConfig::definitions ( UDQVarType  var_type) const

Retrieve defining expressions and evaluation objects for all known UDQs of a particular category.

Parameters
[in]var_typeUDQ category.
Returns
Defining expressions and evaluation objects for all known UDQs of category var_type.

◆ eval()

void Opm::UDQConfig::eval ( std::size_t  report_step,
const WellMatcher &  wm,
const GroupOrder go,
SegmentMatcherFactory  create_segment_matcher,
RegionSetMatcherFactory  create_region_matcher,
SummaryState st,
UDQState udq_state 
) const

Compute new values for all UDQs.

Uses both assignment and defining expressions as applicable. Assigns new UDQ values to both the summary and UDQ state objects.

Parameters
[in]report_stepCurrent report step.
[in]wmWell name pattern matcher.
[in]goGroup name pattern matcher.
[in]create_segment_matcherFactory function for constructing segment set matchers.
[in]create_region_matcherFactory function for constructing region set matchers.
[in,out]stSummary vectors. For output and evaluating ACTION condition purposes. Values pertaining to UDQs being assigned here will be updated.
[in,out]udq_stateDynamic values for all known UDQs. Values pertaining to UDQs being assigned here will be updated.

◆ eval_assign()

void Opm::UDQConfig::eval_assign ( const WellMatcher &  wm,
const GroupOrder go,
SegmentMatcherFactory  create_segment_matcher,
SummaryState st,
UDQState udq_state 
) const

Apply all pending assignments.

Assigns new UDQ values to both the summary and UDQ state objects.

Parameters
[in]wmWell name pattern matcher.
[in]goGroup name pattern matcher.
[in]create_segment_matcherFactory function for constructing segment set matchers.
[in,out]stSummary vectors. For output and evaluating ACTION condition purposes. Values pertaining to UDQs being assigned here will be updated.
[in,out]udq_stateDynamic values for all known UDQs. Values pertaining to UDQs being assigned here will be updated.

◆ exportTypeCount()

void Opm::UDQConfig::exportTypeCount ( std::array< int, static_cast< std::size_t >(UDQVarType::NumTypes)> &  count) const

Export count of all known UDQ categories in the current run.

Parameters
[out]countCount of all active UDQs of all categories in the current run.

◆ has_keyword()

bool Opm::UDQConfig::has_keyword ( const std::string &  keyword) const

Query whether or not a particular UDQ exists in the collection.

Parameters
[in]keywordUnqualified UDQ name such as FUNNY, GUITAR, WURST, or SUSHI.
Returns
whether or not keyword exists in the collection.

◆ has_unit()

bool Opm::UDQConfig::has_unit ( const std::string &  keyword) const

Query whether or not a particular UDQ has an associated unit string.

Parameters
[in]keywordUnqualified UDQ name such as FUNNY, GUITAR, WURST, or SUSHI.
Returns
Whether or not there is a unit string associated to keyword.

◆ input()

std::vector< UDQInput > Opm::UDQConfig::input ( ) const

Retrieve unprocessed input objects for all UDQs.

Needed for restart file output purposes.

◆ operator==()

bool Opm::UDQConfig::operator== ( const UDQConfig config) const

Equality predicate.

Parameters
[in]configObject against which
*this
will be tested for equality.
Returns
Whether or not
*this
is the same as config.

◆ operator[]() [1/2]

UDQInput Opm::UDQConfig::operator[] ( const std::string &  keyword) const

Unprocessed input object for named quantity.

Parameters
[in]keywordUnqualified UDQ name such as FUNNY, GUITAR, WURST, or SUSHI.
Returns
Unprocessed input object for keyword. Throws an exception of type
std::invalid_argument
if no such named UDQ exists.

◆ operator[]() [2/2]

UDQInput Opm::UDQConfig::operator[] ( std::size_t  insert_index) const

Unprocessed input object for enumerated quantity.

Parameters
[in]insert_indexLinear index in order of appearance for an individual UDQ.
Returns
Unprocessed input object for keyword. Throws an exception of type
std::invalid_argument
if no such numbered UDQ exists.

◆ required_summary()

void Opm::UDQConfig::required_summary ( std::unordered_set< std::string > &  summary_keys) const

Export all summary vectors needed to compute values for the current collection of user defined quantities.

Parameters
[in,out]summary_keysNamed summary vectors. Upon completion, any additional summary vectors needed to evaluate the current set of user defined quantities will be included in this set.

◆ serializeOp()

template<class Serializer >
void Opm::UDQConfig::serializeOp ( Serializer serializer)
inline

Convert between byte array and object representation.

Template Parameters
SerializerByte array conversion protocol.
Parameters
[in,out]serializerByte array conversion object.

◆ size()

std::size_t Opm::UDQConfig::size ( ) const

Total number of active DEFINE and ASSIGN statements.

Corresponds to the length of the vector returned from input().

◆ unit()

const std::string & Opm::UDQConfig::unit ( const std::string &  key) const

Retrieve unit string for a particular UDQ.

Parameters
[in]keyUnqualified UDQ name such as FUNNY, GUITAR, WURST, or SUSHI.
Returns
Input unit string associated to key. Throws an exception of type
std::invalid_argument
if no unit string exists for key.

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