Actual source code: cgtype.c
2: #include <../src/ksp/ksp/impls/cg/cgimpl.h>
4: /*@
5: KSPCGSetType - Sets the variant of the conjugate gradient method to
6: use for solving a linear system with a complex coefficient matrix.
7: This option is irrelevant when solving a real system.
9: Logically Collective
11: Input Parameters:
12: + ksp - the iterative context
13: - type - the variant of CG to use, one of
14: .vb
15: KSP_CG_HERMITIAN - complex, Hermitian matrix (default)
16: KSP_CG_SYMMETRIC - complex, symmetric matrix
17: .ve
19: Level: intermediate
21: Options Database Keys:
22: + -ksp_cg_hermitian - Indicates Hermitian matrix
23: - -ksp_cg_symmetric - Indicates symmetric matrix
25: Note:
26: By default, the matrix is assumed to be complex, Hermitian.
28: .seealso: [](chapter_ksp), `KSP`, `KSPCG`
29: @*/
30: PetscErrorCode KSPCGSetType(KSP ksp, KSPCGType type)
31: {
32: PetscFunctionBegin;
34: PetscTryMethod(ksp, "KSPCGSetType_C", (KSP, KSPCGType), (ksp, type));
35: PetscFunctionReturn(PETSC_SUCCESS);
36: }
38: /*@
39: KSPCGUseSingleReduction - Merge the two inner products needed in `KSPCG` into a single `MPI_Allreduce()` call.
41: Logically Collective
43: Input Parameters:
44: + ksp - the iterative context
45: - flg - turn on or off the single reduction
47: Options Database Key:
48: . -ksp_cg_single_reduction <bool> - Merge inner products into single `MPI_Allreduce()`
50: Level: intermediate
52: Notes:
53: The algorithm used in this case is described as Method 1 in [1]. V. Eijkhout credits the algorithm initially to Chronopoulos and Gear.
55: It requires two extra work vectors than the conventional implementation in PETSc.
57: See also `KSPPIPECG`, `KSPPIPECR`, and `KSPGROPPCG` that use non-blocking reductions. [](sec_pipelineksp),
59: References:
60: . [1] - Lapack Working Note 56, "Conjugate Gradient Algorithms with Reduced Synchronization Overhead
61: Distributed Memory Multiprocessors", by E. F. D'Azevedo, V. L. Eijkhout, and C. H. Romine, December 3, 1999.
63: .seealso: [](chapter_ksp), [](sec_pipelineksp), `KSP`, `KSPCG`, `KSPGMRES`, `KSPPIPECG`, `KSPPIPECR`, and `KSPGROPPCG`
64: @*/
65: PetscErrorCode KSPCGUseSingleReduction(KSP ksp, PetscBool flg)
66: {
67: PetscFunctionBegin;
70: PetscTryMethod(ksp, "KSPCGUseSingleReduction_C", (KSP, PetscBool), (ksp, flg));
71: PetscFunctionReturn(PETSC_SUCCESS);
72: }
74: /*@
75: KSPCGSetRadius - Sets the radius of the trust region when the solver is used inside `SNESNEWTONTR`
77: Logically Collective
79: Input Parameters:
80: + ksp - the iterative context
81: - radius - the trust region radius (Infinity is the default)
83: Level: advanced
85: .seealso: [](chapter_ksp), `KSP`, `KSPCG`, `KSPNASH`, `KSPSTCG`, `KSPGLTR`
86: @*/
87: PetscErrorCode KSPCGSetRadius(KSP ksp, PetscReal radius)
88: {
89: PetscFunctionBegin;
92: PetscTryMethod(ksp, "KSPCGSetRadius_C", (KSP, PetscReal), (ksp, radius));
93: PetscFunctionReturn(PETSC_SUCCESS);
94: }
96: /*@
97: KSPCGGetNormD - Got norm of the direction when the solver is used inside `SNESNEWTONTR`
99: Collective
101: Input Parameters:
102: + ksp - the iterative context
103: - norm_d - the norm of the direction
105: Level: advanced
107: .seealso: [](chapter_ksp), `KSP`, `KSPCG`, `KSPNASH`, `KSPSTCG`, `KSPGLTR`
108: @*/
109: PetscErrorCode KSPCGGetNormD(KSP ksp, PetscReal *norm_d)
110: {
111: PetscFunctionBegin;
113: PetscUseMethod(ksp, "KSPCGGetNormD_C", (KSP, PetscReal *), (ksp, norm_d));
114: PetscFunctionReturn(PETSC_SUCCESS);
115: }
117: /*@
118: KSPCGGetObjFcn - Get objective function value when the solver is used inside `SNESNEWTONTR`
120: Collective
122: Input Parameters:
123: + ksp - the iterative context
124: - o_fcn - the objective function value
126: Level: advanced
128: .seealso: [](chapter_ksp), `KSP`, `KSPCG`, `KSPNASH`, `KSPSTCG`, `KSPGLTR`
129: @*/
130: PetscErrorCode KSPCGGetObjFcn(KSP ksp, PetscReal *o_fcn)
131: {
132: PetscFunctionBegin;
134: PetscUseMethod(ksp, "KSPCGGetObjFcn_C", (KSP, PetscReal *), (ksp, o_fcn));
135: PetscFunctionReturn(PETSC_SUCCESS);
136: }