Chi-Tech
chi_math::PETScUtils Namespace Reference

Data Structures

struct  GhostVecLocalRaw
 
struct  PETScSolverSetup
 

Functions

Vec CreateVector (int64_t local_size, int64_t global_size)
 
void CreateVector (Vec &x, int64_t local_size, int64_t global_size)
 
Vec CreateVectorWithGhosts (int64_t local_size, int64_t global_size, int64_t nghosts, const std::vector< int64_t > &ghost_indices)
 
Mat CreateSquareMatrix (int64_t local_size, int64_t global_size)
 
void CreateSquareMatrix (Mat &A, int64_t local_size, int64_t global_size)
 
void InitMatrixSparsity (Mat &A, const std::vector< int64_t > &nodal_nnz_in_diag, const std::vector< int64_t > &nodal_nnz_off_diag)
 
void InitMatrixSparsity (Mat &A, int64_t nodal_nnz_in_diag, int64_t nodal_nnz_off_diag)
 
PETScSolverSetup CreateCommonKrylovSolverSetup (Mat ref_matrix, const std::string &in_solver_name="KSPSolver", const std::string &in_solver_type=KSPGMRES, const std::string &in_preconditioner_type=PCNONE, double in_relative_residual_tolerance=1.0e-6, int64_t in_maximum_iterations=100)
 
PetscErrorCode RelativeResidualConvergenceTest (KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *convergedReason, void *monitordestroy)
 
PetscErrorCode KSPMonitorRelativeToRHS (KSP ksp, PetscInt n, PetscReal rnorm, void *)
 
PetscErrorCode KSPMonitorStraight (KSP ksp, PetscInt n, PetscReal rnorm, void *)
 
void CopyVecToSTLvector (Vec x, std::vector< double > &data, size_t N, bool resize_STL=true)
 
void CopyVecToSTLvectorWithGhosts (Vec x, std::vector< double > &data, size_t N, bool resize_STL=true)
 
void CopySTLvectorToVec (const std::vector< double > &data, Vec x, size_t N)
 
void CopyParallelVectorToVec (const ParallelVector &y, Vec x)
 
void CopyGlobalVecToSTLvector (Vec x, const std::vector< int64_t > &global_indices, std::vector< double > &data)
 
void CommunicateGhostEntries (Vec x)
 
GhostVecLocalRaw GetGhostVectorLocalViewRead (Vec x)
 
void RestoreGhostVectorLocalViewRead (Vec x, GhostVecLocalRaw &local_data)
 

Function Documentation

◆ CommunicateGhostEntries()

void chi_math::PETScUtils::CommunicateGhostEntries ( Vec  x)

Communicates ghost entries of a ghost vector. This operation is suitable when only a single vector is communicated. When more than vector is communicated it would be more efficient to "Begin" all the vectors followed by and "End" of each vector.

Definition at line 134 of file petsc_utils_04_vecmanip.cc.

◆ CopyGlobalVecToSTLvector()

void chi_math::PETScUtils::CopyGlobalVecToSTLvector ( Vec  x,
const std::vector< int64_t > &  global_indices,
std::vector< double > &  data 
)

Copies global values from a PETSc vector to a STL vector.

Definition at line 81 of file petsc_utils_04_vecmanip.cc.

◆ CopyParallelVectorToVec()

void chi_math::PETScUtils::CopyParallelVectorToVec ( const ParallelVector y,
Vec  x 
)

Definition at line 69 of file petsc_utils_04_vecmanip.cc.

◆ CopySTLvectorToVec()

void chi_math::PETScUtils::CopySTLvectorToVec ( const std::vector< double > &  data,
Vec  x,
size_t  N 
)

Definition at line 57 of file petsc_utils_04_vecmanip.cc.

◆ CopyVecToSTLvector()

void chi_math::PETScUtils::CopyVecToSTLvector ( Vec  x,
std::vector< double > &  data,
size_t  N,
bool  resize_STL = true 
)

Copies a PETSc vector to a STL vector. Only the local portion is copied.

Definition at line 10 of file petsc_utils_04_vecmanip.cc.

◆ CopyVecToSTLvectorWithGhosts()

