KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
List of all members
Kratos::TrilinosSpace< TMatrixType, TVectorType > Class Template Reference

The space adapted for Trilinos vectors and matrices. More...

#include <trilinos_space.h>

Collaboration diagram for Kratos::TrilinosSpace< TMatrixType, TVectorType >:

Public Member Functions

Life Cycle
 TrilinosSpace ()
 Default constructor. More...
 
virtual ~TrilinosSpace ()
 Destructor. More...
 

Type Definitions

using DataType = double
 Definition of the data type. More...
 
using MatrixType = TMatrixType
 Definition of the matrix type. More...
 
using VectorType = TVectorType
 Definition of the vector type. More...
 
using IndexType = std::size_t
 Definition of the index type. More...
 
using SizeType = std::size_t
 Definition of the size type. More...
 
using MatrixPointerType = typename Kratos::shared_ptr< TMatrixType >
 Definition of the pointer types. More...
 
using VectorPointerType = typename Kratos::shared_ptr< TVectorType >
 
using DofUpdaterType = TrilinosDofUpdater< TrilinosSpace< TMatrixType, TVectorType > >
 Some other definitions. More...
 
using DofUpdaterPointerType = typename DofUpdater< TrilinosSpace< TMatrixType, TVectorType > >::UniquePointer
 
 KRATOS_CLASS_POINTER_DEFINITION (TrilinosSpace)
 Pointer definition of TrilinosSpace. More...
 

Operations

MatrixPointerType ReadMatrixMarket (const std::string FileName, Epetra_MpiComm &rComm)
 Read a matrix from a MatrixMarket file. More...
 
VectorPointerType ReadMatrixMarketVector (const std::string &rFileName, Epetra_MpiComm &rComm, const int N)
 Read a vector from a MatrixMarket file. More...
 
static MatrixPointerType CreateEmptyMatrixPointer ()
 This method creates an empty pointer to a matrix. More...
 
static VectorPointerType CreateEmptyVectorPointer ()
 This method creates an empty pointer to a vector. More...
 
static MatrixPointerType CreateEmptyMatrixPointer (Epetra_MpiComm &rComm)
 This method creates an empty pointer to a matrix using epetra communicator. More...
 
static VectorPointerType CreateEmptyVectorPointer (Epetra_MpiComm &rComm)
 This method creates an empty pointer to a vector using epetra communicator. More...
 
static IndexType Size (const VectorType &rV)
 Returns size of vector rV. More...
 
static IndexType Size1 (MatrixType const &rM)
 Returns number of rows of rM. More...
 
static IndexType Size2 (MatrixType const &rM)
 Returns number of columns of rM. More...
 
static void GetColumn (const unsigned int j, const MatrixType &rM, VectorType &rX)
 Returns the column of the matrix in the given position. More...
 
static void Copy (const MatrixType &rX, MatrixType &rY)
 Returns a copy of the matrix rX. More...
 
static void Copy (const VectorType &rX, VectorType &rY)
 Returns a copy of the vector rX. More...
 
static double Dot (const VectorType &rX, const VectorType &rY)
 Returns the product of two vectors. More...
 
static double Max (const VectorType &rX)
 Returns the maximum value of the vector rX. More...
 
static double Min (const VectorType &rX)
 Returns the minimum value of the vector rX. More...
 
static double TwoNorm (const VectorType &rX)
 Returns the norm of the vector rX. More...
 
static double TwoNorm (const MatrixType &rA)
 Returns the Frobenius norm of the matrix rX. More...
 
static void Mult (const MatrixType &rA, const VectorType &rX, VectorType &rY)
 Returns the multiplication of a matrix by a vector. More...
 
static void Mult (const MatrixType &rA, const MatrixType &rB, MatrixType &rC, const bool CallFillCompleteOnResult=true, const bool KeepAllHardZeros=false)
 Returns the multiplication matrix-matrix. More...
 
static void TransposeMult (const MatrixType &rA, const VectorType &rX, VectorType &rY)
 Returns the transpose multiplication of a matrix by a vector. More...
 
static void TransposeMult (const MatrixType &rA, const MatrixType &rB, MatrixType &rC, const std::pair< bool, bool > TransposeFlag={false, false}, const bool CallFillCompleteOnResult=true, const bool KeepAllHardZeros=false)
 Returns the transpose multiplication matrix-matrix. More...
 
