Chi-Tech
petsc_utils.h
Go to the documentation of this file.
1#ifndef CHI_MATH_PETSC_UTILS_H
2#define CHI_MATH_PETSC_UTILS_H
3
4#include <petscksp.h>
5#include <vector>
6
7namespace chi_math
8{
9class ParallelVector;
10}
11
13{
14 /**Generalized solver structure.*/
16 {
17 KSP ksp;
18 PC pc;
19
20 std::string in_solver_name = "KSPSolver";
21
22 std::string solver_type = KSPGMRES;
23 std::string preconditioner_type = PCNONE;
24
25 double relative_residual_tol = 1.0e-6;
27 };
28
29 //01
30 Vec CreateVector(int64_t local_size, int64_t global_size);
31 void CreateVector(Vec& x, int64_t local_size, int64_t global_size);
32
33 Vec CreateVectorWithGhosts(int64_t local_size, int64_t global_size,
34 int64_t nghosts,
35 const std::vector<int64_t>& ghost_indices);
36
37 //02
38 Mat CreateSquareMatrix(int64_t local_size, int64_t global_size);
39 void CreateSquareMatrix(Mat& A, int64_t local_size, int64_t global_size);
40 void InitMatrixSparsity(Mat &A,
41 const std::vector<int64_t>& nodal_nnz_in_diag,
42 const std::vector<int64_t>& nodal_nnz_off_diag);
43 void InitMatrixSparsity(Mat &A,
44 int64_t nodal_nnz_in_diag,
45 int64_t nodal_nnz_off_diag);
46
47 //03
49 Mat ref_matrix,
50 const std::string& in_solver_name = "KSPSolver",
51 const std::string& in_solver_type = KSPGMRES,
52 const std::string& in_preconditioner_type = PCNONE,
53 double in_relative_residual_tolerance = 1.0e-6,
54 int64_t in_maximum_iterations = 100);
55
57 KSP ksp, PetscInt n,
58 PetscReal rnorm,
59 KSPConvergedReason* convergedReason,
60 void *monitordestroy);
61
62 PetscErrorCode
63 KSPMonitorRelativeToRHS(KSP ksp, PetscInt n,
64 PetscReal rnorm, void *);
65 PetscErrorCode
66 KSPMonitorStraight(KSP ksp, PetscInt n,
67 PetscReal rnorm, void *);
68
69 //04
70 void CopyVecToSTLvector(Vec x, std::vector<double>& data, size_t N, bool resize_STL = true);
71 void CopyVecToSTLvectorWithGhosts(Vec x, std::vector<double>& data, size_t N, bool resize_STL = true);
72
73 void CopySTLvectorToVec(const std::vector<double>& data, Vec x, size_t N);
75
77 Vec x,
78 const std::vector<int64_t>& global_indices,
79 std::vector<double>& data);
80
82
83 /**Simple data structure to keep track of a ghost vector's
84 * localized views.*/
86 {
87 Vec x_localized = nullptr;
88 double* x_localized_raw = nullptr;
89
90 /**Returns a copy of the value at the specified index.*/
91 double operator[](int index)
92 {return x_localized_raw[index];}
93
94 /**Returns a reference of the value at the specified index.*/
95 double& operator()(int index)
96 {return x_localized_raw[index];}
97 };
98
99 GhostVecLocalRaw GetGhostVectorLocalViewRead(Vec x);
100 void RestoreGhostVectorLocalViewRead(Vec x,GhostVecLocalRaw& local_data);
101
102
103}//namespace chi_math::PETScUtils
104
105#endif
void RestoreGhostVectorLocalViewRead(Vec x, GhostVecLocalRaw &local_data)
GhostVecLocalRaw GetGhostVectorLocalViewRead(Vec x)
void InitMatrixSparsity(Mat &A, const std::vector< int64_t > &nodal_nnz_in_diag, const std::vector< int64_t > &nodal_nnz_off_diag)
void CopyVecToSTLvectorWithGhosts(Vec x, std::vector< double > &data, size_t N, bool resize_STL=true)
PetscErrorCode KSPMonitorStraight(KSP ksp, PetscInt n, PetscReal rnorm, void *)
void CopyGlobalVecToSTLvector(Vec x, const std::vector< int64_t > &global_indices, std::vector< double > &data)
PetscErrorCode RelativeResidualConvergenceTest(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *convergedReason, void *monitordestroy)
void CopyParallelVectorToVec(const ParallelVector &y, Vec x)
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)
void CopyVecToSTLvector(Vec x, std::vector< double > &data, size_t N, bool resize_STL=true)
PetscErrorCode KSPMonitorRelativeToRHS(KSP ksp, PetscInt n, PetscReal rnorm, void *)
Vec CreateVectorWithGhosts(int64_t local_size, int64_t global_size, int64_t nghosts, const std::vector< int64_t > &ghost_indices)
void CopySTLvectorToVec(const std::vector< double > &data, Vec x, size_t N)
Mat CreateSquareMatrix(int64_t local_size, int64_t global_size)
Vec CreateVector(int64_t local_size, int64_t global_size)
struct _p_Mat * Mat
struct _p_Vec * Vec