void chi_math::PETScUtils::CopyVecToSTLvectorWithGhosts ( Vec  x,
std::vector< double > &  data,
size_t  N,
bool  resize_STL = true 
)

Copies a PETSc vector to a STL vector. Only the local portion is copied.

Definition at line 36 of file petsc_utils_04_vecmanip.cc.

◆ CreateCommonKrylovSolverSetup()

chi_math::PETScUtils::PETScSolverSetup chi_math::PETScUtils::CreateCommonKrylovSolverSetup ( Mat  ref_matrix,
const std::string &  in_solver_name = "KSPSolver",
const std::string &  in_solver_type = KSPGMRES,
const std::string &  in_preconditioner_type = PCNONE,
double  in_relative_residual_tolerance = 1.0e-6,
int64_t  in_maximum_iterations = 100 
)

Creates a common Krylov-solver setup.

This is a macro for:

KSPCreate(PETSC_COMM_WORLD,&setup.ksp);
KSPSetOperators(setup.ksp,ref_matrix,ref_matrix);
KSPSetType(setup.ksp,in_solver_type.c_str());
KSPSetOptionsPrefix(setup.ksp,in_solver_name.c_str());
KSPGetPC(setup.ksp,&setup.pc);
PCSetType(setup.pc,in_preconditioner_type.c_str());
KSPSetTolerances(setup.ksp,1.e-50,
in_relative_residual_tolerance,1.0e50,
in_maximum_iterations);
KSPSetInitialGuessNonzero(setup.ksp,PETSC_TRUE);
KSPMonitorSet(setup.ksp,&KSPMonitorRelativeToRHS,NULL,NULL);
KSPSetConvergenceTest(setup.ksp,&RelativeResidualConvergenceTest,NULL,NULL);
return setup;
PetscErrorCode RelativeResidualConvergenceTest(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *convergedReason, void *monitordestroy)
PetscErrorCode KSPMonitorRelativeToRHS(KSP ksp, PetscInt n, PetscReal rnorm, void *)

Definition at line 35 of file petsc_utils_03_createKSP.cc.

◆ CreateSquareMatrix() [1/2]

Mat chi_math::PETScUtils::CreateSquareMatrix ( int64_t  local_size,
int64_t  global_size 
)

Creates a general square matrix.

This is a macro for:

Mat A;
MatCreate(PETSC_COMM_WORLD,&A);
MatSetType(A,MATMPIAIJ);
MatSetSizes(A,local_size, local_size,
global_size, global_size);
MatMPIAIJSetPreallocation(A,1, nullptr,
0, nullptr);
MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
MatSetOption(A, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);
return A;
struct _p_Mat * Mat

Definition at line 25 of file petsc_utils_02_createMat.cc.

◆ CreateSquareMatrix() [2/2]

void chi_math::PETScUtils::CreateSquareMatrix ( Mat A,
int64_t  local_size,
int64_t  global_size 
)

Creates a general square matrix.

This is a macro for:

MatCreate(PETSC_COMM_WORLD,&A);
MatSetType(A,MATMPIAIJ);
MatSetSizes(A,local_size, local_size,
global_size, global_size);
MatMPIAIJSetPreallocation(A,1, nullptr,
0, nullptr);
MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
MatSetOption(A, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);

Definition at line 58 of file petsc_utils_02_createMat.cc.

◆ CreateVector() [1/2]

Vec chi_math::PETScUtils::CreateVector ( int64_t  local_size,
int64_t  global_size 
)

Creates a general vector.

This is a macro for:

Vec x;
VecCreate(PETSC_COMM_WORLD,&x);
VecSetType(x,VECMPI);
VecSetSizes(x, local_size, global_size);
VecSetOption(x,VEC_IGNORE_NEGATIVE_INDICES,PETSC_TRUE);
return x;
struct _p_Vec * Vec

Definition at line 18 of file petsc_utils_01_createVec.cc.

◆ CreateVector() [2/2]

void chi_math::PETScUtils::CreateVector ( Vec x,
int64_t  local_size,
int64_t  global_size 
)

Creates a general vector.

This is a macro for:

