Skip to content

Namespace ae108::cpppetsc

Namespace List > ae108 > cpppetsc

Manage a distributed data mesh. More...

Classes

Type Name
class Context <class Policy>
struct Exception
struct GeneralizedMeshBoundaryCondition <class MeshType>
Defines a boundary condition that sets the target value (specified by a vertex view and a degree of freedom index) to the sum of scaled source values (specified by a factor, a vertex index and a degree of freedom index).
struct InsertionProxy <Mode>
This class is used to provide a "matrix"-like interface (i.e. access elements via operator()) to the Matrix wrapper. PETSc requires that calls to INSERT and ADD are not mixed; therefore specializations of this class only provide one of these methods.
struct InsertionProxy< ADD_VALUES > <>
This class redirects calls to += to the internal functor.
struct InsertionProxy< INSERT_VALUES > <>
This class redirects calls to = to the internal functor.
struct InvalidParametersException
class LeastSquaresSolver <class Policy>
class LinearSolver <class Policy>
This class can be used to solve linear systems Ax = b.
struct LinearSolverDivergedException
class LocalElementView <class MeshType_>
class LocalVertexView <class MeshType_>
class Matrix <class Policy>
class Mesh <class Policy>
struct MeshBoundaryCondition <class MeshType>
class MeshDataProvider <class MeshType_>
class NonlinearSolver <class Policy>
struct NonlinearSolverDivergedException
struct ParallelComputePolicy
The policies configure the wrapper classes (error handling, communicator).
struct SequentialComputePolicy
The policies configure the wrapper classes (error handling, communicator).
class SharedEntity <typename T, typename Policy>
class TAOSolver <class Policy>
struct TAOSolverDivergedException
class TaggedEntity <class Entity, class Tag>
A wrapper that adds a tag to a given class.
class Vector <class Policy>
class Viewer <class Policy>

Public Types

Type Name
enum std::uint8_t OwnershipType
typedef std::unique_ptr< std::remove_pointer_t< T >, std::function< void(T)> > UniqueEntity
typedef TaggedEntity< Entity, DistributedTag > distributed
A tag to annotate a "distributed" vector. This vector contains the data for all vertices, but the data is distributed between MPI nodes.
typedef TaggedEntity< Entity, GlobalTag > global
A tag to annotate a "global" vector. This vector contains the data for all vertices, and all data is available locally.
typedef TaggedEntity< Entity, LocalTag > local
A tag to annotate a "local" vector. This vector only contains the data for all local vertices.

Public Functions

