xref: /petsc/src/ksp/ksp/impls/symmlq/symmlq.c (revision 445ca0904fedcd9c19e8b5d649a2e70cfb1cd372)
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 = 0.0,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   if (ksp->normtype != KSP_NORM_NONE) {
88     ierr  = VecNorm(Z,NORM_2,&np);CHKERRQ(ierr);     /*   np <- ||z||        */
89     KSPCheckNorm(ksp,np);
90   }
91   ierr = KSPLogResidualHistory(ksp,np);CHKERRQ(ierr);
92   ierr       = KSPMonitor(ksp,0,np);CHKERRQ(ierr);
93   ksp->rnorm = np;
94   ierr       = (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); /* test for convergence */
95   if (ksp->reason) PetscFunctionReturn(0);
96 
97   i = 0; ceta = 0.;
98   do {
99     ksp->its = i+1;
100 
101     /*    Update    */
102     if (ksp->its > 1) {
103       ierr = VecCopy(V,VOLD);CHKERRQ(ierr);  /* v_old <- v; */
104       ierr = VecCopy(U,UOLD);CHKERRQ(ierr);  /* u_old <- u; */
105 
106       ierr = VecCopy(R,V);CHKERRQ(ierr);
107       ierr = VecScale(V,1.0/beta);CHKERRQ(ierr); /* v <- ibeta*r; */
108       ierr = VecCopy(Z,U);CHKERRQ(ierr);
109       ierr = VecScale(U,1.0/beta);CHKERRQ(ierr); /* u <- ibeta*z; */
110 
111       ierr = VecCopy(Wbar,W);CHKERRQ(ierr);
112       ierr = VecScale(W,c);CHKERRQ(ierr);
113       ierr = VecAXPY(W,s,U);CHKERRQ(ierr);   /* w  <- c*w_bar + s*u;    (w_k) */
114       ierr = VecScale(Wbar,-s);CHKERRQ(ierr);
115       ierr = VecAXPY(Wbar,c,U);CHKERRQ(ierr); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
116       ierr = VecAXPY(X,ceta,W);CHKERRQ(ierr); /* x <- x + ceta * w;       (xL_k)  */
117 
118       ceta_oold = ceta_old;
119       ceta_old  = ceta;
120     }
121 
122     /*   Lanczos  */
123     ierr = KSP_MatMult(ksp,Amat,U,R);CHKERRQ(ierr);   /*  r     <- Amat*u; */
124     ierr = VecDot(U,R,&alpha);CHKERRQ(ierr);          /*  alpha <- u'*r;   */
125     ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /*      z <- B*r;    */
126 
127     ierr    = VecAXPY(R,-alpha,V);CHKERRQ(ierr);   /*  r <- r - alpha* v;  */
128     ierr    = VecAXPY(Z,-alpha,U);CHKERRQ(ierr);   /*  z <- z - alpha* u;  */
129     ierr    = VecAXPY(R,-beta,VOLD);CHKERRQ(ierr); /*  r <- r - beta * v_old; */
130     ierr    = VecAXPY(Z,-beta,UOLD);CHKERRQ(ierr); /*  z <- z - beta * u_old; */
131     betaold = beta;                                /* beta_k                  */
132     ierr    = VecDot(R,Z,&dp);CHKERRQ(ierr);       /* dp <- r'*z;             */
133     KSPCheckDot(ksp,dp);
134     if (PetscAbsScalar(dp) < symmlq->haptol) {
135       ierr = PetscInfo2(ksp,"Detected happy breakdown %g tolerance %g\n",(double)PetscAbsScalar(dp),(double)symmlq->haptol);CHKERRQ(ierr);
136       dp   = 0.0;
137     }
138 
139 #if !defined(PETSC_USE_COMPLEX)
140     if (dp < 0.0) {
141       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
142       break;
143     }
144 #endif
145     beta = PetscSqrtScalar(dp);                    /*  beta = sqrt(dp); */
146 
147     /*    QR factorization    */
148     coold = cold; cold = c; soold = sold; sold = s;
149     rho0  = cold * alpha - coold * sold * betaold;   /* gamma_bar */
150     rho1  = PetscSqrtScalar(rho0*rho0 + beta*beta);  /* gamma     */
151     rho2  = sold * alpha + coold * cold * betaold;   /* delta     */
152     rho3  = soold * betaold;                         /* epsilon   */
153 
154     /* Givens rotation: [c -s; s c] (different from the Reference!) */
155     c = rho0 / rho1; s = beta / rho1;
156 
157     if (ksp->its==1) ceta = beta1/rho1;
158     else ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
159 
160     s_prod = s_prod*PetscAbsScalar(s);
161     if (c == 0.0) np = s_prod*1.e16;
162     else np = s_prod/PetscAbsScalar(c);       /* residual norm for xc_k (CGNORM) */
163 
164     if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = np;
165     else ksp->rnorm = 0.0;
166     ierr = KSPLogResidualHistory(ksp,ksp->rnorm);CHKERRQ(ierr);
167     ierr = KSPMonitor(ksp,i+1,ksp->rnorm);CHKERRQ(ierr);
168     ierr = (*ksp->converged)(ksp,i+1,ksp->rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); /* test for convergence */
169     if (ksp->reason) break;
170     i++;
171   } while (i<ksp->max_it);
172 
173   /* move to the CG point: xc_(k+1) */
174   if (c == 0.0) ceta_bar = ceta*1.e15;
175   else ceta_bar = ceta/c;
176 
177   ierr = VecAXPY(X,ceta_bar,Wbar);CHKERRQ(ierr); /* x <- x + ceta_bar*w_bar */
178 
179   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
180   PetscFunctionReturn(0);
181 }
182 
183 /*MC
184      KSPSYMMLQ -  This code implements the SYMMLQ method.
185 
186    Options Database Keys:
187 .   see KSPSolve()
188 
189    Level: beginner
190 
191    Notes:
192     The operator and the preconditioner must be symmetric for this method. The
193           preconditioner must be POSITIVE-DEFINITE.
194 
195           Supports only left preconditioning.
196 
197    Reference: Paige & Saunders, 1975.
198 
199 .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
200 M*/
201 PETSC_EXTERN PetscErrorCode KSPCreate_SYMMLQ(KSP ksp)
202 {
203   KSP_SYMMLQ     *symmlq;
204   PetscErrorCode ierr;
205 
206   PetscFunctionBegin;
207   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);CHKERRQ(ierr);
208   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);CHKERRQ(ierr);
209 
210   ierr           = PetscNewLog(ksp,&symmlq);CHKERRQ(ierr);
211   symmlq->haptol = 1.e-18;
212   ksp->data      = (void*)symmlq;
213 
214   /*
215        Sets the functions that are associated with this data structure
216        (in C++ this is the same as defining virtual functions)
217   */
218   ksp->ops->setup          = KSPSetUp_SYMMLQ;
219   ksp->ops->solve          = KSPSolve_SYMMLQ;
220   ksp->ops->destroy        = KSPDestroyDefault;
221   ksp->ops->setfromoptions = NULL;
222   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
223   ksp->ops->buildresidual  = KSPBuildResidualDefault;
224   PetscFunctionReturn(0);
225 }
226