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