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 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: {
 33:   PetscTryMethod(ksp, "KSPCGSetType_C", (KSP, KSPCGType), (ksp, type));
 34:   return 0;
 35: }

 37: /*@
 38:     KSPCGUseSingleReduction - Merge the two inner products needed in CG into a single MPI_Allreduce() call.

 40:     Logically Collective on ksp

 42:     Input Parameters:
 43: +   ksp - the iterative context
 44: -   flg - turn on or off the single reduction

 46:     Options Database:
 47: .   -ksp_cg_single_reduction <bool> - Merge inner products into single MPI_Allreduce

 49:     Level: intermediate

 51:      The algorithm used in this case is described as Method 1 in Lapack Working Note 56, "Conjugate Gradient Algorithms with Reduced Synchronization Overhead
 52:      Distributed Memory Multiprocessors", by E. F. D'Azevedo, V. L. Eijkhout, and C. H. Romine, December 3, 1999. V. Eijkhout credits the algorithm
 53:      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.

 59: .seealso: `KSP`, `KSPCG`, `KSPGMRES`
 60: @*/
 61: PetscErrorCode KSPCGUseSingleReduction(KSP ksp, PetscBool flg)
 62: {
 65:   PetscTryMethod(ksp, "KSPCGUseSingleReduction_C", (KSP, PetscBool), (ksp, flg));
 66:   return 0;
 67: }

 69: /*@
 70:     KSPCGSetRadius - Sets the radius of the trust region.

 72:     Logically Collective on ksp

 74:     Input Parameters:
 75: +   ksp    - the iterative context
 76: -   radius - the trust region radius (Infinity is the default)

 78:     Level: advanced

 80: .seealso: `KSP`, `KSPCG`, `KSPNASH`, `KSPSTCG`, `KSPGLTR`
 81: @*/
 82: PetscErrorCode KSPCGSetRadius(KSP ksp, PetscReal radius)
 83: {
 86:   PetscTryMethod(ksp, "KSPCGSetRadius_C", (KSP, PetscReal), (ksp, radius));
 87:   return 0;
 88: }

 90: /*@
 91:     KSPCGGetNormD - Got norm of the direction.

 93:     Collective on ksp

 95:     Input Parameters:
 96: +   ksp    - the iterative context
 97: -   norm_d - the norm of the direction

 99:     Level: advanced

101: .seealso: `KSP`, `KSPCG`, `KSPNASH`, `KSPSTCG`, `KSPGLTR`
102: @*/
103: PetscErrorCode KSPCGGetNormD(KSP ksp, PetscReal *norm_d)
104: {
106:   PetscUseMethod(ksp, "KSPCGGetNormD_C", (KSP, PetscReal *), (ksp, norm_d));
107:   return 0;
108: }

110: /*@
111:     KSPCGGetObjFcn - Get objective function value.

113:     Collective on ksp

115:     Input Parameters:
116: +   ksp   - the iterative context
117: -   o_fcn - the objective function value

119:     Level: advanced

121: .seealso: `KSP`, `KSPCG`, `KSPNASH`, `KSPSTCG`, `KSPGLTR`
122: @*/
123: PetscErrorCode KSPCGGetObjFcn(KSP ksp, PetscReal *o_fcn)
124: {
126:   PetscUseMethod(ksp, "KSPCGGetObjFcn_C", (KSP, PetscReal *), (ksp, o_fcn));
127:   return 0;
128: }