Type Name
Matrix< Policy > asHermitianTransposedMatrix (const Matrix< Policy > * matrix)
Returns a matrix that behaves like the Hermitian transpose.
template Matrix< SequentialComputePolicy > asHermitianTransposedMatrix (const Matrix< SequentialComputePolicy > * matrix)
template Matrix< ParallelComputePolicy > asHermitianTransposedMatrix (const Matrix< ParallelComputePolicy > * matrix)
Matrix< Policy > asInverseMatrix (const Matrix< Policy > * matrix)
Returns a matrix that behaves like the inverse of the provided matrix.
template Matrix< SequentialComputePolicy > asInverseMatrix (const Matrix< SequentialComputePolicy > *)
template Matrix< ParallelComputePolicy > asInverseMatrix (const Matrix< ParallelComputePolicy > *)
Matrix< Policy > asSchurComplement (const Matrix< Policy > * matrix_00, const Matrix< Policy > * matrix_01, const Matrix< Policy > * matrix_10, const Matrix< Policy > * matrix_11)
Returns a matrix that behaves like the Schur complement: M_11 - M_10 * M_00^-1 * M_01.
template Matrix< SequentialComputePolicy > asSchurComplement (const Matrix< SequentialComputePolicy > * matrix_00, const Matrix< SequentialComputePolicy > * matrix_01, const Matrix< SequentialComputePolicy > * matrix_10, const Matrix< SequentialComputePolicy > * matrix_11)
template Matrix< ParallelComputePolicy > asSchurComplement (const Matrix< ParallelComputePolicy > * matrix_00, const Matrix< ParallelComputePolicy > * matrix_01, const Matrix< ParallelComputePolicy > * matrix_10, const Matrix< ParallelComputePolicy > * matrix_11)
Matrix< Policy > asSchurComplement (const Matrix< Policy > * matrix, const std::vector< typename Matrix< Policy >::size_type > & indices)
Returns a matrix that behaves like the Schur complement: M_11 - M_10 * M_00^-1 * M_01.
template Matrix< SequentialComputePolicy > asSchurComplement (const Matrix< SequentialComputePolicy > *, const std::vector< typename Matrix< SequentialComputePolicy >::size_type > &)
template Matrix< ParallelComputePolicy > asSchurComplement (const Matrix< ParallelComputePolicy > *, const std::vector< typename Matrix< ParallelComputePolicy >::size_type > &)
Matrix< Policy > asThAT (const Matrix< Policy > * matrix, const Matrix< Policy > * transform)
Returns a matrix that behaves like the transformed matrix T^H * A * T, where ^H denotes the Hermitian transpose.
template Matrix< SequentialComputePolicy > asThAT (const Matrix< SequentialComputePolicy > *, const Matrix< SequentialComputePolicy > *)
template Matrix< ParallelComputePolicy > asThAT (const Matrix< ParallelComputePolicy > *, const Matrix< ParallelComputePolicy > *)
Matrix< Policy > asTransposedMatrix (const Matrix< Policy > * matrix)
Returns a matrix that behaves like the transpose.
template Matrix< SequentialComputePolicy > asTransposedMatrix (const Matrix< SequentialComputePolicy > * matrix)
template Matrix< ParallelComputePolicy > asTransposedMatrix (const Matrix< ParallelComputePolicy > * matrix)
TaggedEntity< Vector< Policy >, Tag > clone (const TaggedEntity< Vector< Policy >, Tag > & vector)
Creates a clone of the vector.
Matrix< Policy > clone (const Matrix< Policy > & matrix)
Creates a clone of the matrix.
template TaggedEntity< Vector< SequentialComputePolicy >, LocalTag > clone (const TaggedEntity< Vector< SequentialComputePolicy >, LocalTag > &)
template TaggedEntity< Vector< SequentialComputePolicy >, DistributedTag > clone (const TaggedEntity< Vector< SequentialComputePolicy >, DistributedTag > &)
template TaggedEntity< Vector< SequentialComputePolicy >, GlobalTag > clone (const TaggedEntity< Vector< SequentialComputePolicy >, GlobalTag > &)
template TaggedEntity< Vector< ParallelComputePolicy >, LocalTag > clone (const TaggedEntity< Vector< ParallelComputePolicy >, LocalTag > &)
template TaggedEntity< Vector< ParallelComputePolicy >, DistributedTag > clone (const TaggedEntity< Vector< ParallelComputePolicy >, DistributedTag > &)
template TaggedEntity< Vector< ParallelComputePolicy >, GlobalTag > clone (const TaggedEntity< Vector< ParallelComputePolicy >, GlobalTag > &)
template Matrix< SequentialComputePolicy > clone (const Matrix< SequentialComputePolicy > & matrix)
template Matrix< ParallelComputePolicy > clone (const Matrix< ParallelComputePolicy > & matrix)
Matrix< Policy > computeElementsOfMatrix (const Matrix< Policy > & matrix)
Computes the elements of the given matrix by evaluating matrix-vector products.
template Matrix< SequentialComputePolicy > computeElementsOfMatrix (const Matrix< SequentialComputePolicy > &)
template Matrix< ParallelComputePolicy > computeElementsOfMatrix (const Matrix< ParallelComputePolicy > &)
cpppetsc::Matrix< Policy > computeElementsOfMatrix (const cpppetsc::Matrix< Policy > & matrix)
Computes the elements of the given matrix by evaluating matrix-vector products.
Matrix< Policy > computeSchurComplement (const Matrix< Policy > & matrix)
Explicitly computes the elements of the Schur complement represented by matrix .
template Matrix< SequentialComputePolicy > computeSchurComplement (const Matrix< SequentialComputePolicy > & matrix)
template Matrix< ParallelComputePolicy > computeSchurComplement (const Matrix< ParallelComputePolicy > & matrix)
void copy (const TaggedEntity< Vector< Policy >, Tag > & from, TaggedEntity< Vector< Policy >, Tag > * to)
Copies the contents of vector from to vectorto . The layouts of the vectors must be identical.
template void copy (const TaggedEntity< Vector< SequentialComputePolicy >, LocalTag > &, TaggedEntity< Vector< SequentialComputePolicy >, LocalTag > *)
template void copy (const TaggedEntity< Vector< SequentialComputePolicy >, DistributedTag > &, TaggedEntity< Vector< SequentialComputePolicy >, DistributedTag > *)
template void copy (const TaggedEntity< Vector< SequentialComputePolicy >, GlobalTag > &, TaggedEntity< Vector< SequentialComputePolicy >, GlobalTag > *)
template void copy (const TaggedEntity< Vector< ParallelComputePolicy >, LocalTag > &, TaggedEntity< Vector< ParallelComputePolicy >, LocalTag > *)
template void copy (const TaggedEntity< Vector< ParallelComputePolicy >, DistributedTag > &, TaggedEntity< Vector< ParallelComputePolicy >, DistributedTag > *)
template void copy (const TaggedEntity< Vector< ParallelComputePolicy >, GlobalTag > &, TaggedEntity< Vector< ParallelComputePolicy >, GlobalTag > *)
MeshDataProvider< MeshType > createDataProviderFromMesh (const MeshType *const mesh)
Creates a data provider from the given mesh.
Matrix< Policy > createLhsTransform (const Matrix< Policy > & matrix, const typename Matrix< Policy >::size_type rows)
Creates a matrix P that can be multiplied from the left hand side with the provided Matrix A, i.e. P * A is well-formed.
Matrix< Policy > createLhsTransform (const Matrix< Policy > & matrix, const typename Matrix< Policy >::LocalRows localRows, const typename Matrix< Policy >::GlobalRows globalRows)
Creates a matrix P that can be multiplied from the left hand side with the provided Matrix A, i.e. P * A is well-formed.
template Matrix< SequentialComputePolicy > createLhsTransform (const Matrix< SequentialComputePolicy > &, typename Matrix< SequentialComputePolicy >::size_type)
template Matrix< ParallelComputePolicy > createLhsTransform (const Matrix< ParallelComputePolicy > &, typename Matrix< ParallelComputePolicy >::size_type)
template Matrix< SequentialComputePolicy > createLhsTransform (const Matrix< SequentialComputePolicy > &, const typename Matrix< SequentialComputePolicy >::LocalRows, const typename Matrix< SequentialComputePolicy >::GlobalRows)
template Matrix< ParallelComputePolicy > createLhsTransform (const Matrix< ParallelComputePolicy > &, const typename Matrix< ParallelComputePolicy >::LocalRows, const typename Matrix< ParallelComputePolicy >::GlobalRows)
Matrix< Policy > createRhsTransform (const Matrix< Policy > & matrix, const typename Matrix< Policy >::size_type columns)
Creates a matrix P that can be multiplied from the right hand side with the provided Matrix A, i.e. A * P is well-formed.
Matrix< Policy > createRhsTransform (const Matrix< Policy > & matrix, const typename Matrix< Policy >::LocalCols localColumns, const typename Matrix< Policy >::GlobalCols globalColumns)
Creates a matrix P that can be multiplied from the right hand side with the provided Matrix A, i.e. A * P is well-formed.
template Matrix< SequentialComputePolicy > createRhsTransform (const Matrix< SequentialComputePolicy > &, typename Matrix< SequentialComputePolicy >::size_type)
template Matrix< ParallelComputePolicy > createRhsTransform (const Matrix< ParallelComputePolicy > &, typename Matrix< ParallelComputePolicy >::size_type)
template Matrix< SequentialComputePolicy > createRhsTransform (const Matrix< SequentialComputePolicy > &, const typename Matrix< SequentialComputePolicy >::LocalCols, const typename Matrix< SequentialComputePolicy >::GlobalCols)
template Matrix< ParallelComputePolicy > createRhsTransform (const Matrix< ParallelComputePolicy > &, const typename Matrix< ParallelComputePolicy >::LocalCols, const typename Matrix< ParallelComputePolicy >::GlobalCols)
distributed< Vector< Policy > > createTransformInput (const Matrix< Policy > & transform)
Creates a vector that can be multiplied with the provided matrix. It is initialized with zeroes.
template distributed< Vector< SequentialComputePolicy > > createTransformInput (const Matrix< SequentialComputePolicy > &)
template distributed< Vector< ParallelComputePolicy > > createTransformInput (const Matrix< ParallelComputePolicy > &)
distributed< Vector< Policy > > createTransformOutput (const Matrix< Policy > & transform)
Creates a vector that can store the result of a multiplication with the provided matrix. It is initialized with zeroes.
template distributed< Vector< SequentialComputePolicy > > createTransformOutput (const Matrix< SequentialComputePolicy > &)
template distributed< Vector< ParallelComputePolicy > > createTransformOutput (const Matrix< ParallelComputePolicy > &)
distributed< typename Mesh< Policy >::vector_type > createVectorFromSource (const Mesh< Policy > & mesh, const typename Mesh< Policy >::size_type dofs, std::function< void(typename Mesh< Policy >::size_type, typename Mesh< Policy >::value_type *)> writeVertexData)
Creates a global vector with dofs degrees of freedom per vertex in the provided in the mesh. The vector is filled with data copied (not added) from the source.
template distributed< typename Mesh< SequentialComputePolicy >::vector_type > createVectorFromSource (const Mesh< SequentialComputePolicy > &, Mesh< SequentialComputePolicy >::size_type, std::function< void(typename Mesh< SequentialComputePolicy >::size_type, typename Mesh< SequentialComputePolicy >::value_type *)>)
template distributed< typename Mesh< ParallelComputePolicy >::vector_type > createVectorFromSource (const Mesh< ParallelComputePolicy > &, Mesh< ParallelComputePolicy >::size_type, std::function< void(typename Mesh< ParallelComputePolicy >::size_type, typename Mesh< ParallelComputePolicy >::value_type *)>)
const char * getName (const distributed< Vector< Policy > > & vector)
Returns the name of the provided vector.
template const char * getName< ParallelComputePolicy > (const distributed< Vector< ParallelComputePolicy > > &)
template const char * getName< SequentialComputePolicy > (const distributed< Vector< SequentialComputePolicy > > &)
SharedEntity< T, Policy > makeSharedEntity (T t, const OwnershipType own=OwnershipType::Take)
Create a SharedEntity .
UniqueEntity< T > makeUniqueEntity (const T & t)
Create a UniqueEntity with the default destructor.
distributed< Vector< Policy > > multiply (const Matrix< Policy > & matrix, const distributed< Vector< Policy > > & vector)
Multiplies the given matrix by the given vector and returns the result.
template distributed< Vector< SequentialComputePolicy > > multiply (const Matrix< SequentialComputePolicy > & matrix, const distributed< Vector< SequentialComputePolicy > > & vector)
template distributed< Vector< ParallelComputePolicy > > multiply (const Matrix< ParallelComputePolicy > & matrix, const distributed< Vector< ParallelComputePolicy > > & vector)
Matrix< Policy > nestMatrices (const std::vector< std::vector< Matrix< Policy > * > > & matrices)
Creates a matrix that consists of all the provided matrices combined. The matrices are not actually copied; instead the returned matrix contains a reference to the provided matrices.
template Matrix< SequentialComputePolicy > nestMatrices (const std::vector< std::vector< Matrix< SequentialComputePolicy > * > > &)
template Matrix< ParallelComputePolicy > nestMatrices (const std::vector< std::vector< Matrix< ParallelComputePolicy > * > > &)
distributed< Vector< Policy > > nestVectors (const std::vector< cpppetsc::distributed< Vector< Policy > > * > & vectors)
Creates a vector that consists of all the provided vectors stacked on top of each other. The vectors are not actually copied; instead the returned vector contains a reference to the provided vectors.
template distributed< Vector< SequentialComputePolicy > > nestVectors (const std::vector< cpppetsc::distributed< Vector< SequentialComputePolicy > > * > &)
template distributed< Vector< ParallelComputePolicy > > nestVectors (const std::vector< cpppetsc::distributed< Vector< ParallelComputePolicy > > * > &)
std::ostream & operator<< (std::ostream & stream, const Matrix< Policy > & matrix)
Writes the local rows of the matrix to the stream.
template std::ostream & operator<< (std::ostream &, const Matrix< SequentialComputePolicy > &)
template std::ostream & operator<< (std::ostream &, const Matrix< ParallelComputePolicy > &)
std::ostream & operator<< (std::ostream & stream, const Vector< Policy > & vector)
Writes the local rows of the vector to the stream.
template std::ostream & operator<< (std::ostream &, const Vector< SequentialComputePolicy > &)
template std::ostream & operator<< (std::ostream &, const Vector< ParallelComputePolicy > &)
void readFromViewer (const Viewer< Policy > & viewer, distributed< Vector< Policy > > *const vector)
template void readFromViewer (const Viewer< SequentialComputePolicy > & viewer, distributed< Vector< SequentialComputePolicy > > *const vector)
template void readFromViewer (const Viewer< ParallelComputePolicy > & viewer, distributed< Vector< ParallelComputePolicy > > *const vector)
Vector< Policy >::value_type scalarProduct (const TaggedEntity< Vector< Policy >, Tag > & x, const TaggedEntity< Vector< Policy >, Tag > & y)
Computes the scalar product of the two vectors x andy .
void setName (const char * name, distributed< Vector< Policy > > *const vector)
Sets the name of the provided vector.
template void setName< ParallelComputePolicy > (const char *, distributed< Vector< ParallelComputePolicy > > *)
template void setName< SequentialComputePolicy > (const char *, distributed< Vector< SequentialComputePolicy > > *)
constexpr TaggedEntity< Entity, Tag > tag (Entity && entity)
Wraps the given entity in a TaggedEntity with the provided tag.
std::vector< typename cpppetsc::Mesh< Policy >::size_type > vertexDataOffsets (const cpppetsc::Mesh< Policy > & mesh)
Returns a vector of size N , whereN is the number of vertices in the mesh. Elementi of this vector is the index of the first degree of freedom associated with vertexi in a distributed/global vector.
template std::vector< typename cpppetsc::Mesh< SequentialComputePolicy >::size_type > vertexDataOffsets (const cpppetsc::Mesh< SequentialComputePolicy > &)
template std::vector< typename cpppetsc::Mesh< ParallelComputePolicy >::size_type > vertexDataOffsets (const cpppetsc::Mesh< ParallelComputePolicy > &)
void writeToViewer (const Mesh< Policy > & mesh, const distributed< Vector< Policy > > & coordinates, Viewer< Policy > * viewer)
Writes a mesh to a Viewer .
void writeToViewer (const distributed< Vector< Policy > > & data, Viewer< Policy > * viewer)
Writes a vector to a Viewer .
template void writeToViewer< ParallelComputePolicy > (const Mesh< ParallelComputePolicy > &, const distributed< Vector< ParallelComputePolicy > > &, Viewer< ParallelComputePolicy > *)
template void writeToViewer< ParallelComputePolicy > (const distributed< Vector< ParallelComputePolicy > > &, Viewer< ParallelComputePolicy > *)
template void writeToViewer< SequentialComputePolicy > (const Mesh< SequentialComputePolicy > &, const distributed< Vector< SequentialComputePolicy > > &, Viewer< SequentialComputePolicy > *)
template void writeToViewer< SequentialComputePolicy > (const distributed< Vector< SequentialComputePolicy > > &, Viewer< SequentialComputePolicy > *)