static void BtDBProductOperation (MatrixType &rA, const MatrixType &rD, const MatrixType &rB, const bool CallFillCompleteOnResult=true, const bool KeepAllHardZeros=false, const bool EnforceInitialGraph=false)
 Calculates the product operation B'DB. More...
 
static void BDBtProductOperation (MatrixType &rA, const MatrixType &rD, const MatrixType &rB, const bool CallFillCompleteOnResult=true, const bool KeepAllHardZeros=false, const bool EnforceInitialGraph=false)
 Calculates the product operation BDB'. More...
 
static void InplaceMult (VectorType &rX, const double A)
 Returns the multiplication of a vector by a scalar. More...
 
static void Assign (VectorType &rX, const double A, const VectorType &rY)
 Returns the multiplication of a vector by a scalar. More...
 
static void UnaliasedAdd (VectorType &rX, const double A, const VectorType &rY)
 Returns the unaliased addition of a vector by a scalar times a vector. More...
 
static void ScaleAndAdd (const double A, const VectorType &rX, const double B, const VectorType &rY, VectorType &rZ)
 Returns the unaliased addition of two vectors by a scalar. More...
 
static void ScaleAndAdd (const double A, const VectorType &rX, const double B, VectorType &rY)
 Returns the unaliased addition of two vectors by a scalar. More...
 
static void SetValue (VectorType &rX, IndexType i, const double value)
 Sets a value in a vector. More...
 
static void Set (VectorType &rX, const DataType A)
 assigns a scalar to a vector More...
 
static void Resize (MatrixType &rA, const SizeType m, const SizeType n)
 Resizes a matrix. More...
 
static void Resize (VectorType &rX, const SizeType n)
 Resizes a vector. More...
 
static void Resize (VectorPointerType &pX, const SizeType n)
 Resizes a vector. More...
 
static void Clear (MatrixPointerType &pA)
 Clears a matrix. More...
 
static void Clear (VectorPointerType &pX)
 Clears a vector. More...
 
static void SetToZero (MatrixType &rA)
 Sets a matrix to zero. More...
 
static void SetToZero (VectorType &rX)
 Sets a vector to zero. More...
 
static void AssembleLHS (MatrixType &rA, const Matrix &rLHSContribution, const std::vector< std::size_t > &rEquationId)
 TODO: creating the the calculating reaction version. More...
 
static void AssembleRHS (VectorType &rb, const Vector &rRHSContribution, const std::vector< std::size_t > &rEquationId)
 TODO: creating the the calculating reaction version. More...
 
static constexpr bool IsDistributed ()
 This function returns if we are in a distributed system. More...
 
static double GetValue (const VectorType &rX, const std::size_t I)
 This function returns a value from a given vector according to a given index. More...
 
static void GatherValues (const VectorType &rX, const std::vector< int > &IndexArray, double *pValues)
 This function gathers the values of a given vector according to a given index array. More...
 
static Epetra_CrsGraph CombineMatricesGraphs (const MatrixType &rA, const MatrixType &rB)
 Generates a graph combining the graphs of two matrices. More...
 
static void CopyMatrixValues (MatrixType &rA, const MatrixType &rB)
 Copy values from one matrix to another. More...
 
static double CheckAndCorrectZeroDiagonalValues (const ProcessInfo &rProcessInfo, MatrixType &rA, VectorType &rb, const SCALING_DIAGONAL ScalingDiagonal=SCALING_DIAGONAL::NO_SCALING)
 This method checks and corrects the zero diagonal values. More...
 
static double GetScaleNorm (const ProcessInfo &rProcessInfo, const MatrixType &rA, const SCALING_DIAGONAL ScalingDiagonal=SCALING_DIAGONAL::NO_SCALING)
 This method returns the scale norm considering for scaling the diagonal. More...
 
static double GetDiagonalNorm (const MatrixType &rA)
 This method returns the diagonal norm considering for scaling the diagonal. More...
 
static double GetAveragevalueDiagonal (const MatrixType &rA)
 This method returns the diagonal max value. More...
 
static double GetMaxDiagonal (const MatrixType &rA)
 This method returns the diagonal max value. More...
 
static double GetMinDiagonal (const MatrixType &rA)
 This method returns the diagonal min value. More...
 
