preCICE v3.1.2
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Attributes | List of all members
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.
 
Eigen::VectorXd & values ()
 Returns a reference to the data values.
 
const Eigen::VectorXd & values () const
 Returns a const reference to the data values.
 
Eigen::MatrixXd & gradients ()
 Returns a reference to the gradient data values.
 
const Eigen::MatrixXd & gradients () const
 Returns a const reference to the gradient data values.
 
time::Samplesample ()
 Returns a reference to the _sample.
 
const time::Samplesample () const
 Returns a const reference to the _sample.
 
Eigen::VectorXd sampleAtTime (double time) const
 Samples _waveform at given time.
 
int getWaveformDegree () const
 get degree of _waveform.
 
time::StoragetimeStepsStorage ()
 Returns a reference to the _timeStepsStorage of _waveform.
 
void moveToNextWindow ()
 
auto stamples () const
 Returns a the stamples from _timeStepsStorage.
 
void setSampleAtTime (double time, const time::Sample &sample)
 Add sample at given time to _timeStepsStorage.
 
const std::stringgetName () 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).
 
void toZero ()
 Sets all values to zero.
 
bool hasGradient () const
 Returns if the data contains gradient 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
 Waveform wrapping this Data.
 
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 29 of file Data.hpp.

Member Enumeration Documentation

◆ typeName

Enumerator
SCALAR 
VECTOR 

Definition at line 32 of file Data.hpp.

Constructor & Destructor Documentation

◆ Data()

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.

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 121 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 116 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 93 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 88 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 159 of file Data.cpp.

◆ getWaveformDegree()

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

get degree of _waveform.

Returns
int degree of _waveform

Definition at line 63 of file Data.cpp.

Here is the call graph for this function:

◆ gradients() [1/2]

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

Returns a reference to the gradient data values.

Definition at line 38 of file Data.cpp.

◆ gradients() [2/2]

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

Returns a const reference to the gradient data values.

Definition at line 43 of file Data.cpp.

◆ hasGradient()

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

Returns if the data contains gradient data.

Definition at line 106 of file Data.cpp.

◆ moveToNextWindow()

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

Definition at line 73 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 111 of file Data.cpp.

◆ sample() [1/2]

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

Returns a reference to the _sample.

Definition at line 48 of file Data.cpp.

◆ sample() [2/2]

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

Returns a const reference to the _sample.

Definition at line 53 of file Data.cpp.

◆ sampleAtTime()

Eigen::VectorXd 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 58 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 _timeStepsStorage.

Definition at line 82 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 _timeStepsStorage.

Definition at line 100 of file Data.hpp.

Here is the call graph for this function:

◆ timeStepsStorage()

time::Storage & precice::mesh::Data::timeStepsStorage ( )

Returns a reference to the _timeStepsStorage of _waveform.

Definition at line 68 of file Data.cpp.

Here is the call graph for this function:

◆ toZero()

void precice::mesh::Data::toZero ( )

Sets all values to zero.

Definition at line 98 of file Data.cpp.

◆ values() [1/2]

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

Returns a reference to the data values.

Definition at line 28 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 33 of file Data.cpp.

Member Data Documentation

◆ _dimensions

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

Dimensionality of one data value.

Definition at line 149 of file Data.hpp.

◆ _hasGradient

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

Whether gradient data is available or not.

Definition at line 155 of file Data.hpp.

◆ _id

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

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

Definition at line 146 of file Data.hpp.

◆ _log

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

Definition at line 137 of file Data.hpp.

◆ _name

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

Name of the data set.

Definition at line 143 of file Data.hpp.

◆ _sample

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

Definition at line 157 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 152 of file Data.hpp.

◆ _waveform

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

Waveform wrapping this Data.

Definition at line 140 of file Data.hpp.


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