Detailed Description

Library for creating a FEM mesh and managing the degrees of freedom associated with the elements and their vertices (based on PETSc’s DMPLEX). It also provides wrappers around PETSc’s solvers (KSP, SNES, and TAO).

Public Types Documentation

enum OwnershipType

enum ae108::cpppetsc::OwnershipType {
    Take,
    Share
};

typedef UniqueEntity

using ae108::cpppetsc::UniqueEntity = typedef std::unique_ptr<std::remove_pointer_t<T>, std::function<void(T)> >;

typedef distributed

A tag to annotate a "distributed" vector. This vector contains the data for all vertices, but the data is distributed between MPI nodes.

using ae108::cpppetsc::distributed = typedef TaggedEntity<Entity, DistributedTag>;


typedef global

A tag to annotate a "global" vector. This vector contains the data for all vertices, and all data is available locally.

using ae108::cpppetsc::global = typedef TaggedEntity<Entity, GlobalTag>;


typedef local

A tag to annotate a "local" vector. This vector only contains the data for all local vertices.

using ae108::cpppetsc::local = typedef TaggedEntity<Entity, LocalTag>;


Public Functions Documentation

function asHermitianTransposedMatrix

Returns a matrix that behaves like the Hermitian transpose.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::asHermitianTransposedMatrix (
    const Matrix < Policy > * matrix
) 

