Actual source code: cgtype.c
petsc-3.9.4 2018-09-11
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: .keywords: CG, conjugate gradient, Hermitian, symmetric, set, type
30: .seealso: KSP, KSPCG
31: @*/
32: PetscErrorCode KSPCGSetType(KSP ksp,KSPCGType type)
33: {
38: PetscTryMethod(ksp,"KSPCGSetType_C",(KSP,KSPCGType),(ksp,type));
39: return(0);
40: }
42: /*@
43: KSPCGUseSingleReduction - Merge the two inner products needed in CG into a single MPI_Allreduce() call.
45: Logically Collective on KSP
47: Input Parameters:
48: + ksp - the iterative context
49: - flg - turn on or off the single reduction
51: Options Database:
52: . -ksp_cg_single_reduction
54: Level: intermediate
56: The algorithm used in this case is described as Method 1 in Lapack Working Note 56, "Conjugate Gradient Algorithms with Reduced Synchronization Overhead
57: Distributed Memory Multiprocessors", by E. F. D'Azevedo, V. L. Eijkhout, and C. H. Romine, December 3, 1999. V. Eijkhout creates the algorithm
58: initially to Chronopoulos and Gear.
60: It requires two extra work vectors than the conventional implementation in PETSc.
62: See also KSPPIPECG, KSPPIPECR, and KSPGROPPCG that use non-blocking reductions.
64: .keywords: CG, conjugate gradient, Hermitian, symmetric, set, type
66: .seealso: KSP, KSPCG, KSPGMRES
67: @*/
68: PetscErrorCode KSPCGUseSingleReduction(KSP ksp,PetscBool flg)
69: {
75: PetscTryMethod(ksp,"KSPCGUseSingleReduction_C",(KSP,PetscBool),(ksp,flg));
76: return(0);
77: }
79: /*@
80: KSPCGSetRadius - Sets the radius of the trust region.
82: Logically Collective on KSP
84: Input Parameters:
85: + ksp - the iterative context
86: - radius - the trust region radius (Infinity is the default)
88: Level: advanced
90: .keywords: set, trust region radius
92: .seealso: KSP, KSPCG, KSPNASH, KSPSTCG, KSPGLTR
93: @*/
94: PetscErrorCode KSPCGSetRadius(KSP ksp, PetscReal radius)
95: {
100: if (radius < 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE, "Radius negative");
102: PetscTryMethod(ksp,"KSPCGSetRadius_C",(KSP,PetscReal),(ksp,radius));
103: return(0);
104: }
106: /*@
107: KSPCGGetNormD - Got norm of the direction.
109: Collective on KSP
111: Input Parameters:
112: + ksp - the iterative context
113: - norm_d - the norm of the direction
115: Level: advanced
117: .keywords: get, norm direction
119: .seealso: KSP, KSPCG, KSPNASH, KSPSTCG, KSPGLTR
120: @*/
121: PetscErrorCode KSPCGGetNormD(KSP ksp, PetscReal *norm_d)
122: {
127: PetscUseMethod(ksp,"KSPCGGetNormD_C",(KSP,PetscReal*),(ksp,norm_d));
128: return(0);
129: }
131: /*@
132: KSPCGGetObjFcn - Get objective function value.
134: Collective on KSP
136: Input Parameters:
137: + ksp - the iterative context
138: - o_fcn - the objective function value
140: Level: advanced
142: .keywords: get, objective function
144: .seealso: KSP, KSPCG, KSPNASH, KSPSTCG, KSPGLTR
145: @*/
146: PetscErrorCode KSPCGGetObjFcn(KSP ksp, PetscReal *o_fcn)
147: {
152: PetscUseMethod(ksp,"KSPCGGetObjFcn_C",(KSP,PetscReal*),(ksp,o_fcn));
153: return(0);
154: }