preCICE
Loading...
Searching...
No Matches
precice::mesh::Data Class Reference

Describes a set of data values belonging to the vertices of a mesh. More...

#include <Data.hpp>

Collaboration diagram for precice::mesh::Data:
[legend]

Public Types

enum  typeName { SCALAR , VECTOR }

Public Member Functions

 Data (std::string name, DataID id, int dimension, int spatialDimensions=-1, int waveformDegree=time::Time::DEFAULT_WAVEFORM_DEGREE)
 Constructor.
 Data (std::string name, DataID id, int dimension, int spatialDimensions, int waveformDegree, std::vector< std::optional< double > > lowerBound, std::vector< std::optional< double > > upperBound)
 Constructor with data bounds.
Eigen::VectorXd & values ()
 Returns a reference to the data values.
const Eigen::VectorXd & values () const
 Returns a const reference to the data values.
const Eigen::MatrixXd & gradients () const
 Returns a const reference to the gradient data values.
const time::Samplesample () const
 Returns a const reference to the _sample.
time::SampleResult sampleAtTime (double time) const
 Samples _waveform at given time.
int getWaveformDegree () const
 get degree of _waveform.
time::Waveformwaveform ()
 Returns a reference to the waveform.
std::vector< std::optional< double > > getLowerBound () const
std::vector< std::optional< double > > getUpperBound () const
void moveToNextWindow ()
auto stamples () const
 Returns a the stamples from the waveform.
void setSampleAtTime (double time, const time::Sample &sample)
 Add sample at given time to the waveform.
void setGlobalSample (const time::Sample &sample)
 Set _sample.
void emplaceSampleAtTime (double time)
 Creates an empty sample at given time.
void emplaceSampleAtTime (double time, std::initializer_list< double > values)
 Creates a sample at given time with given values.
void emplaceSampleAtTime (double time, std::initializer_list< double > values, std::initializer_list< double > gradients)
 Creates a sample at given time with given values and gradients.
const std::string & getName () const
 Returns the name of the data set, as set in the config file.
DataID getID () const
 Returns the ID of the data set (supposed to be unique).
bool hasGradient () const
 Returns if the data contains gradient data.
bool hasSamples () const
 Returns if there are sample of this data.
void requireDataGradient ()
 Set the additional requirement of gradient data.
int getSpatialDimensions () const
 Returns the mesh dimension (i.e., number of rows) of one gradient data value .
int getDimensions () const
 Returns the dimension (i.e., number of components) of one data value (i.e number of columns of one gradient data value).
void allocateValues (int expectedCount)
 Allocates memory for the data values and corresponding gradient values.

Private Attributes

logging::Logger _log {"mesh::Data"}
time::Waveform _waveform
 Sample storage of this Data.
std::vector< std::optional< double > > _lowerBound
 Lower bound for data values. This vector with optional elements has size of the data dimension. The elements only have values if the user has configured relevant attributes.
std::vector< std::optional< double > > _upperBound
 Upper bound for data values. This vector with optional elements has size of the data dimension. The elements only have values if the user has configured relevant attributes.
std::string _name
 Name of the data set.
DataID _id
 ID of the data set (supposed to be unique).
int _dimensions
 Dimensionality of one data value.
int _spatialDimensions
 Spatial Dimension of one element -> number of rows (only 2, 3 allowed for 2D, 3D).
bool _hasGradient = false
 Whether gradient data is available or not.
time::Sample _sample

Detailed Description

Describes a set of data values belonging to the vertices of a mesh.

Definition at line 26 of file Data.hpp.

Member Enumeration Documentation

◆ typeName

Enumerator
SCALAR 
VECTOR 

Definition at line 29 of file Data.hpp.

Constructor & Destructor Documentation

◆ Data() [1/2]

precice::mesh::Data::Data ( std::string name,
DataID id,
int dimension,
int spatialDimensions = -1,
int waveformDegree = time::Time::DEFAULT_WAVEFORM_DEGREE )

Constructor.

Definition at line 12 of file Data.cpp.

◆ Data() [2/2]

precice::mesh::Data::Data ( std::string name,
DataID id,
int dimension,
int spatialDimensions,
int waveformDegree,
std::vector< std::optional< double > > lowerBound,
std::vector< std::optional< double > > upperBound )

Constructor with data bounds.

Definition at line 30 of file Data.cpp.

Member Function Documentation

◆ allocateValues()

void precice::mesh::Data::allocateValues ( int expectedCount)

Allocates memory for the data values and corresponding gradient values.

Parameters
expectedCountexpected number of values count (i.e. number of mesh vertices)

Definition at line 169 of file Data.cpp.

Here is the call graph for this function:

◆ emplaceSampleAtTime() [1/3]

void precice::mesh::Data::emplaceSampleAtTime ( double time)

Creates an empty sample at given time.

Definition at line 119 of file Data.cpp.

Here is the call graph for this function:

◆ emplaceSampleAtTime() [2/3]

void precice::mesh::Data::emplaceSampleAtTime ( double time,
std::initializer_list< double > values )

Creates a sample at given time with given values.

Definition at line 124 of file Data.cpp.

Here is the call graph for this function:

◆ emplaceSampleAtTime() [3/3]

void precice::mesh::Data::emplaceSampleAtTime ( double time,
std::initializer_list< double > values,
std::initializer_list< double > gradients )

Creates a sample at given time with given values and gradients.

Definition at line 130 of file Data.cpp.