static constexpr bool IsDistributedSpace ()
 Check if the TrilinosSpace is distributed. More...
 

Input and output

virtual std::string Info () const
 Turn back information as a string. More...
 
virtual void PrintInfo (std::ostream &rOStream) const
 Print information about this object. More...
 
virtual void PrintData (std::ostream &rOStream) const
 Print object's data. More...
 
template<class TOtherMatrixType >
static bool WriteMatrixMarketMatrix (const char *pFileName, const TOtherMatrixType &rM, const bool Symmetric)
 Writes a matrix to a file in MatrixMarket format. More...
 
template<class VectorType >
static bool WriteMatrixMarketVector (const char *pFileName, const VectorType &rV)
 Writes a vector to a file in MatrixMarket format. More...
 
static DofUpdaterPointerType CreateDofUpdater ()
 Creates a new dof updater. More...
 

Detailed Description

template<class TMatrixType, class TVectorType>
class Kratos::TrilinosSpace< TMatrixType, TVectorType >

The space adapted for Trilinos vectors and matrices.

Author
Riccardo Rossi
Template Parameters
TMatrixTypeThe matrix type considered
TVectorTypethe vector type considered

Member Typedef Documentation

◆ DataType

template<class TMatrixType , class TVectorType >
using Kratos::TrilinosSpace< TMatrixType, TVectorType >::DataType = double

Definition of the data type.

◆ DofUpdaterPointerType

template<class TMatrixType , class TVectorType >
using Kratos::TrilinosSpace< TMatrixType, TVectorType >::DofUpdaterPointerType = typename DofUpdater<TrilinosSpace<TMatrixType,TVectorType> >::UniquePointer

◆ DofUpdaterType

template<class TMatrixType , class TVectorType >
using Kratos::TrilinosSpace< TMatrixType, TVectorType >::DofUpdaterType = TrilinosDofUpdater< TrilinosSpace<TMatrixType,TVectorType> >

Some other definitions.

◆ IndexType

template<class TMatrixType , class TVectorType >
using Kratos::TrilinosSpace< TMatrixType, TVectorType >::IndexType = std::size_t

Definition of the index type.

◆ MatrixPointerType

template<class TMatrixType , class TVectorType >
using Kratos::TrilinosSpace< TMatrixType, TVectorType >::MatrixPointerType = typename Kratos::shared_ptr<TMatrixType>

Definition of the pointer types.

◆ MatrixType

template<class TMatrixType , class TVectorType >
using Kratos::TrilinosSpace< TMatrixType, TVectorType >::MatrixType = TMatrixType

Definition of the matrix type.

◆ SizeType

template<class TMatrixType , class TVectorType >
using Kratos::TrilinosSpace< TMatrixType, TVectorType >::SizeType = std::size_t

Definition of the size type.

◆ VectorPointerType

template<class TMatrixType , class TVectorType >
using Kratos::TrilinosSpace< TMatrixType, TVectorType >::VectorPointerType = typename Kratos::shared_ptr<TVectorType>

◆ VectorType

template<class TMatrixType , class TVectorType >
using Kratos::TrilinosSpace< TMatrixType, TVectorType >::VectorType = TVectorType

Definition of the vector type.

Constructor & Destructor Documentation

◆ TrilinosSpace()

template<class TMatrixType , class TVectorType >
Kratos::TrilinosSpace< TMatrixType, TVectorType >::TrilinosSpace ( )
inline

Default constructor.

◆ ~TrilinosSpace()

template<class TMatrixType , class TVectorType >
virtual Kratos::TrilinosSpace< TMatrixType, TVectorType >::~TrilinosSpace ( )
inlinevirtual

Destructor.

Member Function Documentation

◆ AssembleLHS()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::AssembleLHS ( MatrixType rA,
const Matrix rLHSContribution,
const std::vector< std::size_t > &  rEquationId 
)
inlinestatic

TODO: creating the the calculating reaction version.

Assembles the LHS of the system

Parameters
rAThe LHS matrix
rLHSContributionThe contribution to the LHS
rEquationIdThe equation ids

◆ AssembleRHS()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::AssembleRHS ( VectorType rb,
const Vector rRHSContribution,
const std::vector< std::size_t > &  rEquationId 
)
inlinestatic

TODO: creating the the calculating reaction version.

