Actual source code: cgs.c
1: #define PETSCKSP_DLL
3: /*
5: Note that for the complex numbers version, the VecDot() arguments
6: within the code MUST remain in the order given for correct computation
7: of inner products.
8: */
9: #include private/kspimpl.h
13: static PetscErrorCode KSPSetUp_CGS(KSP ksp)
14: {
18: if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPCGS");
19: KSPDefaultGetWork(ksp,7);
20: return(0);
21: }
25: static PetscErrorCode KSPSolve_CGS(KSP ksp)
26: {
28: PetscInt i;
29: PetscScalar rho,rhoold,a,s,b;
30: Vec X,B,V,P,R,RP,T,Q,U,AUQ;
31: PetscReal dp = 0.0;
32: PetscTruth diagonalscale;
35: /* not sure what residual norm it does use, should use for right preconditioning */
37: PCDiagonalScale(ksp->pc,&diagonalscale);
38: if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
40: X = ksp->vec_sol;
41: B = ksp->vec_rhs;
42: R = ksp->work[0];
43: RP = ksp->work[1];
44: V = ksp->work[2];
45: T = ksp->work[3];
46: Q = ksp->work[4];
47: P = ksp->work[5];
48: U = ksp->work[6];
49: AUQ = V;
51: /* Compute initial preconditioned residual */
52: KSPInitialResidual(ksp,X,V,T,R,B);
54: /* Test for nothing to do */
55: VecNorm(R,NORM_2,&dp);
56: if (ksp->normtype == KSP_NORM_NATURAL) {
57: dp *= dp;
58: }
59: PetscObjectTakeAccess(ksp);
60: ksp->its = 0;
61: ksp->rnorm = dp;
62: PetscObjectGrantAccess(ksp);
63: KSPLogResidualHistory(ksp,dp);
64: KSPMonitor(ksp,0,dp);
65: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
66: if (ksp->reason) return(0);
68: /* Make the initial Rp == R */
69: VecCopy(R,RP);
70: /* added for Fidap */
71: /* Penalize Startup - Isaac Hasbani Trick for CGS
72: Since most initial conditions result in a mostly 0 residual,
73: we change all the 0 values in the vector RP to the maximum.
74: */
75: if (ksp->normtype == KSP_NORM_NATURAL) {
76: PetscReal vr0max;
77: PetscScalar *tmp_RP=0;
78: PetscInt numnp=0, *max_pos=0;
79: VecMax(RP, max_pos, &vr0max);
80: VecGetArray(RP, &tmp_RP);
81: VecGetLocalSize(RP, &numnp);
82: for (i=0; i<numnp; i++) {
83: if (tmp_RP[i] == 0.0) tmp_RP[i] = vr0max;
84: }
85: VecRestoreArray(RP, &tmp_RP);
86: }
87: /* end of addition for Fidap */
89: /* Set the initial conditions */
90: VecDot(R,RP,&rhoold); /* rhoold = (r,rp) */
91: VecCopy(R,U);
92: VecCopy(R,P);
93: KSP_PCApplyBAorAB(ksp,P,V,T);
95: i = 0;
96: do {
98: VecDot(V,RP,&s); /* s <- (v,rp) */
99: a = rhoold / s; /* a <- rho / s */
100: VecWAXPY(Q,-a,V,U); /* q <- u - a v */
101: VecWAXPY(T,1.0,U,Q); /* t <- u + q */
102: VecAXPY(X,a,T); /* x <- x + a (u + q) */
103: KSP_PCApplyBAorAB(ksp,T,AUQ,U);
104: VecAXPY(R,-a,AUQ); /* r <- r - a K (u + q) */
105: VecDot(R,RP,&rho); /* rho <- (r,rp) */
106: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
107: VecNorm(R,NORM_2,&dp);
108: } else if (ksp->normtype == KSP_NORM_NATURAL) {
109: dp = PetscAbsScalar(rho);
110: }
112: PetscObjectTakeAccess(ksp);
113: ksp->its++;
114: ksp->rnorm = dp;
115: PetscObjectGrantAccess(ksp);
116: KSPLogResidualHistory(ksp,dp);
117: KSPMonitor(ksp,i+1,dp);
118: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
119: if (ksp->reason) break;
121: b = rho / rhoold; /* b <- rho / rhoold */
122: VecWAXPY(U,b,Q,R); /* u <- r + b q */
123: VecAXPY(Q,b,P);
124: VecWAXPY(P,b,Q,U); /* p <- u + b(q + b p) */
125: KSP_PCApplyBAorAB(ksp,P,V,Q); /* v <- K p */
126: rhoold = rho;
127: i++;
128: } while (i<ksp->max_it);
129: if (i >= ksp->max_it) {
130: ksp->reason = KSP_DIVERGED_ITS;
131: }
133: KSPUnwindPreconditioner(ksp,X,T);
134: return(0);
135: }
137: /*MC
138: KSPCGS - This code implements the CGS (Conjugate Gradient Squared) method.
140: Options Database Keys:
141: . see KSPSolve()
143: Level: beginner
145: Notes: Reference: Sonneveld, 1989.
147: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBCGS
148: M*/
152: PetscErrorCode KSPCreate_CGS(KSP ksp)
153: {
155: ksp->data = (void*)0;
156: ksp->pc_side = PC_LEFT;
157: ksp->ops->setup = KSPSetUp_CGS;
158: ksp->ops->solve = KSPSolve_CGS;
159: ksp->ops->destroy = KSPDefaultDestroy;
160: ksp->ops->buildsolution = KSPDefaultBuildSolution;
161: ksp->ops->buildresidual = KSPDefaultBuildResidual;
162: ksp->ops->setfromoptions = 0;
163: ksp->ops->view = 0;
164: return(0);
165: }