Note that the matrix is not actually computed, but the returned matrix stores a reference to the given matrix to compute the result of operations on demand.

Parameters:

  • matrix Valid nonzero pointer.

function asHermitianTransposedMatrix

template Matrix < SequentialComputePolicy > ae108::cpppetsc::asHermitianTransposedMatrix (
    const Matrix < SequentialComputePolicy > * matrix
) 

function asHermitianTransposedMatrix

template Matrix < ParallelComputePolicy > ae108::cpppetsc::asHermitianTransposedMatrix (
    const Matrix < ParallelComputePolicy > * matrix
) 

function asInverseMatrix

Returns a matrix that behaves like the inverse of the provided matrix.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::asInverseMatrix (
    const Matrix < Policy > * matrix
) 

Note that the inverse matrix is not actually computed, but the returned matrix stores a reference to the given matrix to compute the result of operations on demand.

Parameters:

  • matrix Valid nonzero pointer.

function asInverseMatrix

template Matrix < SequentialComputePolicy > ae108::cpppetsc::asInverseMatrix (
    const Matrix < SequentialComputePolicy > *
) 

function asInverseMatrix

template Matrix < ParallelComputePolicy > ae108::cpppetsc::asInverseMatrix (
    const Matrix < ParallelComputePolicy > *
) 

function asSchurComplement

Returns a matrix that behaves like the Schur complement: M_11 - M_10 * M_00^-1 * M_01.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::asSchurComplement (
    const Matrix < Policy > * matrix_00,
    const Matrix < Policy > * matrix_01,
    const Matrix < Policy > * matrix_10,
    const Matrix < Policy > * matrix_11
) 

Note that the matrix is not actually computed, but the returned matrix stores a reference to the given matrices to compute the result of operations on demand.

Parameters:

  • matrix_00 Valid nonzero pointer.
  • matrix_01 Valid nonzero pointer.
  • matrix_10 Valid nonzero pointer.
  • matrix_11 Valid nonzero pointer.

function asSchurComplement

template Matrix < SequentialComputePolicy > ae108::cpppetsc::asSchurComplement (
    const Matrix < SequentialComputePolicy > * matrix_00,
    const Matrix < SequentialComputePolicy > * matrix_01,
    const Matrix < SequentialComputePolicy > * matrix_10,
    const Matrix < SequentialComputePolicy > * matrix_11
) 

