Actual source code: cr.c
1: /*$Id: cr.c,v 1.63 2001/03/23 23:23:33 balay Exp $*/
3: /*
4: This implements Preconditioned Conjugate Residuals.
5: */
6: #include src/sles/ksp/kspimpl.h
8: static int KSPSetUp_CR(KSP ksp)
9: {
13: if (ksp->pc_side == PC_RIGHT) {SETERRQ(2,"no right preconditioning for KSPCR");}
14: else if (ksp->pc_side == PC_SYMMETRIC) {SETERRQ(2,"no symmetric preconditioning for KSPCR");}
15: KSPDefaultGetWork(ksp,9);
16: return(0);
17: }
19: static int KSPSolve_CR(KSP ksp,int *its)
20: {
21: int i,maxit,pres,ierr;
22: MatStructure pflag;
23: PetscReal dp;
24: Scalar lambda,alpha0,alpha1,btop,bbot,bbotold,tmp,zero = 0.0,mone = -1.0;
25: Vec X,B,R,Pm1,P,Pp1,Sm1,S,Qm1,Q,Qp1,T,Tmp;
26: Mat Amat,Pmat;
27: PetscTruth diagonalscale;
30: ierr = PCDiagonalScale(ksp->B,&diagonalscale);
31: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
33: pres = ksp->use_pres;
34: maxit = ksp->max_it;
35: X = ksp->vec_sol;
36: B = ksp->vec_rhs;
37: R = ksp->work[0];
38: Pm1 = ksp->work[1];
39: P = ksp->work[2];
40: Pp1 = ksp->work[3];
41: Qm1 = ksp->work[4];
42: Q = ksp->work[5];
43: Qp1 = T = ksp->work[6];
44: Sm1 = ksp->work[7];
45: S = ksp->work[8];
47: ierr = PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
48: bbotold = 1.0; /* a hack */
49: if (!ksp->guess_zero) {
50: KSP_MatMult(ksp,Amat,X,R); /* r <- b - Ax */
51: VecAYPX(&mone,B,R);
52: } else {
53: VecCopy(B,R); /* r <- b (x is 0) */
54: }
55: VecSet(&zero,Pm1); /* pm1 <- 0 */
56: VecSet(&zero,Sm1); /* sm1 <- 0 */
57: VecSet(&zero,Qm1); /* Qm1 <- 0 */
58: KSP_PCApply(ksp,ksp->B,R,P); /* p <- Br */
59: if (!ksp->avoidnorms) {
60: if (pres) {
61: VecNorm(P,NORM_2,&dp);/* dp <- z'*z */
62: } else {
63: VecNorm(R,NORM_2,&dp);/* dp <- r'*r */
64: }
65: }
66: PetscObjectTakeAccess(ksp);
67: ksp->its = 0;
68: ksp->rnorm = dp;
69: PetscObjectGrantAccess(ksp);
70: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
71: if (ksp->reason) {*its = 0; return(0);}
72: KSPLogResidualHistory(ksp,dp);
73: KSPMonitor(ksp,0,dp);
74: KSP_MatMult(ksp,Amat,P,Q); /* q <- A p */
76: for (i=0; i<maxit; i++) {
78: ierr = KSP_PCApply(ksp,ksp->B,Q,S); /* s <- Bq */
79: ierr = VecDot(R,S,&btop); /* */
80: ierr = VecDot(Q,S,&bbot); /* lambda = */
81: lambda = btop/bbot;
82: ierr = VecAXPY(&lambda,P,X); /* x <- x + lambda p */
83: tmp = -lambda;
84: ierr = VecAXPY(&tmp,Q,R); /* r <- r - lambda q */
85: if (!ksp->avoidnorms) {
86: ierr = VecNorm(R,NORM_2,&dp); /* dp <- r'*r */
87: } else { dp = 0.0; }
88: PetscObjectTakeAccess(ksp);
89: ksp->its++;
90: ksp->rnorm = dp;
91: PetscObjectGrantAccess(ksp);
92: KSPLogResidualHistory(ksp,dp);
93: KSPMonitor(ksp,i+1,dp);
94: ierr = (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
95: if (ksp->reason) break;
96: ierr = KSP_MatMult(ksp,Amat,S,T); /* T <- As */
97: ierr = VecDot(T,S,&btop);
98: alpha0 = btop/bbot;
99: ierr = VecDot(T,Sm1,&btop);
100: alpha1 = btop/bbotold;
102: tmp = -alpha0; VecWAXPY(&tmp,P,S,Pp1);
103: tmp = -alpha1; VecAXPY(&tmp,Pm1,Pp1);
104: /* MM(Pp1,Qp1); use 3 term recurrence relation instead */
105: tmp = -alpha0; VecAXPY(&tmp,Q,Qp1);
106: tmp = -alpha1; VecAXPY(&tmp,Qm1,Qp1);
107: /* scale the search direction !! Not mentioned in any reference */
108: VecNorm(Pp1,NORM_2,&dp);
109: tmp = 1.0/dp; VecScale(&tmp,Pp1);
110: VecScale(&tmp,Qp1);
111: /* rotate work vectors */
112: Tmp = Sm1; Sm1 = S; S = Tmp;
113: Tmp = Pm1; Pm1 = P; P = Pp1; Pp1 = Tmp;
114: Tmp = Qm1; Qm1 = Q; Q = Qp1; Qp1 = T = Tmp;
115: bbotold = bbot;
116: }
117: if (i == maxit) {ksp->reason = KSP_DIVERGED_ITS;}
118: *its = ksp->its;
119: return(0);
120: }
122: EXTERN_C_BEGIN
123: int KSPCreate_CR(KSP ksp)
124: {
126: ksp->pc_side = PC_LEFT;
127: ksp->calc_res = PETSC_TRUE;
128: ksp->ops->setup = KSPSetUp_CR;
129: ksp->ops->solve = KSPSolve_CR;
130: ksp->ops->destroy = KSPDefaultDestroy;
131: ksp->ops->buildsolution = KSPDefaultBuildSolution;
132: ksp->ops->buildresidual = KSPDefaultBuildResidual;
133: ksp->ops->setfromoptions = 0;
134: ksp->ops->view = 0;
135: return(0);
136: }
137: EXTERN_C_END