Here is the call graph for this function:

◆ getDimensions()

int precice::mesh::Data::getDimensions ( ) const

Returns the dimension (i.e., number of components) of one data value (i.e number of columns of one gradient data value).

Definition at line 164 of file Data.cpp.

◆ getID()

DataID precice::mesh::Data::getID ( ) const

Returns the ID of the data set (supposed to be unique).

Definition at line 144 of file Data.cpp.

◆ getLowerBound()

std::vector< std::optional< double > > precice::mesh::Data::getLowerBound ( ) const

Definition at line 87 of file Data.cpp.

◆ getName()

const std::string & precice::mesh::Data::getName ( ) const

Returns the name of the data set, as set in the config file.

Definition at line 139 of file Data.cpp.

◆ getSpatialDimensions()

int precice::mesh::Data::getSpatialDimensions ( ) const

Returns the mesh dimension (i.e., number of rows) of one gradient data value .

Definition at line 207 of file Data.cpp.

◆ getUpperBound()

std::vector< std::optional< double > > precice::mesh::Data::getUpperBound ( ) const

Definition at line 92 of file Data.cpp.

◆ getWaveformDegree()

int precice::mesh::Data::getWaveformDegree ( ) const

get degree of _waveform.

Returns
int degree of _waveform

Definition at line 77 of file Data.cpp.

◆ gradients()

const Eigen::MatrixXd & precice::mesh::Data::gradients ( ) const

Returns a const reference to the gradient data values.

Definition at line 62 of file Data.cpp.

◆ hasGradient()

bool precice::mesh::Data::hasGradient ( ) const

Returns if the data contains gradient data.

Definition at line 149 of file Data.cpp.

◆ hasSamples()

bool precice::mesh::Data::hasSamples ( ) const

Returns if there are sample of this data.

Definition at line 154 of file Data.cpp.

◆ moveToNextWindow()

void precice::mesh::Data::moveToNextWindow ( )

Definition at line 97 of file Data.cpp.

Here is the call graph for this function:

◆ requireDataGradient()

void precice::mesh::Data::requireDataGradient ( )

Set the additional requirement of gradient data.

Definition at line 159 of file Data.cpp.

◆ sample()

const time::Sample & precice::mesh::Data::sample ( ) const

Returns a const reference to the _sample.

Definition at line 67 of file Data.cpp.

◆ sampleAtTime()

time::SampleResult precice::mesh::Data::sampleAtTime ( double time) const

Samples _waveform at given time.

Parameters
timeTime where the sampling happens.
Returns
Value of _waveform at time time.

Definition at line 72 of file Data.cpp.

◆ setGlobalSample()

void precice::mesh::Data::setGlobalSample ( const time::Sample & sample)

Set _sample.

Definition at line 113 of file Data.cpp.

Here is the call graph for this function:

◆ setSampleAtTime()

void precice::mesh::Data::setSampleAtTime ( double time,
const time::Sample & sample )

Add sample at given time to the waveform.

Definition at line 106 of file Data.cpp.

Here is the call graph for this function:

◆ stamples()

auto precice::mesh::Data::stamples ( ) const
inline

Returns a the stamples from the waveform.

Definition at line 109 of file Data.hpp.

◆ values() [1/2]

Eigen::VectorXd & precice::mesh::Data::values ( )

Returns a reference to the data values.

Definition at line 52 of file Data.cpp.

◆ values() [2/2]

const Eigen::VectorXd & precice::mesh::Data::values ( ) const

Returns a const reference to the data values.

Definition at line 57 of file Data.cpp.

◆ waveform()

time::Waveform & precice::mesh::Data::waveform ( )

Returns a reference to the waveform.

Definition at line 82 of file Data.cpp.

Member Data Documentation

◆ _dimensions

int precice::mesh::Data::_dimensions
private

Dimensionality of one data value.

Definition at line 176 of file Data.hpp.

◆ _hasGradient

bool precice::mesh::Data::_hasGradient = false
private

Whether gradient data is available or not.

Definition at line 182 of file Data.hpp.

◆ _id

DataID precice::mesh::Data::_id
private

ID of the data set (supposed to be unique).

Definition at line 173 of file Data.hpp.

◆ _log

logging::Logger precice::mesh::Data::_log {"mesh::Data"}
private

Definition at line 158 of file Data.hpp.

◆ _lowerBound

std::vector<std::optional<double> > precice::mesh::Data::_lowerBound
private

Lower bound for data values. This vector with optional elements has size of the data dimension. The elements only have values if the user has configured relevant attributes.

Definition at line 164 of file Data.hpp.

◆ _name

std::string precice::mesh::Data::_name
private

Name of the data set.

Definition at line 170 of file Data.hpp.

◆ _sample

time::Sample precice::mesh::Data::_sample
private

Definition at line 184 of file Data.hpp.

◆ _spatialDimensions

int precice::mesh::Data::_spatialDimensions
private

Spatial Dimension of one element -> number of rows (only 2, 3 allowed for 2D, 3D).

Definition at line 179 of file Data.hpp.

◆ _upperBound

std::vector<std::optional<double> > precice::mesh::Data::_upperBound
private

Upper bound for data values. This vector with optional elements has size of the data dimension. The elements only have values if the user has configured relevant attributes.

Definition at line 167 of file Data.hpp.

◆ _waveform

time::Waveform precice::mesh::Data::_waveform
private

Sample storage of this Data.

Definition at line 161 of file Data.hpp.


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