function asSchurComplement

template Matrix < ParallelComputePolicy > ae108::cpppetsc::asSchurComplement (
    const Matrix < ParallelComputePolicy > * matrix_00,
    const Matrix < ParallelComputePolicy > * matrix_01,
    const Matrix < ParallelComputePolicy > * matrix_10,
    const Matrix < ParallelComputePolicy > * matrix_11
) 

function asSchurComplement

Returns a matrix that behaves like the Schur complement: M_11 - M_10 * M_00^-1 * M_01.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::asSchurComplement (
    const Matrix < Policy > * matrix,
    const std::vector< typename Matrix < Policy >::size_type > & indices
) 

Note that the matrix is not actually computed, but the returned matrix stores a reference to the given matrices to compute the result of operations on demand.

Parameters:

  • matrix Valid nonzero pointer.
  • indices The row/column indices that define the matrix M_11. Only indices in the local row range of the matrix need to be provided.

function asSchurComplement

template Matrix < SequentialComputePolicy > ae108::cpppetsc::asSchurComplement (
    const Matrix < SequentialComputePolicy > *,
    const std::vector< typename Matrix < SequentialComputePolicy >::size_type > &
) 

function asSchurComplement

template Matrix < ParallelComputePolicy > ae108::cpppetsc::asSchurComplement (
    const Matrix < ParallelComputePolicy > *,
    const std::vector< typename Matrix < ParallelComputePolicy >::size_type > &
) 

function asThAT

Returns a matrix that behaves like the transformed matrix T^H * A * T, where ^H denotes the Hermitian transpose.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::asThAT (
    const Matrix < Policy > * matrix,
    const Matrix < Policy > * transform
) 

Note that the matrix is not actually computed, but the return matrix stores a reference to the given matrices to compute the result of operations on demand.

Parameters:

  • matrix Valid nonzero pointer.
  • transform Valid nonzero pointer.

function asThAT

template Matrix < SequentialComputePolicy > ae108::cpppetsc::asThAT (
    const Matrix < SequentialComputePolicy > *,
    const Matrix < SequentialComputePolicy > *
) 

function asThAT

template Matrix < ParallelComputePolicy > ae108::cpppetsc::asThAT (
    const Matrix < ParallelComputePolicy > *,
    const Matrix < ParallelComputePolicy > *
) 

function asTransposedMatrix

Returns a matrix that behaves like the transpose.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::asTransposedMatrix (
    const Matrix < Policy > * matrix
) 

Note that the matrix is not actually computed, but the returned matrix stores a reference to the given matrix to compute the result of operations on demand.

Parameters:

  • matrix Valid nonzero pointer.

function asTransposedMatrix

template Matrix < SequentialComputePolicy > ae108::cpppetsc::asTransposedMatrix (
    const Matrix < SequentialComputePolicy > * matrix
) 

function asTransposedMatrix

template Matrix < ParallelComputePolicy > ae108::cpppetsc::asTransposedMatrix (
    const Matrix < ParallelComputePolicy > * matrix
) 

function clone

Creates a clone of the vector.

template<class Policy, class Tag>
TaggedEntity < Vector < Policy >, Tag > ae108::cpppetsc::clone (
    const TaggedEntity < Vector < Policy >, Tag > & vector
) 


function clone

Creates a clone of the matrix.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::clone (
    const Matrix < Policy > & matrix
) 


function clone

template TaggedEntity < Vector < SequentialComputePolicy >, LocalTag > ae108::cpppetsc::clone (
    const TaggedEntity < Vector < SequentialComputePolicy >, LocalTag > &
) 

function clone

template TaggedEntity < Vector < SequentialComputePolicy >, DistributedTag > ae108::cpppetsc::clone (
    const TaggedEntity < Vector < SequentialComputePolicy >, DistributedTag > &
) 

function clone

template TaggedEntity < Vector < SequentialComputePolicy >, GlobalTag > ae108::cpppetsc::clone (
    const TaggedEntity < Vector < SequentialComputePolicy >, GlobalTag > &
) 

function clone

template TaggedEntity < Vector < ParallelComputePolicy >, LocalTag > ae108::cpppetsc::clone (
    const TaggedEntity < Vector < ParallelComputePolicy >, LocalTag > &
) 

function clone

template TaggedEntity < Vector < ParallelComputePolicy >, DistributedTag > ae108::cpppetsc::clone (
    const TaggedEntity < Vector < ParallelComputePolicy >, DistributedTag > &
) 

function clone

template TaggedEntity < Vector < ParallelComputePolicy >, GlobalTag > ae108::cpppetsc::clone (
    const TaggedEntity < Vector < ParallelComputePolicy >, GlobalTag > &
) 

function clone

template Matrix < SequentialComputePolicy > ae108::cpppetsc::clone (
    const Matrix < SequentialComputePolicy > & matrix
) 

function clone

template Matrix < ParallelComputePolicy > ae108::cpppetsc::clone (
    const Matrix < ParallelComputePolicy > & matrix
) 

function computeElementsOfMatrix

Computes the elements of the given matrix by evaluating matrix-vector products.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::computeElementsOfMatrix (
    const Matrix < Policy > & matrix
) 


function computeElementsOfMatrix

template Matrix < SequentialComputePolicy > ae108::cpppetsc::computeElementsOfMatrix (
    const Matrix < SequentialComputePolicy > &
) 

function computeElementsOfMatrix

template Matrix < ParallelComputePolicy > ae108::cpppetsc::computeElementsOfMatrix (
    const Matrix < ParallelComputePolicy > &
) 

function computeElementsOfMatrix

Computes the elements of the given matrix by evaluating matrix-vector products.

template<class Policy>
cpppetsc::Matrix < Policy > ae108::cpppetsc::computeElementsOfMatrix (
    const cpppetsc::Matrix < Policy > & matrix
) 


