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

Representation of a UDQ ASSIGN statement. More...

#include <UDQAssign.hpp>

Public Types

using VString = std::vector< std::string >
 Type alias for a vector of strings. Simplifies function signatures.
 
using VEnumItems = std::vector< UDQSet::EnumeratedItems >
 Type alias for a vector of enumerated items.
 
using WGNameMatcher = std::function< VString(const VString &)>
 Call-back function type for a well/group name matcher.
 
using ItemMatcher = std::function< VEnumItems(const VString &)>
 Call-back function type for a matcher of enumerated items.
 

Public Member Functions

 UDQAssign ()=default
 Default constructor.
 
 UDQAssign (const std::string &keyword, const VString &selector, double value, std::size_t report_step)
 Constructor.
 
 UDQAssign (const std::string &keyword, const VEnumItems &selector, double value, std::size_t report_step)
 Constructor.
 
 UDQAssign (const std::string &keyword, VEnumItems &&selector, double value, std::size_t report_step)
 Constructor.
 
 UDQAssign (const std::string &keyword, const RestartIO::RstUDQ &assignRst, const std::size_t report_step)
 Constructor.
 
const std::string & keyword () const
 Name of UDQ to which this assignment applies.
 
UDQVarType var_type () const
 Kind of UDQ to which this assignment applies.
 
void add_record (const VString &selector, double value, std::size_t report_step)
 Add new record to existing UDQ assignment.
 
void add_record (const VEnumItems &selector, double value, std::size_t report_step)
 Add new record to existing UDQ assignment.
 
void add_record (VEnumItems &&selector, double value, std::size_t report_step)
 Add new record to existing UDQ assignment.
 
void add_record (const RestartIO::RstUDQ &assignRst, const std::size_t report_step)
 Add new record to existing UDQ assignment.
 
UDQSet eval (const VEnumItems &items) const
 Apply current assignment to a selection of enumerated items.
 
UDQSet eval (const VString &wells) const
 Apply current assignment to a selection of named items.
 
UDQSet eval () const
 Construct scalar UDQ set for a scalar UDQ assignment.
 
UDQSet eval (const VString &wgNames, WGNameMatcher matcher) const
 Apply current assignment to a selection of named items.
 
UDQSet eval (const VEnumItems &items, ItemMatcher matcher) const
 Apply current assignment to a selection of enumerated items.
 
std::size_t report_step () const
 Time at which this assignment happens.
 
bool operator== (const UDQAssign &data) const
 Equality predicate.
 
template<class Serializer >
void serializeOp (Serializer &serializer)
 Convert between byte array and object representation.
 

Static Public Member Functions

static UDQAssign serializationTestObject ()
 Create a serialisation test object.
 

Detailed Description

Representation of a UDQ ASSIGN statement.

Member Typedef Documentation

◆ ItemMatcher

using Opm::UDQAssign::ItemMatcher = std::function<VEnumItems(const VString&)>

Call-back function type for a matcher of enumerated items.

Takes a selector and returns a vector of such items.

◆ VEnumItems

Type alias for a vector of enumerated items.

Simplifies function signatures.

◆ WGNameMatcher

using Opm::UDQAssign::WGNameMatcher = std::function<VString(const VString&)>

Call-back function type for a well/group name matcher.

Takes a selector and returns a vector of matching well/group names.

Constructor & Destructor Documentation

◆ UDQAssign() [1/4]

Opm::UDQAssign::UDQAssign ( const std::string &  keyword,
const VString selector,
double  value,
std::size_t  report_step 
)

Constructor.

Parameters
[in]keywordUDQ name.
[in]selectorCollection of entity names to which this assignment applies. Might, for instance, be a selection of well or group names for a well/group level UDQ, or an empty vector for a scalar/field level UDQ.
[in]valueNumeric UDQ value for the entities named in the selector.
[in]report_stepTime at which this assignment happens. Assignments should be performed exactly once and the time value ensures this behaviour.

◆ UDQAssign() [2/4]

Opm::UDQAssign::UDQAssign ( const std::string &  keyword,
const VEnumItems selector,
double  value,
std::size_t  report_step 
)

Constructor.

Parameters
[in]keywordUDQ name.
[in]selectorCollection of named and numbered entities to which this assignment applies. Might, for instance, be a selection of well segments for a segment level UDQ.
[in]valueNumeric UDQ value for the entities identified in the selector.
[in]report_stepTime at which this assignment happens. Assignments should be performed exactly once and the time value ensures this behaviour.

◆ UDQAssign() [3/4]

Opm::UDQAssign::UDQAssign ( const std::string &  keyword,
VEnumItems &&  selector,
double  value,
std::size_t  report_step 
)

Constructor.

Assumes ownership over the selector.

Parameters
[in]keywordUDQ name.
[in]selectorCollection of named and numbered entities to which this assignment applies. Might, for instance, be a selection of well segments for a segment level UDQ.
[in]valueNumeric UDQ value for the entities identified in the selector.
[in]report_stepTime at which this assignment happens. Assignments should be performed exactly once and the time value ensures this behaviour.

◆ UDQAssign() [4/4]

Opm::UDQAssign::UDQAssign ( const std::string &  keyword,
const RestartIO::RstUDQ &  assignRst,
const std::size_t  report_step 
)

Constructor.

