Actual source code: cgs.c


  2: /*

  4:     Note that for the complex numbers version, the VecDot() arguments
  5:     within the code MUST remain in the order given for correct computation
  6:     of inner products.
  7: */
  8: #include <petsc/private/kspimpl.h>

 10: static PetscErrorCode KSPSetUp_CGS(KSP ksp)
 11: {
 12:   KSPSetWorkVecs(ksp, 7);
 13:   return 0;
 14: }

 16: static PetscErrorCode KSPSolve_CGS(KSP ksp)
 17: {
 18:   PetscInt    i;
 19:   PetscScalar rho, rhoold, a, s, b;
 20:   Vec         X, B, V, P, R, RP, T, Q, U, AUQ;
 21:   PetscReal   dp = 0.0;
 22:   PetscBool   diagonalscale;

 24:   /* not sure what residual norm it does use, should use for right preconditioning */

 26:   PCGetDiagonalScale(ksp->pc, &diagonalscale);

 29:   X   = ksp->vec_sol;
 30:   B   = ksp->vec_rhs;
 31:   R   = ksp->work[0];
 32:   RP  = ksp->work[1];
 33:   V   = ksp->work[2];
 34:   T   = ksp->work[3];
 35:   Q   = ksp->work[4];
 36:   P   = ksp->work[5];
 37:   U   = ksp->work[6];
 38:   AUQ = V;

 40:   /* Compute initial preconditioned residual */
 41:   KSPInitialResidual(ksp, X, V, T, R, B);

 43:   /* Test for nothing to do */
 44:   if (ksp->normtype != KSP_NORM_NONE) {
 45:     VecNorm(R, NORM_2, &dp);
 46:     KSPCheckNorm(ksp, dp);
 47:     if (ksp->normtype == KSP_NORM_NATURAL) dp *= dp;
 48:   } else dp = 0.0;

 50:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
 51:   ksp->its   = 0;
 52:   ksp->rnorm = dp;
 53:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
 54:   KSPLogResidualHistory(ksp, dp);
 55:   KSPMonitor(ksp, 0, dp);
 56:   (*ksp->converged)(ksp, 0, dp, &ksp->reason, ksp->cnvP);
 57:   if (ksp->reason) return 0;

 59:   /* Make the initial Rp == R */
 60:   VecCopy(R, RP);
 61:   /*  added for Fidap */
 62:   /* Penalize Startup - Isaac Hasbani Trick for CGS
 63:      Since most initial conditions result in a mostly 0 residual,
 64:      we change all the 0 values in the vector RP to the maximum.
 65:   */
 66:   if (ksp->normtype == KSP_NORM_NATURAL) {
 67:     PetscReal    vr0max;
 68:     PetscScalar *tmp_RP = NULL;
 69:     PetscInt     numnp = 0, *max_pos = NULL;
 70:     VecMax(RP, max_pos, &vr0max);
 71:     VecGetArray(RP, &tmp_RP);
 72:     VecGetLocalSize(RP, &numnp);
 73:     for (i = 0; i < numnp; i++) {
 74:       if (tmp_RP[i] == 0.0) tmp_RP[i] = vr0max;
 75:     }
 76:     VecRestoreArray(RP, &tmp_RP);
 77:   }
 78:   /*  end of addition for Fidap */

 80:   /* Set the initial conditions */
 81:   VecDot(R, RP, &rhoold); /* rhoold = (r,rp)      */
 82:   VecCopy(R, U);
 83:   VecCopy(R, P);
 84:   KSP_PCApplyBAorAB(ksp, P, V, T);

 86:   i = 0;
 87:   do {
 88:     VecDot(V, RP, &s); /* s <- (v,rp)          */
 89:     KSPCheckDot(ksp, s);
 90:     a = rhoold / s;                    /* a <- rho / s         */
 91:     VecWAXPY(Q, -a, V, U);  /* q <- u - a v         */
 92:     VecWAXPY(T, 1.0, U, Q); /* t <- u + q           */
 93:     VecAXPY(X, a, T);       /* x <- x + a (u + q)   */
 94:     KSP_PCApplyBAorAB(ksp, T, AUQ, U);
 95:     VecAXPY(R, -a, AUQ); /* r <- r - a K (u + q) */
 96:     VecDot(R, RP, &rho); /* rho <- (r,rp)        */
 97:     KSPCheckDot(ksp, rho);
 98:     if (ksp->normtype == KSP_NORM_NATURAL) {
 99:       dp = PetscAbsScalar(rho);
100:     } else if (ksp->normtype != KSP_NORM_NONE) {
101:       VecNorm(R, NORM_2, &dp);
102:       KSPCheckNorm(ksp, dp);
103:     } else dp = 0.0;

105:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
106:     ksp->its++;
107:     ksp->rnorm = dp;
108:     PetscObjectSAWsGrantAccess((PetscObject)ksp);
109:     KSPLogResidualHistory(ksp, dp);
110:     KSPMonitor(ksp, i + 1, dp);
111:     (*ksp->converged)(ksp, i + 1, dp, &ksp->reason, ksp->cnvP);
112:     if (ksp->reason) break;

114:     b = rho / rhoold;                /* b <- rho / rhoold    */
115:     VecWAXPY(U, b, Q, R); /* u <- r + b q         */
116:     VecAXPY(Q, b, P);
117:     VecWAXPY(P, b, Q, U);            /* p <- u + b(q + b p)  */
118:     KSP_PCApplyBAorAB(ksp, P, V, Q); /* v <- K p    */
119:     rhoold = rho;
120:     i++;
121:   } while (i < ksp->max_it);
122:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;

124:   KSPUnwindPreconditioner(ksp, X, T);
125:   return 0;
126: }

128: /*MC
129:      KSPCGS - This code implements the CGS (Conjugate Gradient Squared) method.

131:    Options Database Keys:
132:     see KSPSolve()

134:    Level: beginner

136:    References:
137: .  * - Sonneveld, 1989.

139:    Notes:
140:     Does not require a symmetric matrix. Does not apply transpose of the matrix.
141:         Supports left and right preconditioning, but not symmetric.

143:    Developer Notes:
144:     Has this weird support for doing the convergence test with the natural norm, I assume this works only with
145:       no preconditioning and symmetric positive definite operator.

147: .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPBCGS`, `KSPSetPCSide()`
148: M*/
149: PETSC_EXTERN PetscErrorCode KSPCreate_CGS(KSP ksp)
150: {
151:   ksp->data = (void *)0;

153:   KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3);
154:   KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2);
155:   KSPSetSupportedNorm(ksp, KSP_NORM_NATURAL, PC_LEFT, 2);
156:   KSPSetSupportedNorm(ksp, KSP_NORM_NATURAL, PC_RIGHT, 2);
157:   KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1);
158:   KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1);

160:   ksp->ops->setup          = KSPSetUp_CGS;
161:   ksp->ops->solve          = KSPSolve_CGS;
162:   ksp->ops->destroy        = KSPDestroyDefault;
163:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
164:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
165:   ksp->ops->setfromoptions = NULL;
166:   ksp->ops->view           = NULL;
167:   return 0;
168: }