function computeSchurComplement

Explicitly computes the elements of the Schur complement represented by matrix .

template<class Policy>
Matrix < Policy > ae108::cpppetsc::computeSchurComplement (
    const Matrix < Policy > & matrix
) 

Parameters:

  • matrix A matrix created by asSchurComplement.

function computeSchurComplement

template Matrix < SequentialComputePolicy > ae108::cpppetsc::computeSchurComplement (
    const Matrix < SequentialComputePolicy > & matrix
) 

function computeSchurComplement

template Matrix < ParallelComputePolicy > ae108::cpppetsc::computeSchurComplement (
    const Matrix < ParallelComputePolicy > & matrix
) 

function copy

Copies the contents of vector from to vectorto . The layouts of the vectors must be identical.

template<class Policy, class Tag>
void ae108::cpppetsc::copy (
    const TaggedEntity < Vector < Policy >, Tag > & from,
    TaggedEntity < Vector < Policy >, Tag > * to
) 


function copy

template void ae108::cpppetsc::copy (
    const TaggedEntity < Vector < SequentialComputePolicy >, LocalTag > &,
    TaggedEntity < Vector < SequentialComputePolicy >, LocalTag > *
) 

function copy

template void ae108::cpppetsc::copy (
    const TaggedEntity < Vector < SequentialComputePolicy >, DistributedTag > &,
    TaggedEntity < Vector < SequentialComputePolicy >, DistributedTag > *
) 

function copy

template void ae108::cpppetsc::copy (
    const TaggedEntity < Vector < SequentialComputePolicy >, GlobalTag > &,
    TaggedEntity < Vector < SequentialComputePolicy >, GlobalTag > *
) 

function copy

template void ae108::cpppetsc::copy (
    const TaggedEntity < Vector < ParallelComputePolicy >, LocalTag > &,
    TaggedEntity < Vector < ParallelComputePolicy >, LocalTag > *
) 

function copy

template void ae108::cpppetsc::copy (
    const TaggedEntity < Vector < ParallelComputePolicy >, DistributedTag > &,
    TaggedEntity < Vector < ParallelComputePolicy >, DistributedTag > *
) 

function copy

template void ae108::cpppetsc::copy (
    const TaggedEntity < Vector < ParallelComputePolicy >, GlobalTag > &,
    TaggedEntity < Vector < ParallelComputePolicy >, GlobalTag > *
) 

function createDataProviderFromMesh

Creates a data provider from the given mesh.

template<class MeshType>
MeshDataProvider < MeshType > ae108::cpppetsc::createDataProviderFromMesh (
    const MeshType *const mesh
) 

Parameters:

  • mesh Valid nonzero pointer.

function createLhsTransform

Creates a matrix P that can be multiplied from the left hand side with the provided Matrix A, i.e. P * A is well-formed.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::createLhsTransform (
    const Matrix < Policy > & matrix,
    const typename Matrix < Policy >::size_type rows
) 

Parameters:

  • rows Global number of rows of the created matrix.

function createLhsTransform

Creates a matrix P that can be multiplied from the left hand side with the provided Matrix A, i.e. P * A is well-formed.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::createLhsTransform (
    const Matrix < Policy > & matrix,
    const typename Matrix < Policy >::LocalRows localRows,
    const typename Matrix < Policy >::GlobalRows globalRows
) 


function createLhsTransform

template Matrix < SequentialComputePolicy > ae108::cpppetsc::createLhsTransform (
    const Matrix < SequentialComputePolicy > &,
    typename Matrix < SequentialComputePolicy >::size_type
) 

function createLhsTransform

template Matrix < ParallelComputePolicy > ae108::cpppetsc::createLhsTransform (
    const Matrix < ParallelComputePolicy > &,
    typename Matrix < ParallelComputePolicy >::size_type
) 

function createLhsTransform

template Matrix < SequentialComputePolicy > ae108::cpppetsc::createLhsTransform (
    const Matrix < SequentialComputePolicy > &,
    const typename Matrix < SequentialComputePolicy >::LocalRows,
    const typename Matrix < SequentialComputePolicy >::GlobalRows
) 

function createLhsTransform

template Matrix < ParallelComputePolicy > ae108::cpppetsc::createLhsTransform (
    const Matrix < ParallelComputePolicy > &,
    const typename Matrix < ParallelComputePolicy >::LocalRows,
    const typename Matrix < ParallelComputePolicy >::GlobalRows
) 

function createRhsTransform

Creates a matrix P that can be multiplied from the right hand side with the provided Matrix A, i.e. A * P is well-formed.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::createRhsTransform (
    const Matrix < Policy > & matrix,
    const typename Matrix < Policy >::size_type columns
) 

Parameters:

  • columns Global number of columns of the created matrix.

function createRhsTransform

Creates a matrix P that can be multiplied from the right hand side with the provided Matrix A, i.e. A * P is well-formed.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::createRhsTransform (
    const Matrix < Policy > & matrix,
    const typename Matrix < Policy >::LocalCols localColumns,
    const typename Matrix < Policy >::GlobalCols globalColumns
) 


function createRhsTransform

template Matrix < SequentialComputePolicy > ae108::cpppetsc::createRhsTransform (
    const Matrix < SequentialComputePolicy > &,
    typename Matrix < SequentialComputePolicy >::size_type
) 

function createRhsTransform

template Matrix < ParallelComputePolicy > ae108::cpppetsc::createRhsTransform (
    const Matrix < ParallelComputePolicy > &,
    typename Matrix < ParallelComputePolicy >::size_type
) 

function createRhsTransform

template Matrix < SequentialComputePolicy > ae108::cpppetsc::createRhsTransform (
    const Matrix < SequentialComputePolicy > &,
    const typename Matrix < SequentialComputePolicy >::LocalCols,
    const typename Matrix < SequentialComputePolicy >::GlobalCols
) 

function createRhsTransform

