Actual source code: cgs.c
petsc-3.5.4 2015-05-23
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>
12: static PetscErrorCode KSPSetUp_CGS(KSP ksp)
13: {
17: KSPSetWorkVecs(ksp,7);
18: return(0);
19: }
23: static PetscErrorCode KSPSolve_CGS(KSP ksp)
24: {
26: PetscInt i;
27: PetscScalar rho,rhoold,a,s,b;
28: Vec X,B,V,P,R,RP,T,Q,U,AUQ;
29: PetscReal dp = 0.0;
30: PetscBool diagonalscale;
33: /* not sure what residual norm it does use, should use for right preconditioning */
35: PCGetDiagonalScale(ksp->pc,&diagonalscale);
36: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
38: X = ksp->vec_sol;
39: B = ksp->vec_rhs;
40: R = ksp->work[0];
41: RP = ksp->work[1];
42: V = ksp->work[2];
43: T = ksp->work[3];
44: Q = ksp->work[4];
45: P = ksp->work[5];
46: U = ksp->work[6];
47: AUQ = V;
49: /* Compute initial preconditioned residual */
50: KSPInitialResidual(ksp,X,V,T,R,B);
52: /* Test for nothing to do */
53: VecNorm(R,NORM_2,&dp);
54: if (ksp->normtype == KSP_NORM_NATURAL) dp *= dp;
55: PetscObjectSAWsTakeAccess((PetscObject)ksp);
56: ksp->its = 0;
57: ksp->rnorm = dp;
58: PetscObjectSAWsGrantAccess((PetscObject)ksp);
59: KSPLogResidualHistory(ksp,dp);
60: KSPMonitor(ksp,0,dp);
61: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
62: if (ksp->reason) return(0);
64: /* Make the initial Rp == R */
65: VecCopy(R,RP);
66: /* added for Fidap */
67: /* Penalize Startup - Isaac Hasbani Trick for CGS
68: Since most initial conditions result in a mostly 0 residual,
69: we change all the 0 values in the vector RP to the maximum.
70: */
71: if (ksp->normtype == KSP_NORM_NATURAL) {
72: PetscReal vr0max;
73: PetscScalar *tmp_RP=0;
74: PetscInt numnp =0, *max_pos=0;
75: VecMax(RP, max_pos, &vr0max);
76: VecGetArray(RP, &tmp_RP);
77: VecGetLocalSize(RP, &numnp);
78: for (i=0; i<numnp; i++) {
79: if (tmp_RP[i] == 0.0) tmp_RP[i] = vr0max;
80: }
81: VecRestoreArray(RP, &tmp_RP);
82: }
83: /* end of addition for Fidap */
85: /* Set the initial conditions */
86: VecDot(R,RP,&rhoold); /* rhoold = (r,rp) */
87: VecCopy(R,U);
88: VecCopy(R,P);
89: KSP_PCApplyBAorAB(ksp,P,V,T);
91: i = 0;
92: do {
94: VecDot(V,RP,&s); /* s <- (v,rp) */
95: a = rhoold / s; /* a <- rho / s */
96: VecWAXPY(Q,-a,V,U); /* q <- u - a v */
97: VecWAXPY(T,1.0,U,Q); /* t <- u + q */
98: VecAXPY(X,a,T); /* x <- x + a (u + q) */
99: KSP_PCApplyBAorAB(ksp,T,AUQ,U);
100: VecAXPY(R,-a,AUQ); /* r <- r - a K (u + q) */
101: VecDot(R,RP,&rho); /* rho <- (r,rp) */
102: if (ksp->normtype == KSP_NORM_NATURAL) {
103: dp = PetscAbsScalar(rho);
104: } else {
105: VecNorm(R,NORM_2,&dp);
106: }
108: PetscObjectSAWsTakeAccess((PetscObject)ksp);
109: ksp->its++;
110: ksp->rnorm = dp;
111: PetscObjectSAWsGrantAccess((PetscObject)ksp);
112: KSPLogResidualHistory(ksp,dp);
113: KSPMonitor(ksp,i+1,dp);
114: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
115: if (ksp->reason) break;
117: b = rho / rhoold; /* b <- rho / rhoold */
118: VecWAXPY(U,b,Q,R); /* u <- r + b q */
119: VecAXPY(Q,b,P);
120: VecWAXPY(P,b,Q,U); /* p <- u + b(q + b p) */
121: KSP_PCApplyBAorAB(ksp,P,V,Q); /* v <- K p */
122: rhoold = rho;
123: i++;
124: } while (i<ksp->max_it);
125: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
127: KSPUnwindPreconditioner(ksp,X,T);
128: return(0);
129: }
131: /*MC
132: KSPCGS - This code implements the CGS (Conjugate Gradient Squared) method.
134: Options Database Keys:
135: . see KSPSolve()
137: Level: beginner
139: References: Sonneveld, 1989.
141: Notes: Does not require a symmetric matrix. Does not apply transpose of the matrix.
142: Supports left and right preconditioning, but not symmetric.
144: Developer Notes: 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 (for list of available types), KSP, KSPBCGS, KSPSetPCSide()
148: M*/
151: PETSC_EXTERN PetscErrorCode KSPCreate_CGS(KSP ksp)
152: {
156: ksp->data = (void*)0;
158: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
159: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
160: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
161: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_RIGHT,2);
163: ksp->ops->setup = KSPSetUp_CGS;
164: ksp->ops->solve = KSPSolve_CGS;
165: ksp->ops->destroy = KSPDestroyDefault;
166: ksp->ops->buildsolution = KSPBuildSolutionDefault;
167: ksp->ops->buildresidual = KSPBuildResidualDefault;
168: ksp->ops->setfromoptions = 0;
169: ksp->ops->view = 0;
170: return(0);
171: }