xref: /petsc/src/ksp/ksp/impls/symmlq/symmlq.c (revision c0e053deb26eb7702ba63f4e0953a2eb3ba4d64c)
1 
2 #include <petsc/private/kspimpl.h>
3 
4 typedef struct {
5   PetscReal haptol;
6 } KSP_SYMMLQ;
7 
8 PetscErrorCode KSPSetUp_SYMMLQ(KSP ksp)
9 {
10   PetscErrorCode ierr;
11 
12   PetscFunctionBegin;
13   ierr = KSPSetWorkVecs(ksp,9);CHKERRQ(ierr);
14   PetscFunctionReturn(0);
15 }
16 
17 PetscErrorCode  KSPSolve_SYMMLQ(KSP ksp)
18 {
19   PetscErrorCode ierr;
20   PetscInt       i;
21   PetscScalar    alpha,beta,ibeta,betaold,beta1,ceta = 0,ceta_oold = 0.0, ceta_old = 0.0,ceta_bar;
22   PetscScalar    c  = 1.0,cold=1.0,s=0.0,sold=0.0,coold,soold,rho0,rho1,rho2,rho3;
23   PetscScalar    dp = 0.0;
24   PetscReal      np,s_prod;
25   Vec            X,B,R,Z,U,V,W,UOLD,VOLD,Wbar;
26   Mat            Amat,Pmat;
27   KSP_SYMMLQ     *symmlq = (KSP_SYMMLQ*)ksp->data;
28   PetscBool      diagonalscale;
29 
30   PetscFunctionBegin;
31   ierr = PCGetDiagonalScale(ksp->pc,&diagonalscale);CHKERRQ(ierr);
32   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
33 
34   X    = ksp->vec_sol;
35   B    = ksp->vec_rhs;
36   R    = ksp->work[0];
37   Z    = ksp->work[1];
38   U    = ksp->work[2];
39   V    = ksp->work[3];
40   W    = ksp->work[4];
41   UOLD = ksp->work[5];
42   VOLD = ksp->work[6];
43   Wbar = ksp->work[7];
44 
45   ierr = PCGetOperators(ksp->pc,&Amat,&Pmat);CHKERRQ(ierr);
46 
47   ksp->its = 0;
48 
49   ierr = VecSet(UOLD,0.0);CHKERRQ(ierr);           /* u_old <- zeros;  */
50   ierr = VecCopy(UOLD,VOLD);CHKERRQ(ierr);          /* v_old <- u_old;  */
51   ierr = VecCopy(UOLD,W);CHKERRQ(ierr);             /* w     <- u_old;  */
52   ierr = VecCopy(UOLD,Wbar);CHKERRQ(ierr);          /* w_bar <- u_old;  */
53   if (!ksp->guess_zero) {
54     ierr = KSP_MatMult(ksp,Amat,X,R);CHKERRQ(ierr); /*     r <- b - A*x */
55     ierr = VecAYPX(R,-1.0,B);CHKERRQ(ierr);
56   } else {
57     ierr = VecCopy(B,R);CHKERRQ(ierr);              /*     r <- b (x is 0) */
58   }
59 
60   ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z  <- B*r       */
61   ierr = VecDot(R,Z,&dp);CHKERRQ(ierr);             /* dp = r'*z;      */
62   KSPCheckDot(ksp,dp);
63   if (PetscAbsScalar(dp) < symmlq->haptol) {
64     ierr        = PetscInfo2(ksp,"Detected happy breakdown %g tolerance %g\n",(double)PetscAbsScalar(dp),(double)symmlq->haptol);CHKERRQ(ierr);
65     ksp->rnorm  = 0.0;  /* what should we really put here? */
66     ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;  /* bugfix proposed by Lourens (lourens.vanzanen@shell.com) */
67     PetscFunctionReturn(0);
68   }
69 
70 #if !defined(PETSC_USE_COMPLEX)
71   if (dp < 0.0) {
72     ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
73     PetscFunctionReturn(0);
74   }
75 #endif
76   dp     = PetscSqrtScalar(dp);
77   beta   = dp;                         /*  beta <- sqrt(r'*z)  */
78   beta1  = beta;
79   s_prod = PetscAbsScalar(beta1);
80 
81   ierr  = VecCopy(R,V);CHKERRQ(ierr); /* v <- r; */
82   ierr  = VecCopy(Z,U);CHKERRQ(ierr); /* u <- z; */
83   ibeta = 1.0 / beta;
84   ierr  = VecScale(V,ibeta);CHKERRQ(ierr);    /* v <- ibeta*v; */
85   ierr  = VecScale(U,ibeta);CHKERRQ(ierr);    /* u <- ibeta*u; */
86   ierr  = VecCopy(U,Wbar);CHKERRQ(ierr);       /* w_bar <- u;   */
87   ierr  = VecNorm(Z,NORM_2,&np);CHKERRQ(ierr);     /*   np <- ||z||        */
88   KSPCheckNorm(ksp,np);
89   ierr = KSPLogResidualHistory(ksp,np);CHKERRQ(ierr);
90   ierr       = KSPMonitor(ksp,0,np);CHKERRQ(ierr);
91   ksp->rnorm = np;
92   ierr       = (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); /* test for convergence */
93   if (ksp->reason) PetscFunctionReturn(0);
94 
95   i = 0; ceta = 0.;
96   do {
97     ksp->its = i+1;
98 
99     /*    Update    */
100     if (ksp->its > 1) {
101       ierr = VecCopy(V,VOLD);CHKERRQ(ierr);  /* v_old <- v; */
102       ierr = VecCopy(U,UOLD);CHKERRQ(ierr);  /* u_old <- u; */
103 
104       ierr = VecCopy(R,V);CHKERRQ(ierr);
105       ierr = VecScale(V,1.0/beta);CHKERRQ(ierr); /* v <- ibeta*r; */
106       ierr = VecCopy(Z,U);CHKERRQ(ierr);
107       ierr = VecScale(U,1.0/beta);CHKERRQ(ierr); /* u <- ibeta*z; */
108 
109       ierr = VecCopy(Wbar,W);CHKERRQ(ierr);
110       ierr = VecScale(W,c);CHKERRQ(ierr);
111       ierr = VecAXPY(W,s,U);CHKERRQ(ierr);   /* w  <- c*w_bar + s*u;    (w_k) */
112       ierr = VecScale(Wbar,-s);CHKERRQ(ierr);
113       ierr = VecAXPY(Wbar,c,U);CHKERRQ(ierr); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
114       ierr = VecAXPY(X,ceta,W);CHKERRQ(ierr); /* x <- x + ceta * w;       (xL_k)  */
115 
116       ceta_oold = ceta_old;
117       ceta_old  = ceta;
118     }
119 
120     /*   Lanczos  */
121     ierr = KSP_MatMult(ksp,Amat,U,R);CHKERRQ(ierr);   /*  r     <- Amat*u; */
122     ierr = VecDot(U,R,&alpha);CHKERRQ(ierr);          /*  alpha <- u'*r;   */
123     ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /*      z <- B*r;    */
124 
125     ierr    = VecAXPY(R,-alpha,V);CHKERRQ(ierr);   /*  r <- r - alpha* v;  */
126     ierr    = VecAXPY(Z,-alpha,U);CHKERRQ(ierr);   /*  z <- z - alpha* u;  */
127     ierr    = VecAXPY(R,-beta,VOLD);CHKERRQ(ierr); /*  r <- r - beta * v_old; */
128     ierr    = VecAXPY(Z,-beta,UOLD);CHKERRQ(ierr); /*  z <- z - beta * u_old; */
129     betaold = beta;                                /* beta_k                  */
130     ierr    = VecDot(R,Z,&dp);CHKERRQ(ierr);       /* dp <- r'*z;             */
131     KSPCheckDot(ksp,dp);
132     if (PetscAbsScalar(dp) < symmlq->haptol) {
133       ierr = PetscInfo2(ksp,"Detected happy breakdown %g tolerance %g\n",(double)PetscAbsScalar(dp),(double)symmlq->haptol);CHKERRQ(ierr);
134       dp   = 0.0;
135     }
136 
137 #if !defined(PETSC_USE_COMPLEX)
138     if (dp < 0.0) {
139       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
140       break;
141     }
142 #endif
143     beta = PetscSqrtScalar(dp);                    /*  beta = sqrt(dp); */
144 
145     /*    QR factorization    */
146     coold = cold; cold = c; soold = sold; sold = s;
147     rho0  = cold * alpha - coold * sold * betaold;   /* gamma_bar */
148     rho1  = PetscSqrtScalar(rho0*rho0 + beta*beta);  /* gamma     */
149     rho2  = sold * alpha + coold * cold * betaold;   /* delta     */
150     rho3  = soold * betaold;                         /* epsilon   */
151 
152     /* Givens rotation: [c -s; s c] (different from the Reference!) */
153     c = rho0 / rho1; s = beta / rho1;
154 
155     if (ksp->its==1) ceta = beta1/rho1;
156     else ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
157 
158     s_prod = s_prod*PetscAbsScalar(s);
159     if (c == 0.0) np = s_prod*1.e16;
160     else np = s_prod/PetscAbsScalar(c);       /* residual norm for xc_k (CGNORM) */
161 
162     ksp->rnorm = np;
163     ierr = KSPLogResidualHistory(ksp,np);CHKERRQ(ierr);
164     ierr = KSPMonitor(ksp,i+1,np);CHKERRQ(ierr);
165     ierr = (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); /* test for convergence */
166     if (ksp->reason) break;
167     i++;
168   } while (i<ksp->max_it);
169 
170   /* move to the CG point: xc_(k+1) */
171   if (c == 0.0) ceta_bar = ceta*1.e15;
172   else ceta_bar = ceta/c;
173 
174   ierr = VecAXPY(X,ceta_bar,Wbar);CHKERRQ(ierr); /* x <- x + ceta_bar*w_bar */
175 
176   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
177   PetscFunctionReturn(0);
178 }
179 
180 /*MC
181      KSPSYMMLQ -  This code implements the SYMMLQ method.
182 
183    Options Database Keys:
184 .   see KSPSolve()
185 
186    Level: beginner
187 
188    Notes:
189     The operator and the preconditioner must be symmetric for this method. The
190           preconditioner must be POSITIVE-DEFINITE.
191 
192           Supports only left preconditioning.
193 
194    Reference: Paige & Saunders, 1975.
195 
196 .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
197 M*/
198 PETSC_EXTERN PetscErrorCode KSPCreate_SYMMLQ(KSP ksp)
199 {
200   KSP_SYMMLQ     *symmlq;
201   PetscErrorCode ierr;
202 
203   PetscFunctionBegin;
204   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);CHKERRQ(ierr);
205 
206   ierr           = PetscNewLog(ksp,&symmlq);CHKERRQ(ierr);
207   symmlq->haptol = 1.e-18;
208   ksp->data      = (void*)symmlq;
209 
210   /*
211        Sets the functions that are associated with this data structure
212        (in C++ this is the same as defining virtual functions)
213   */
214   ksp->ops->setup          = KSPSetUp_SYMMLQ;
215   ksp->ops->solve          = KSPSolve_SYMMLQ;
216   ksp->ops->destroy        = KSPDestroyDefault;
217   ksp->ops->setfromoptions = 0;
218   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
219   ksp->ops->buildresidual  = KSPBuildResidualDefault;
220   PetscFunctionReturn(0);
221 }
222 
223 
224 
225 
226 
227