Assembles the RHS of the system

Parameters
rbThe RHS vector
rRHSContributionThe RHS contribution
rEquationIdThe equation ids

◆ Assign()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::Assign ( VectorType rX,
const double  A,
const VectorType rY 
)
inlinestatic

Returns the multiplication of a vector by a scalar.

x = A*y Checks if a multiplication is needed and tries to do otherwise

Note
ATTENTION it is assumed no aliasing between rX and rY
Parameters
rXThe resulting vector considered
AThe scalar considered
rYThe multiplied vector considered

◆ BDBtProductOperation()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::BDBtProductOperation ( MatrixType rA,
const MatrixType rD,
const MatrixType rB,
const bool  CallFillCompleteOnResult = true,
const bool  KeepAllHardZeros = false,
const bool  EnforceInitialGraph = false 
)
inlinestatic

Calculates the product operation BDB'.

Parameters
rAThe resulting matrix
rDThe "center" matrix
rBThe matrices to be transposed
CallFillCompleteOnResultOptional argument, defaults to true. Power users may specify this argument to be false if they DON'T want this function to call C.FillComplete. (It is often useful to allow this function to call C.FillComplete, in cases where one or both of the input matrices are rectangular and it is not trivial to know which maps to use for the domain- and range-maps.)
KeepAllHardZerosOptional argument, defaults to false. If true, Multiply, keeps all entries in C corresponding to hard zeros. If false, the following happens by case: A*B^T, A^T*B^T - Does not store entries caused by hard zeros in C. A^T*B (unoptimized) - Hard zeros are always stored (this option has no effect) A*B, A^T*B (optimized) - Hard zeros in corresponding to hard zeros in A are not stored, There are certain cases involving reuse of C, where this can be useful.
Template Parameters
TEnforceInitialGraphIf the initial graph is enforced, or a new one is generated

◆ BtDBProductOperation()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::BtDBProductOperation ( MatrixType rA,
const MatrixType rD,
const MatrixType rB,
const bool  CallFillCompleteOnResult = true,
const bool  KeepAllHardZeros = false,
const bool  EnforceInitialGraph = false 
)
inlinestatic

Calculates the product operation B'DB.

Parameters
rAThe resulting matrix
rDThe "center" matrix
rBThe matrices to be transposed
CallFillCompleteOnResultOptional argument, defaults to true. Power users may specify this argument to be false if they DON'T want this function to call C.FillComplete. (It is often useful to allow this function to call C.FillComplete, in cases where one or both of the input matrices are rectangular and it is not trivial to know which maps to use for the domain- and range-maps.)
KeepAllHardZerosOptional argument, defaults to false. If true, Multiply, keeps all entries in C corresponding to hard zeros. If false, the following happens by case: A*B^T, A^T*B^T - Does not store entries caused by hard zeros in C. A^T*B (unoptimized) - Hard zeros are always stored (this option has no effect) A*B, A^T*B (optimized) - Hard zeros in corresponding to hard zeros in A are not stored, There are certain cases involving reuse of C, where this can be useful.
Template Parameters
TEnforceInitialGraphIf the initial graph is enforced, or a new one is generated

◆ CheckAndCorrectZeroDiagonalValues()

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::CheckAndCorrectZeroDiagonalValues ( const ProcessInfo rProcessInfo,
MatrixType rA,
VectorType rb,
const SCALING_DIAGONAL  ScalingDiagonal = SCALING_DIAGONAL::NO_SCALING 
)
inlinestatic

This method checks and corrects the zero diagonal values.

This method returns the scale norm considering for scaling the diagonal

Parameters
rProcessInfoThe problem process info
rAThe LHS matrix
rbThe RHS vector
ScalingDiagonalThe type of caling diagonal considered
Returns
The scale norm