Reconstitutes an assignment from restart file information

Parameters
[in]keywordUDQ name.
[in]assignRstAggregate UDQ assignment information restored from restart file information.
[in]report_stepTime at which this assignment happens. Assignments should be performed exactly once and the time value ensures this behaviour.

Member Function Documentation

◆ add_record() [1/4]

void Opm::UDQAssign::add_record ( const RestartIO::RstUDQ &  assignRst,
const std::size_t  report_step 
)

Add new record to existing UDQ assignment.

Reconstitutes an assignment from restart file information. Mostly needed for interface compatibility in generic code.

Parameters
[in]assignRstAggregate UDQ assignment information restored from restart file information.
[in]report_stepTime at which this assignment happens. Assignments should be performed exactly once and the time value ensures this behaviour.

◆ add_record() [2/4]

void Opm::UDQAssign::add_record ( const VEnumItems selector,
double  value,
std::size_t  report_step 
)

Add new record to existing UDQ assignment.

Parameters
[in]selectorCollection of named and numbered entities to which this assignment applies. Might, for instance, be a selection of well segments for a segment level UDQ.
[in]valueNumeric UDQ value for the entities identified in the selector.
[in]report_stepTime at which this assignment happens. Assignments should be performed exactly once and the time value ensures this behaviour.

◆ add_record() [3/4]

void Opm::UDQAssign::add_record ( const VString selector,
double  value,
std::size_t  report_step 
)

Add new record to existing UDQ assignment.

Parameters
[in]selectorCollection of entity names to which this assignment applies. Might, for instance, be a selection of well or group names for a well/group level UDQ, or an empty vector for a scalar/field level UDQ.
[in]valueNumeric UDQ value for the entities named in the selector.
[in]report_stepTime at which this assignment happens. Assignments should be performed exactly once and the time value ensures this behaviour.

◆ add_record() [4/4]

void Opm::UDQAssign::add_record ( VEnumItems &&  selector,
double  value,
std::size_t  report_step 
)

Add new record to existing UDQ assignment.

Assumes ownership over the selector.

Parameters
[in]selectorCollection of named and numbered entities to which this assignment applies. Might, for instance, be a selection of well segments for a segment level UDQ.
[in]valueNumeric UDQ value for the entities identified in the selector.
[in]report_stepTime at which this assignment happens. Assignments should be performed exactly once and the time value ensures this behaviour.

◆ eval() [1/5]

UDQSet Opm::UDQAssign::eval ( ) const

Construct scalar UDQ set for a scalar UDQ assignment.

Throws an exception of type

std::invalid_argument

if this assignment statement does not pertain to a scalar or field level UDQ.

Returns
UDQ set of the current assignment of a scalar UDQ.

◆ eval() [2/5]

UDQSet Opm::UDQAssign::eval ( const VEnumItems items) const

Apply current assignment to a selection of enumerated items.

Parameters
[in]itemsCollection of named and numbered entities to which apply this assignment. Might, for instance, be a selection of well segments for a segment level UDQ.
Returns
UDQ set of the current assignment applied to items. Items known at construction time, or defined in subsequent calls to member function add_record(), will have a defined value in the resulting UDQ set while unrecognised items will have an undefined value.

◆ eval() [3/5]

UDQSet Opm::UDQAssign::eval ( const VEnumItems items,
ItemMatcher  matcher 
) const

Apply current assignment to a selection of enumerated items.

Parameters
[in]itemsCollection of named and numbered entities to which apply this assignment. Might, for instance, be a selection of well segments for a segment level UDQ. Return value will be sized according to this sequence.
[in]matcherCall-back for identifying items matching a selector.
Returns
UDQ set of the current assignment applied to items. Those items that are identified by the matcher, will have a defined value in the resulting UDQ set while unrecognised items will have an undefined value.

◆ eval() [4/5]

UDQSet Opm::UDQAssign::eval ( const VString wells) const

Apply current assignment to a selection of named items.

Parameters
[in]wellsCollection of named entities to which apply this assignment. Might, for instance, be a selection of well names for a well level UDQ.
Returns
UDQ set of the current assignment applied to wells. Named items known at construction time, or defined in subsequent calls to member function add_record(), will have a defined value in the resulting UDQ set while unrecognised items will have an undefined value.

◆ eval() [5/5]

UDQSet Opm::UDQAssign::eval ( const VString wgNames,
WGNameMatcher  matcher 
) const

Apply current assignment to a selection of named items.

Parameters
[in]wgNamesBacking sequence of wells/groups for which to create a UDQ well/group set. Return value will be sized according to this sequence.
[in]matcherCall-back for identifying wells/groups matching a selector. Might for instance be a wrapper around
WellMatcher::wells(pattern)
.
Returns
UDQ set of the current assignment applied to wgNames. Wells/groups identified by the matcher, will have a defined value in the resulting UDQ set while unrecognised wells/groups will have an undefined value.

◆ operator==()

bool Opm::UDQAssign::operator== ( const UDQAssign data) const

Equality predicate.

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

◆ report_step()

std::size_t Opm::UDQAssign::report_step ( ) const

Time at which this assignment happens.

Assignments should be performed exactly once and the time value ensures this behaviour.

Returns
Constructor argument of the same name.

◆ serializeOp()

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

Convert between byte array and object representation.

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

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