C++ API#

Basic C++ API of the Problem Class and Cell-Centered Uniform Mesh Class.

Problem Class#

A pressio-demoapps C++ problem class meets the following API

class Problem#
type scalar_type#

The scalar type: this is by default double.

type independent_variable_type#

This represents “time” and is the same as the scalar_type.

type state_type#

Data structure type storing the state: this depends on which specific implementation you use when you instantiate the problem using the create_problem_<> API (see Supported Backends And Types for more details about backends and corresponding types). For example, if you use create_problem_eigen, the state_type will be an Eigen vector.

type right_hand_side_type#

Data structure type to store the RHS: this depends on which specific implementation you use when you instantiate the problem using the create_problem_<> API (see Supported Backends And Types for more details about backends and corresponding types). For example, if you use create_problem_eigen, the right_hand_side_type will be an Eigen vector.

type rhs_type#

Shortcut for right_hand_side_type

type jacobian_type#

Data structure type to store the Jacobian: this depends on which specific implementation you use when you instantiate the problem using the create_problem_<> API (see Supported Backends And Types for more details about backends and corresponding types). For example, if you use create_problem_eigen, the jacobian_type will be an Eigen CRS matrix.

state_type initialCondition()#

Constructs and returns an instance of the initial condition for the target problem.

right_hand_side_type createRightHandSide()#

Constructs and returns an instance of the right hand side.

right_hand_side_type createRhs()#

Shortcut for createRightHandSide().

jacobian_type createJacobian()#

Constructs and returns an instance of the jacobian.

void operator()(const state_type &y, scalar_type time, right_hand_side_type &v)#

Given a state \(y\) and time \(time\), evaluates the RHS of the system overwriting \(v\).

void operator()(const state_type &y, scalar_type time, right_hand_side_type &v, jacobian_type &J, bool computeJac)#

Given a state \(y\) and time \(time\), evaluates the RHS of the system overwriting \(v\) and if computeJac == true, its Jacobian and stores it into \(J\).

void rhsAndJacobian(const state_type &y, scalar_type time, rhs_type &v, std::optional<jacobian_type*> J)#

Given a state \(y\) and time \(time\), evaluates the RHS of the system overwriting \(v\) and if J == true, its Jacobian and stores it into *(J.value()).

auto totalDofSampleMesh()#

Returns the total number of degrees of freedom on the sample mesh. Note that, in general, this is not the same as the number of sample mesh cells. When you have multiple dofs/cell (for example Euler equations), then the total # of dofs on sample mesh = # dofs/cell times the # of sample mesh cells.

auto totalDofStencilMesh()#

Returns the total number of degrees of freedom on the stencil mesh. Note that, in general, this is not the same as the number of stencil mesh cells. When you have multiple dofs/cell (for example Euler equations), then the total # of dofs on stencil mesh = # dofs/cell times the # of stencil mesh cells.

C++ Backends and Corresponding Types#

Backend

Type alias

Eigen

using state_type = Eigen::Matrix<scalar_type, Eigen::Dynamic, 1>

using right_hand_side_type = Eigen::Matrix<scalar_type, Eigen::Dynamic, 1>

using jacobian_type = Eigen::SparseMatrix<scalar_type, Eigen::RowMajor, int32_t>;

Cell-Centered Uniform Mesh Class#

A pressio-demoapps C++ cell-centered mesh class meets the following API

class CellCenteredUniformMesh#
type scalar_type#

This is by default double.

type index_type#

This is the type used for all integers for indexing (e.g. cell IDs, etc) so basically for all ordinals. Defaults to int32_t, which is big enough for most current cases.

int dimensionality() const#

Returns the dimensionality of this mesh: returns 1 for 1d problem, 2 for 2d, etc.

int stencilSize() const#

Returns the size of the stencil (connectivity) of this mesh object.

index_type stencilMeshSize() const#

Returns the number of stencil cells in the mesh. This corresponds to all cells where the state is defined.

index_type sampleMeshSize() const#

Returns the number of sample cells in the mesh. This corresponds to all cells where the RHS is defined.

scalar_type dx() const#

Returns the cell size along the x axis.

scalar_type dy() const#

Returns the cell size along the y axis. This is applicable only if the dimensionality is >= 2.

scalar_type dz() const#

Returns the cell size along the z axis. This is applicable only if the dimensionality is == 3.

auto viewX() const#

Returns a reference to the vector of x-coordinates of all stencil mesh cells. The type of container returned by reference depends on the backend used.

auto viewY() const#

Returns a reference to the vector of y-coordinates of all stencil mesh cells The type of container returned by reference depends on the backend used.

auto viewZ() const#

Returns a reference to the vector of z-coordinates of all stencil mesh cells The type of container returned by reference depends on the backend used.