template Matrix < ParallelComputePolicy > ae108::cpppetsc::createRhsTransform (
    const Matrix < ParallelComputePolicy > &,
    const typename Matrix < ParallelComputePolicy >::LocalCols,
    const typename Matrix < ParallelComputePolicy >::GlobalCols
) 

function createTransformInput

Creates a vector that can be multiplied with the provided matrix. It is initialized with zeroes.

template<class Policy>
distributed < Vector < Policy > > ae108::cpppetsc::createTransformInput (
    const Matrix < Policy > & transform
) 


function createTransformInput

template distributed < Vector < SequentialComputePolicy > > ae108::cpppetsc::createTransformInput (
    const Matrix < SequentialComputePolicy > &
) 

function createTransformInput

template distributed < Vector < ParallelComputePolicy > > ae108::cpppetsc::createTransformInput (
    const Matrix < ParallelComputePolicy > &
) 

function createTransformOutput

Creates a vector that can store the result of a multiplication with the provided matrix. It is initialized with zeroes.

template<class Policy>
distributed < Vector < Policy > > ae108::cpppetsc::createTransformOutput (
    const Matrix < Policy > & transform
) 


function createTransformOutput

template distributed < Vector < SequentialComputePolicy > > ae108::cpppetsc::createTransformOutput (
    const Matrix < SequentialComputePolicy > &
) 

function createTransformOutput

template distributed < Vector < ParallelComputePolicy > > ae108::cpppetsc::createTransformOutput (
    const Matrix < ParallelComputePolicy > &
) 

function createVectorFromSource

Creates a global vector with dofs degrees of freedom per vertex in the provided in the mesh. The vector is filled with data copied (not added) from the source.

template<class Policy>
distributed < typename Mesh < Policy >::vector_type > ae108::cpppetsc::createVectorFromSource (
    const Mesh < Policy > & mesh,
    const typename Mesh < Policy >::size_type dofs,
    std::function< void(typename Mesh < Policy >::size_type, typename Mesh < Policy >::value_type *)> writeVertexData
) 

Parameters:

  • mesh The mesh that defines the geometry, e.g. the number of vertices.
  • dofs The desired number of degrees of freedom per vertex.
  • writeVertexData For every local vertex, the function is called with the vertex index and a buffer of size dofs. The function shall write the data for the given vertex to the buffer.

function createVectorFromSource

template distributed < typename Mesh < SequentialComputePolicy >::vector_type > ae108::cpppetsc::createVectorFromSource (
    const Mesh < SequentialComputePolicy > &,
    Mesh < SequentialComputePolicy >::size_type,
    std::function< void(typename Mesh < SequentialComputePolicy >::size_type, typename Mesh < SequentialComputePolicy >::value_type *)>
) 

function createVectorFromSource

template distributed < typename Mesh < ParallelComputePolicy >::vector_type > ae108::cpppetsc::createVectorFromSource (
    const Mesh < ParallelComputePolicy > &,
    Mesh < ParallelComputePolicy >::size_type,
    std::function< void(typename Mesh < ParallelComputePolicy >::size_type, typename Mesh < ParallelComputePolicy >::value_type *)>
) 

function getName

Returns the name of the provided vector.

template<class Policy>
const char * ae108::cpppetsc::getName (
    const distributed < Vector < Policy > > & vector
) 


function getName< ParallelComputePolicy >

template const char * ae108::cpppetsc::getName< ParallelComputePolicy > (
    const distributed < Vector < ParallelComputePolicy > > &
) 

function getName< SequentialComputePolicy >

template const char * ae108::cpppetsc::getName< SequentialComputePolicy > (
    const distributed < Vector < SequentialComputePolicy > > &
) 

function makeSharedEntity

Create a SharedEntity .

template<typename Policy, typename T>
SharedEntity < T, Policy > ae108::cpppetsc::makeSharedEntity (
    T t,
    const OwnershipType own=OwnershipType::Take
) 

Parameters:

  • own Take: do not increase reference count; Share: increase reference count

function makeUniqueEntity

Create a UniqueEntity with the default destructor.

template<typename Policy, typename T>
UniqueEntity< T > ae108::cpppetsc::makeUniqueEntity (
    const T & t
) 

Template parameters:

  • Policy The policy used to handle errors.

function multiply

Multiplies the given matrix by the given vector and returns the result.

template<class Policy>
distributed < Vector < Policy > > ae108::cpppetsc::multiply (
    const Matrix < Policy > & matrix,
    const distributed < Vector < Policy > > & vector
) 


function multiply

template distributed < Vector < SequentialComputePolicy > > ae108::cpppetsc::multiply (
    const Matrix < SequentialComputePolicy > & matrix,
    const distributed < Vector < SequentialComputePolicy > > & vector
) 

function multiply

template distributed < Vector < ParallelComputePolicy > > ae108::cpppetsc::multiply (
    const Matrix < ParallelComputePolicy > & matrix,
    const distributed < Vector < ParallelComputePolicy > > & vector
) 

function nestMatrices

Creates a matrix that consists of all the provided matrices combined. The matrices are not actually copied; instead the returned matrix contains a reference to the provided matrices.

template<class Policy>
Matrix < Policy > ae108::cpppetsc::nestMatrices (
    const std::vector< std::vector< Matrix < Policy > * > > & matrices
) 

Exception:


function nestMatrices

template Matrix < SequentialComputePolicy > ae108::cpppetsc::nestMatrices (
    const std::vector< std::vector< Matrix < SequentialComputePolicy > * > > &
) 

function nestMatrices

template Matrix < ParallelComputePolicy > ae108::cpppetsc::nestMatrices (
    const std::vector< std::vector< Matrix < ParallelComputePolicy > * > > &
) 

function nestVectors

Creates a vector that consists of all the provided vectors stacked on top of each other. The vectors are not actually copied; instead the returned vector contains a reference to the provided vectors.

