Actual source code: cgtype.c
petsc-3.13.6 2020-09-29
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 on ksp
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: KSP, KSPCG
29: @*/
30: PetscErrorCode KSPCGSetType(KSP ksp,KSPCGType type)
31: {
36: PetscTryMethod(ksp,"KSPCGSetType_C",(KSP,KSPCGType),(ksp,type));
37: return(0);
38: }
40: /*@
41: KSPCGUseSingleReduction - Merge the two inner products needed in CG into a single MPI_Allreduce() call.
43: Logically Collective on ksp
45: Input Parameters:
46: + ksp - the iterative context
47: - flg - turn on or off the single reduction
49: Options Database:
50: . -ksp_cg_single_reduction
52: Level: intermediate
54: The algorithm used in this case is described as Method 1 in Lapack Working Note 56, "Conjugate Gradient Algorithms with Reduced Synchronization Overhead
55: Distributed Memory Multiprocessors", by E. F. D'Azevedo, V. L. Eijkhout, and C. H. Romine, December 3, 1999. V. Eijkhout credits the algorithm
56: initially to Chronopoulos and Gear.
58: It requires two extra work vectors than the conventional implementation in PETSc.
60: See also KSPPIPECG, KSPPIPECR, and KSPGROPPCG that use non-blocking reductions.
62: .seealso: KSP, KSPCG, KSPGMRES
63: @*/
64: PetscErrorCode KSPCGUseSingleReduction(KSP ksp,PetscBool flg)
65: {
71: PetscTryMethod(ksp,"KSPCGUseSingleReduction_C",(KSP,PetscBool),(ksp,flg));
72: return(0);
73: }
75: /*@
76: KSPCGSetRadius - Sets the radius of the trust region.
78: Logically Collective on ksp
80: Input Parameters:
81: + ksp - the iterative context
82: - radius - the trust region radius (Infinity is the default)
84: Level: advanced
86: .seealso: KSP, KSPCG, KSPNASH, KSPSTCG, KSPGLTR
87: @*/
88: PetscErrorCode KSPCGSetRadius(KSP ksp, PetscReal radius)
89: {
94: if (radius < 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE, "Radius negative");
96: PetscTryMethod(ksp,"KSPCGSetRadius_C",(KSP,PetscReal),(ksp,radius));
97: return(0);
98: }
100: /*@
101: KSPCGGetNormD - Got norm of the direction.
103: Collective on ksp
105: Input Parameters:
106: + ksp - the iterative context
107: - norm_d - the norm of the direction
109: Level: advanced
111: .seealso: KSP, KSPCG, KSPNASH, KSPSTCG, KSPGLTR
112: @*/
113: PetscErrorCode KSPCGGetNormD(KSP ksp, PetscReal *norm_d)
114: {
119: PetscUseMethod(ksp,"KSPCGGetNormD_C",(KSP,PetscReal*),(ksp,norm_d));
120: return(0);
121: }
123: /*@
124: KSPCGGetObjFcn - Get objective function value.
126: Collective on ksp
128: Input Parameters:
129: + ksp - the iterative context
130: - o_fcn - the objective function value
132: Level: advanced
134: .seealso: KSP, KSPCG, KSPNASH, KSPSTCG, KSPGLTR
135: @*/
136: PetscErrorCode KSPCGGetObjFcn(KSP ksp, PetscReal *o_fcn)
137: {
142: PetscUseMethod(ksp,"KSPCGGetObjFcn_C",(KSP,PetscReal*),(ksp,o_fcn));
143: return(0);
144: }