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 <private/kspimpl.h>
12: static PetscErrorCode KSPSetUp_CGS(KSP ksp)
13: {
17: KSPDefaultGetWork(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(((PetscObject)ksp)->comm,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) {
55: dp *= dp;
56: }
57: PetscObjectTakeAccess(ksp);
58: ksp->its = 0;
59: ksp->rnorm = dp;
60: PetscObjectGrantAccess(ksp);
61: KSPLogResidualHistory(ksp,dp);
62: KSPMonitor(ksp,0,dp);
63: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
64: if (ksp->reason) return(0);
66: /* Make the initial Rp == R */
67: VecCopy(R,RP);
68: /* added for Fidap */
69: /* Penalize Startup - Isaac Hasbani Trick for CGS
70: Since most initial conditions result in a mostly 0 residual,
71: we change all the 0 values in the vector RP to the maximum.
72: */
73: if (ksp->normtype == KSP_NORM_NATURAL) {
74: PetscReal vr0max;
75: PetscScalar *tmp_RP=0;
76: PetscInt numnp=0, *max_pos=0;
77: VecMax(RP, max_pos, &vr0max);
78: VecGetArray(RP, &tmp_RP);
79: VecGetLocalSize(RP, &numnp);
80: for (i=0; i<numnp; i++) {
81: if (tmp_RP[i] == 0.0) tmp_RP[i] = vr0max;
82: }
83: VecRestoreArray(RP, &tmp_RP);
84: }
85: /* end of addition for Fidap */
87: /* Set the initial conditions */
88: VecDot(R,RP,&rhoold); /* rhoold = (r,rp) */
89: VecCopy(R,U);
90: VecCopy(R,P);
91: KSP_PCApplyBAorAB(ksp,P,V,T);
93: i = 0;
94: do {
96: VecDot(V,RP,&s); /* s <- (v,rp) */
97: a = rhoold / s; /* a <- rho / s */
98: VecWAXPY(Q,-a,V,U); /* q <- u - a v */
99: VecWAXPY(T,1.0,U,Q); /* t <- u + q */
100: VecAXPY(X,a,T); /* x <- x + a (u + q) */
101: KSP_PCApplyBAorAB(ksp,T,AUQ,U);
102: VecAXPY(R,-a,AUQ); /* r <- r - a K (u + q) */
103: VecDot(R,RP,&rho); /* rho <- (r,rp) */
104: if (ksp->normtype == KSP_NORM_NATURAL) {
105: dp = PetscAbsScalar(rho);
106: } else {
107: VecNorm(R,NORM_2,&dp);
108: }
110: PetscObjectTakeAccess(ksp);
111: ksp->its++;
112: ksp->rnorm = dp;
113: PetscObjectGrantAccess(ksp);
114: KSPLogResidualHistory(ksp,dp);
115: KSPMonitor(ksp,i+1,dp);
116: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
117: if (ksp->reason) break;
119: b = rho / rhoold; /* b <- rho / rhoold */
120: VecWAXPY(U,b,Q,R); /* u <- r + b q */
121: VecAXPY(Q,b,P);
122: VecWAXPY(P,b,Q,U); /* p <- u + b(q + b p) */
123: KSP_PCApplyBAorAB(ksp,P,V,Q); /* v <- K p */
124: rhoold = rho;
125: i++;
126: } while (i<ksp->max_it);
127: if (i >= ksp->max_it) {
128: ksp->reason = KSP_DIVERGED_ITS;
129: }
131: KSPUnwindPreconditioner(ksp,X,T);
132: return(0);
133: }
135: /*MC
136: KSPCGS - This code implements the CGS (Conjugate Gradient Squared) method.
138: Options Database Keys:
139: . see KSPSolve()
141: Level: beginner
143: References: Sonneveld, 1989.
145: Notes: Does not require a symmetric matrix. Does not apply transpose of the matrix.
146: Supports left and right preconditioning, but not symmetric.
148: Developer Notes: Has this weird support for doing the convergence test with the natural norm, I assume this works only with
149: no preconditioning and symmetric positive definite operator.
151: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBCGS, KSPSetPCSide()
152: M*/
156: PetscErrorCode KSPCreate_CGS(KSP ksp)
157: {
161: ksp->data = (void*)0;
162: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
163: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);
164: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,1);
165: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_RIGHT,1);
167: ksp->ops->setup = KSPSetUp_CGS;
168: ksp->ops->solve = KSPSolve_CGS;
169: ksp->ops->destroy = KSPDefaultDestroy;
170: ksp->ops->buildsolution = KSPDefaultBuildSolution;
171: ksp->ops->buildresidual = KSPDefaultBuildResidual;
172: ksp->ops->setfromoptions = 0;
173: ksp->ops->view = 0;
174: return(0);
175: }