VecCreate(PETSC_COMM_WORLD,&x);
VecSetType(x,VECMPI);
VecSetSizes(x, local_size, global_size);
VecSetOption(x,VEC_IGNORE_NEGATIVE_INDICES,PETSC_TRUE);

Definition at line 40 of file petsc_utils_01_createVec.cc.

◆ CreateVectorWithGhosts()

Vec chi_math::PETScUtils::CreateVectorWithGhosts ( int64_t  local_size,
int64_t  global_size,
int64_t  nghosts,
const std::vector< int64_t > &  ghost_indices 
)

Creates a general vector with ghost value support.

This is a macro for:

Vec x;
VecCreateGhost(PETSC_COMM_WORLD,
local_size,
global_size,
nghosts,
ghost_indices.data(),
&x);
VecSetOption(x,VEC_IGNORE_NEGATIVE_INDICES,PETSC_TRUE);
return x;

Definition at line 66 of file petsc_utils_01_createVec.cc.

◆ GetGhostVectorLocalViewRead()

chi_math::PETScUtils::GhostVecLocalRaw chi_math::PETScUtils::GetGhostVectorLocalViewRead ( Vec  x)

Gets a local raw view of a ghost vector.

Definition at line 143 of file petsc_utils_04_vecmanip.cc.

◆ InitMatrixSparsity() [1/2]

void chi_math::PETScUtils::InitMatrixSparsity ( Mat A,
const std::vector< int64_t > &  nodal_nnz_in_diag,
const std::vector< int64_t > &  nodal_nnz_off_diag 
)

Initializes the sparsity pattern of a matrix.

This is a macro for:

MatMPIAIJSetPreallocation(A,0,nodal_nnz_in_diag.data(),
0,nodal_nnz_off_diag.data());
MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
MatSetOption(A, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);
MatSetUp(A);

Definition at line 84 of file petsc_utils_02_createMat.cc.

◆ InitMatrixSparsity() [2/2]

void chi_math::PETScUtils::InitMatrixSparsity ( Mat A,
int64_t  nodal_nnz_in_diag,
int64_t  nodal_nnz_off_diag 
)

Initializes the sparsity pattern of a matrix.

This is a macro for:

MatMPIAIJSetPreallocation(A,0,nodal_nnz_in_diag.data(),
0,nodal_nnz_off_diag.data());
MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
MatSetOption(A, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);
MatSetUp(A);

Definition at line 108 of file petsc_utils_02_createMat.cc.

◆ KSPMonitorRelativeToRHS()

PetscErrorCode chi_math::PETScUtils::KSPMonitorRelativeToRHS ( KSP  ksp,
PetscInt  n,
PetscReal  rnorm,
void *   
)

General monitor that print the residual norm relative to the right-hand side norm.

Definition at line 111 of file petsc_utils_03_createKSP.cc.

◆ KSPMonitorStraight()

PetscErrorCode chi_math::PETScUtils::KSPMonitorStraight ( KSP  ksp,
PetscInt  n,
PetscReal  rnorm,
void *   
)

General monitor that print the residual norm relative to the right-hand side norm.

Definition at line 149 of file petsc_utils_03_createKSP.cc.

◆ RelativeResidualConvergenceTest()

PetscErrorCode chi_math::PETScUtils::RelativeResidualConvergenceTest ( KSP  ksp,
PetscInt  n,
PetscReal  rnorm,
KSPConvergedReason *  convergedReason,
void *  monitordestroy 
)

Relative Residual Convergence test. The test uses the L2-norm of the residual ( $ ||b-Ax||_2 $), divided by by the L2-norm of the right hand side ( $ ||b||_2 $) compared to a tolerance, $ \epsilon $.

\[ \frac{||b-Ax||_2}{||b||_2} < \epsilon \]

Definition at line 81 of file petsc_utils_03_createKSP.cc.

◆ RestoreGhostVectorLocalViewRead()

void chi_math::PETScUtils::RestoreGhostVectorLocalViewRead ( Vec  x,
GhostVecLocalRaw local_data 
)

Gets a local raw view of a ghost vector.

Definition at line 160 of file petsc_utils_04_vecmanip.cc.