◆ Clear() [1/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::Clear ( MatrixPointerType pA)
inlinestatic

Clears a matrix.

Parameters
pAThe pointer to the matrix to be cleared

◆ Clear() [2/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::Clear ( VectorPointerType pX)
inlinestatic

Clears a vector.

Parameters
pXThe pointer to the vector to be cleared

◆ CombineMatricesGraphs()

template<class TMatrixType , class TVectorType >
static Epetra_CrsGraph Kratos::TrilinosSpace< TMatrixType, TVectorType >::CombineMatricesGraphs ( const MatrixType rA,
const MatrixType rB 
)
inlinestatic

Generates a graph combining the graphs of two matrices.

Parameters
rAThe first matrix
rBThe second matrix

◆ Copy() [1/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::Copy ( const MatrixType rX,
MatrixType rY 
)
inlinestatic

Returns a copy of the matrix rX.

rY = rX

Parameters
rXThe matrix considered
rYThe copy of the matrix rX

◆ Copy() [2/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::Copy ( const VectorType rX,
VectorType rY 
)
inlinestatic

Returns a copy of the vector rX.

rY = rX

Parameters
rXThe vector considered
rYThe copy of the vector rX

◆ CopyMatrixValues()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::CopyMatrixValues ( MatrixType rA,
const MatrixType rB 
)
inlinestatic

Copy values from one matrix to another.

It is assumed that the sparcity of both matrices is compatible

Parameters
rAThe matrix where assigning values
rBThe matrix to be copied

◆ CreateDofUpdater()

template<class TMatrixType , class TVectorType >
static DofUpdaterPointerType Kratos::TrilinosSpace< TMatrixType, TVectorType >::CreateDofUpdater ( )
inlinestatic

Creates a new dof updater.

Returns
The new dof updater

◆ CreateEmptyMatrixPointer() [1/2]

template<class TMatrixType , class TVectorType >
static MatrixPointerType Kratos::TrilinosSpace< TMatrixType, TVectorType >::CreateEmptyMatrixPointer ( )
inlinestatic

This method creates an empty pointer to a matrix.

Returns
The pointer to the matrix

◆ CreateEmptyMatrixPointer() [2/2]

template<class TMatrixType , class TVectorType >
static MatrixPointerType Kratos::TrilinosSpace< TMatrixType, TVectorType >::CreateEmptyMatrixPointer ( Epetra_MpiComm &  rComm)
inlinestatic

This method creates an empty pointer to a matrix using epetra communicator.

Parameters
rCommThe epetra communicator
Returns
The pointer to the matrix

◆ CreateEmptyVectorPointer() [1/2]

template<class TMatrixType , class TVectorType >
static VectorPointerType Kratos::TrilinosSpace< TMatrixType, TVectorType >::CreateEmptyVectorPointer ( )
inlinestatic

This method creates an empty pointer to a vector.

Returns
The pointer to the vector

◆ CreateEmptyVectorPointer() [2/2]

template<class TMatrixType , class TVectorType >
static VectorPointerType Kratos::TrilinosSpace< TMatrixType, TVectorType >::CreateEmptyVectorPointer ( Epetra_MpiComm &  rComm)
inlinestatic

This method creates an empty pointer to a vector using epetra communicator.

Parameters
rCommThe epetra communicator
Returns
The pointer to the vector

◆ Dot()

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::Dot ( const VectorType rX,
const VectorType rY 
)
inlinestatic

Returns the product of two vectors.

rX * rY

Parameters
rXThe first vector considered
rYThe second vector considered

◆ GatherValues()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::GatherValues ( const VectorType rX,
const std::vector< int > &  IndexArray,
double pValues 
)
inlinestatic

This function gathers the values of a given vector according to a given index array.

Parameters
rXThe vector from which values are to be gathered
IndexArrayThe array containing the indices of the values to be gathered
pValuesThe array containing the gathered values

◆ GetAveragevalueDiagonal()

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::GetAveragevalueDiagonal ( const MatrixType rA)
inlinestatic

This method returns the diagonal max value.

Parameters
rAThe LHS matrix
Returns
The diagonal max value

◆ GetColumn()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::GetColumn ( const unsigned int  j,
const MatrixType rM,
VectorType rX 
)
inlinestatic

Returns the column of the matrix in the given position.

rXi = rMij

Parameters
jThe position of the column
rMThe matrix considered
rXThe column considered
Todo:
Implement this method

◆ GetDiagonalNorm()

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::GetDiagonalNorm ( const MatrixType rA)
inlinestatic

This method returns the diagonal norm considering for scaling the diagonal.

Parameters
rAThe LHS matrix
Returns
The diagonal norm

◆ GetMaxDiagonal()

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::GetMaxDiagonal ( const MatrixType rA)
inlinestatic

This method returns the diagonal max value.

Parameters
rAThe LHS matrix
Returns
The diagonal max value

◆ GetMinDiagonal()

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::GetMinDiagonal ( const MatrixType rA)
inlinestatic

This method returns the diagonal min value.

Parameters
rAThe LHS matrix
Returns
The diagonal min value

◆ GetScaleNorm()

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::GetScaleNorm ( const ProcessInfo rProcessInfo,
const MatrixType rA,
const SCALING_DIAGONAL  ScalingDiagonal = SCALING_DIAGONAL::NO_SCALING 
)
inlinestatic

This method returns the scale norm considering for scaling the diagonal.

Parameters
rProcessInfoThe problem process info
rAThe LHS matrix
ScalingDiagonalThe type of caling diagonal considered
Returns
The scale norm

◆ GetValue()

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::GetValue ( const VectorType rX,
const std::size_t  I 
)
inlinestatic

This function returns a value from a given vector according to a given index.

Parameters
rXThe vector from which values are to be gathered
IThe index of the value to be gathered
Returns
The value of the vector corresponding to the index I

◆ Info()

template<class TMatrixType , class TVectorType >
virtual std::string Kratos::TrilinosSpace< TMatrixType, TVectorType >::Info ( ) const
inlinevirtual

Turn back information as a string.

Returns
Info as a string.

◆ InplaceMult()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::InplaceMult ( VectorType rX,
const double  A 
)
inlinestatic

Returns the multiplication of a vector by a scalar.

y = A*x Checks if a multiplication is needed and tries to do otherwise

Parameters
rXThe vector considered
AThe scalar considered

◆ IsDistributed()

template<class TMatrixType , class TVectorType >
static constexpr bool Kratos::TrilinosSpace< TMatrixType, TVectorType >::IsDistributed ( )
inlinestaticconstexpr

This function returns if we are in a distributed system.

Returns
True if we are in a distributed system, false otherwise (always true in this case)

◆ IsDistributedSpace()

template<class TMatrixType , class TVectorType >
static constexpr bool Kratos::TrilinosSpace< TMatrixType, TVectorType >::IsDistributedSpace ( )
inlinestaticconstexpr

Check if the TrilinosSpace is distributed.

This static member function checks whether the TrilinosSpace is distributed or not.

Returns
True if the space is distributed, false otherwise.

◆ KRATOS_CLASS_POINTER_DEFINITION()

template<class TMatrixType , class TVectorType >
Kratos::TrilinosSpace< TMatrixType, TVectorType >::KRATOS_CLASS_POINTER_DEFINITION ( TrilinosSpace< TMatrixType, TVectorType >  )

Pointer definition of TrilinosSpace.

◆ Max()

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::Max ( const VectorType rX)
inlinestatic

Returns the maximum value of the vector rX.

Parameters
rXThe vector considered
Returns
The maximum value of the vector rX

◆ Min()

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::Min ( const VectorType rX)
inlinestatic

Returns the minimum value of the vector rX.

Parameters
rXThe vector considered
Returns
The minimum value of the vector rX

◆ Mult() [1/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::Mult ( const MatrixType rA,
const MatrixType rB,
MatrixType rC,
const bool  CallFillCompleteOnResult = true,
const bool  KeepAllHardZeros = false 
)
inlinestatic

Returns the multiplication matrix-matrix.

C = A*B

Parameters
rAThe first matrix considered
rBThe second matrix considered
rCThe result of the multiplication
CallFillCompleteOnResultOptional argument, defaults to true. Power users may specify this argument to be false if they DON'T want this function to call C.FillComplete. (It is often useful to allow this function to call C.FillComplete, in cases where one or both of the input matrices are rectangular and it is not trivial to know which maps to use for the domain- and range-maps.)
KeepAllHardZerosOptional argument, defaults to false. If true, Multiply, keeps all entries in C corresponding to hard zeros. If false, the following happens by case: A*B^T, A^T*B^T - Does not store entries caused by hard zeros in C. A^T*B (unoptimized) - Hard zeros are always stored (this option has no effect) A*B, A^T*B (optimized) - Hard zeros in corresponding to hard zeros in A are not stored, There are certain cases involving reuse of C, where this can be useful.

◆ Mult() [2/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::Mult ( const MatrixType rA,
const VectorType rX,
VectorType rY 
)
inlinestatic

Returns the multiplication of a matrix by a vector.

y = A*x

Parameters
rAThe matrix considered
rXThe vector considered
rYThe result of the multiplication

◆ PrintData()

template<class TMatrixType , class TVectorType >
virtual void Kratos::TrilinosSpace< TMatrixType, TVectorType >::PrintData ( std::ostream &  rOStream) const
inlinevirtual

Print object's data.

Parameters
rOStreamThe output stream to print on.

◆ PrintInfo()

template<class TMatrixType , class TVectorType >
virtual void Kratos::TrilinosSpace< TMatrixType, TVectorType >::PrintInfo ( std::ostream &  rOStream) const
inlinevirtual

Print information about this object.

Parameters
rOStreamThe output stream to print on.

◆ ReadMatrixMarket()

template<class TMatrixType , class TVectorType >
MatrixPointerType Kratos::TrilinosSpace< TMatrixType, TVectorType >::ReadMatrixMarket ( const std::string  FileName,
Epetra_MpiComm &  rComm 
)
inline

Read a matrix from a MatrixMarket file.

Parameters
rFileNameThe name of the file to read
rCommThe MPI communicator
Returns
The matrix read from the file

◆ ReadMatrixMarketVector()

template<class TMatrixType , class TVectorType >
VectorPointerType Kratos::TrilinosSpace< TMatrixType, TVectorType >::ReadMatrixMarketVector ( const std::string &  rFileName,
Epetra_MpiComm &  rComm,
const int  N 
)
inline

Read a vector from a MatrixMarket file.

Parameters
rFileNameThe name of the file to read
rCommThe MPI communicator
NThe size of the vector

◆ Resize() [1/3]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::Resize ( MatrixType rA,
const SizeType  m,
const SizeType  n 
)
inlinestatic

Resizes a matrix.

Parameters
rAThe matrix to be resized
mThe new number of rows
nThe new number of columns

◆ Resize() [2/3]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::Resize ( VectorPointerType pX,
const SizeType  n 
)
inlinestatic

Resizes a vector.

Parameters
pAThe pointer to the vector to be resized
nThe new size

◆ Resize() [3/3]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::Resize ( VectorType rX,
const SizeType  n 
)
inlinestatic

Resizes a vector.

Parameters
rXThe vector to be resized
nThe new size

◆ ScaleAndAdd() [1/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::ScaleAndAdd ( const double  A,
const VectorType rX,
const double  B,
const VectorType rY,
VectorType rZ 
)
inlinestatic

Returns the unaliased addition of two vectors by a scalar.

rZ = (A * rX) + (B * rY)

Parameters
AThe scalar considered
rXThe first vector considered
BThe scalar considered
rYThe second vector considered
rZThe resulting vector considered

◆ ScaleAndAdd() [2/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::ScaleAndAdd ( const double  A,
const VectorType rX,
const double  B,
VectorType rY 
)
inlinestatic

Returns the unaliased addition of two vectors by a scalar.

rY = (A * rX) + (B * rY)

Parameters
AThe scalar considered
rXThe first vector considered
BThe scalar considered
rYThe resulting vector considered

◆ Set()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::Set ( VectorType rX,
const DataType  A 
)
inlinestatic

assigns a scalar to a vector

rX = A

Parameters
rXThe vector considered
AThe scalar considered

◆ SetToZero() [1/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::SetToZero ( MatrixType rA)
inlinestatic

Sets a matrix to zero.

Parameters
rXThe matrix to be set

◆ SetToZero() [2/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::SetToZero ( VectorType rX)
inlinestatic

Sets a vector to zero.

Parameters
rXThe vector to be set

◆ SetValue()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::SetValue ( VectorType rX,
IndexType  i,
const double  value 
)
inlinestatic

Sets a value in a vector.

Parameters
rXThe vector considered
iThe index of the value considered
valueThe value considered

◆ Size()

template<class TMatrixType , class TVectorType >
static IndexType Kratos::TrilinosSpace< TMatrixType, TVectorType >::Size ( const VectorType rV)
inlinestatic

Returns size of vector rV.

Parameters
rVThe vector considered
Returns
The size of the vector

◆ Size1()

template<class TMatrixType , class TVectorType >
static IndexType Kratos::TrilinosSpace< TMatrixType, TVectorType >::Size1 ( MatrixType const &  rM)
inlinestatic

Returns number of rows of rM.

Parameters
rMThe matrix considered
Returns
The number of rows of rM

◆ Size2()

template<class TMatrixType , class TVectorType >
static IndexType Kratos::TrilinosSpace< TMatrixType, TVectorType >::Size2 ( MatrixType const &  rM)
inlinestatic

Returns number of columns of rM.

Parameters
rMThe matrix considered
Returns
The number of columns of rM

◆ TransposeMult() [1/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::TransposeMult ( const MatrixType rA,
const MatrixType rB,
MatrixType rC,
const std::pair< bool, bool TransposeFlag = {false, false},
const bool  CallFillCompleteOnResult = true,
const bool  KeepAllHardZeros = false 
)
inlinestatic

Returns the transpose multiplication matrix-matrix.

C = A*B

Parameters
rAThe first matrix considered
rBThe second matrix considered
rCThe result of the multiplication
TransposeFlagFlags to transpose the matrices
CallFillCompleteOnResultOptional argument, defaults to true. Power users may specify this argument to be false if they DON'T want this function to call C.FillComplete. (It is often useful to allow this function to call C.FillComplete, in cases where one or both of the input matrices are rectangular and it is not trivial to know which maps to use for the domain- and range-maps.)
KeepAllHardZerosOptional argument, defaults to false. If true, Multiply, keeps all entries in C corresponding to hard zeros. If false, the following happens by case: A*B^T, A^T*B^T - Does not store entries caused by hard zeros in C. A^T*B (unoptimized) - Hard zeros are always stored (this option has no effect) A*B, A^T*B (optimized) - Hard zeros in corresponding to hard zeros in A are not stored, There are certain cases involving reuse of C, where this can be useful.

◆ TransposeMult() [2/2]

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::TransposeMult ( const MatrixType rA,
const VectorType rX,
VectorType rY 
)
inlinestatic

Returns the transpose multiplication of a matrix by a vector.

y = AT*x

Parameters
rAThe matrix considered
rXThe vector considered
rYThe result of the multiplication

◆ TwoNorm() [1/2]

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::TwoNorm ( const MatrixType rA)
inlinestatic

Returns the Frobenius norm of the matrix rX.

||rA||2

Parameters
rAThe matrix considered
Returns
The Frobenius norm of the matrix rX

◆ TwoNorm() [2/2]

template<class TMatrixType , class TVectorType >
static double Kratos::TrilinosSpace< TMatrixType, TVectorType >::TwoNorm ( const VectorType rX)
inlinestatic

Returns the norm of the vector rX.

||rX||2

Parameters
rXThe vector considered
Returns
The norm of the vector rX

◆ UnaliasedAdd()

template<class TMatrixType , class TVectorType >
static void Kratos::TrilinosSpace< TMatrixType, TVectorType >::UnaliasedAdd ( VectorType rX,
const double  A,
const VectorType rY 
)
inlinestatic

Returns the unaliased addition of a vector by a scalar times a vector.

X += A*y; Checks if a multiplication is needed and tries to do otherwise

Note
ATTENTION it is assumed no aliasing between rX and rY
Parameters
rXThe resulting vector considered
AThe scalar considered
rYThe multiplied vector considered

◆ WriteMatrixMarketMatrix()

template<class TMatrixType , class TVectorType >
template<class TOtherMatrixType >
static bool Kratos::TrilinosSpace< TMatrixType, TVectorType >::WriteMatrixMarketMatrix ( const char *  pFileName,
const TOtherMatrixType &  rM,
const bool  Symmetric 
)
inlinestatic

Writes a matrix to a file in MatrixMarket format.

Parameters
pFileNameThe name of the file to be written
rMThe matrix to be written
SymmetricIf the matrix is symmetric
Returns
True if the file was successfully written, false otherwise

◆ WriteMatrixMarketVector()

template<class TMatrixType , class TVectorType >
template<class VectorType >
static bool Kratos::TrilinosSpace< TMatrixType, TVectorType >::WriteMatrixMarketVector ( const char *  pFileName,
const VectorType rV 
)
inlinestatic

Writes a vector to a file in MatrixMarket format.

Parameters
pFileNameThe name of the file to be written
rVThe vector to be written
Returns
True if the file was successfully written, false otherwise

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