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