template<class Policy>
distributed < Vector < Policy > > ae108::cpppetsc::nestVectors (
    const std::vector< cpppetsc::distributed < Vector < Policy > > * > & vectors
) 


function nestVectors

template distributed < Vector < SequentialComputePolicy > > ae108::cpppetsc::nestVectors (
    const std::vector< cpppetsc::distributed < Vector < SequentialComputePolicy > > * > &
) 

function nestVectors

template distributed < Vector < ParallelComputePolicy > > ae108::cpppetsc::nestVectors (
    const std::vector< cpppetsc::distributed < Vector < ParallelComputePolicy > > * > &
) 

function operator<<

Writes the local rows of the matrix to the stream.

template<class Policy>
std::ostream & ae108::cpppetsc::operator<< (
    std::ostream & stream,
    const Matrix < Policy > & matrix
) 


function operator<<

template std::ostream & ae108::cpppetsc::operator<< (
    std::ostream &,
    const Matrix < SequentialComputePolicy > &
) 

function operator<<

template std::ostream & ae108::cpppetsc::operator<< (
    std::ostream &,
    const Matrix < ParallelComputePolicy > &
) 

function operator<<

Writes the local rows of the vector to the stream.

template<class Policy>
std::ostream & ae108::cpppetsc::operator<< (
    std::ostream & stream,
    const Vector < Policy > & vector
) 


function operator<<

template std::ostream & ae108::cpppetsc::operator<< (
    std::ostream &,
    const Vector < SequentialComputePolicy > &
) 

function operator<<

template std::ostream & ae108::cpppetsc::operator<< (
    std::ostream &,
    const Vector < ParallelComputePolicy > &
) 

function readFromViewer

template<class Policy>
void ae108::cpppetsc::readFromViewer (
    const Viewer < Policy > & viewer,
    distributed < Vector < Policy > > *const vector
) 

function readFromViewer

template void ae108::cpppetsc::readFromViewer (
    const Viewer < SequentialComputePolicy > & viewer,
    distributed < Vector < SequentialComputePolicy > > *const vector
) 

function readFromViewer

template void ae108::cpppetsc::readFromViewer (
    const Viewer < ParallelComputePolicy > & viewer,
    distributed < Vector < ParallelComputePolicy > > *const vector
) 

function scalarProduct

Computes the scalar product of the two vectors x andy .

template<class Policy, class Tag>
Vector < Policy >::value_type ae108::cpppetsc::scalarProduct (
    const TaggedEntity < Vector < Policy >, Tag > & x,
    const TaggedEntity < Vector < Policy >, Tag > & y
) 


function setName

Sets the name of the provided vector.

template<class Policy>
void ae108::cpppetsc::setName (
    const char * name,
    distributed < Vector < Policy > > *const vector
) 


function setName< ParallelComputePolicy >

template void ae108::cpppetsc::setName< ParallelComputePolicy > (
    const char *,
    distributed < Vector < ParallelComputePolicy > > *
) 

function setName< SequentialComputePolicy >

template void ae108::cpppetsc::setName< SequentialComputePolicy > (
    const char *,
    distributed < Vector < SequentialComputePolicy > > *
) 

function tag

Wraps the given entity in a TaggedEntity with the provided tag.

template<class Tag, class Entity>
constexpr TaggedEntity < Entity, Tag > ae108::cpppetsc::tag (
    Entity && entity
) 


function vertexDataOffsets

Returns a vector of size N , whereN is the number of vertices in the mesh. Elementi of this vector is the index of the first degree of freedom associated with vertexi in a distributed/global vector.

template<class Policy>
std::vector< typename cpppetsc::Mesh < Policy >::size_type > ae108::cpppetsc::vertexDataOffsets (
    const cpppetsc::Mesh < Policy > & mesh
) 


function vertexDataOffsets

template std::vector< typename cpppetsc::Mesh < SequentialComputePolicy >::size_type > ae108::cpppetsc::vertexDataOffsets (
    const cpppetsc::Mesh < SequentialComputePolicy > &
) 

function vertexDataOffsets

template std::vector< typename cpppetsc::Mesh < ParallelComputePolicy >::size_type > ae108::cpppetsc::vertexDataOffsets (
    const cpppetsc::Mesh < ParallelComputePolicy > &
) 

function writeToViewer

Writes a mesh to a Viewer .

template<class Policy>
void ae108::cpppetsc::writeToViewer (
    const Mesh < Policy > & mesh,
    const distributed < Vector < Policy > > & coordinates,
    Viewer < Policy > * viewer
) 


function writeToViewer

Writes a vector to a Viewer .

template<class Policy>
void ae108::cpppetsc::writeToViewer (
    const distributed < Vector < Policy > > & data,
    Viewer < Policy > * viewer
) 


function writeToViewer< ParallelComputePolicy >

template void ae108::cpppetsc::writeToViewer< ParallelComputePolicy > (
    const Mesh < ParallelComputePolicy > &,
    const distributed < Vector < ParallelComputePolicy > > &,
    Viewer < ParallelComputePolicy > *
) 

function writeToViewer< ParallelComputePolicy >

template void ae108::cpppetsc::writeToViewer< ParallelComputePolicy > (
    const distributed < Vector < ParallelComputePolicy > > &,
    Viewer < ParallelComputePolicy > *
) 

function writeToViewer< SequentialComputePolicy >

template void ae108::cpppetsc::writeToViewer< SequentialComputePolicy > (
    const Mesh < SequentialComputePolicy > &,
    const distributed < Vector < SequentialComputePolicy > > &,
    Viewer < SequentialComputePolicy > *
) 

function writeToViewer< SequentialComputePolicy >

template void ae108::cpppetsc::writeToViewer< SequentialComputePolicy > (
    const distributed < Vector < SequentialComputePolicy > > &,
    Viewer < SequentialComputePolicy > *
) 


The documentation for this class was generated from the following file cpppetsc/src/include/ae108/cpppetsc/asHermitianTransposedMatrix.h