xref: /petsc/src/ksp/ksp/impls/cr/cr.c (revision 37eeb8152ec6a2cf24186d3591c2c5de5dfd8fa5)
1 
2 #include <petsc/private/kspimpl.h>
3 
4 static PetscErrorCode KSPSetUp_CR(KSP ksp)
5 {
6   PetscErrorCode ierr;
7 
8   PetscFunctionBegin;
9   if (ksp->pc_side == PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"no right preconditioning for KSPCR");
10   else if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"no symmetric preconditioning for KSPCR");
11   ierr = KSPSetWorkVecs(ksp,6);CHKERRQ(ierr);
12   PetscFunctionReturn(0);
13 }
14 
15 static PetscErrorCode  KSPSolve_CR(KSP ksp)
16 {
17   PetscErrorCode ierr;
18   PetscInt       i = 0;
19   PetscReal      dp;
20   PetscScalar    ai, bi;
21   PetscScalar    apq,btop, bbot;
22   Vec            X,B,R,RT,P,AP,ART,Q;
23   Mat            Amat, Pmat;
24 
25   PetscFunctionBegin;
26   X   = ksp->vec_sol;
27   B   = ksp->vec_rhs;
28   R   = ksp->work[0];
29   RT  = ksp->work[1];
30   P   = ksp->work[2];
31   AP  = ksp->work[3];
32   ART = ksp->work[4];
33   Q   = ksp->work[5];
34 
35   /* R is the true residual norm, RT is the preconditioned residual norm */
36   ierr = PCGetOperators(ksp->pc,&Amat,&Pmat);CHKERRQ(ierr);
37   if (!ksp->guess_zero) {
38     ierr = KSP_MatMult(ksp,Amat,X,R);CHKERRQ(ierr);     /*   R <- A*X           */
39     ierr = VecAYPX(R,-1.0,B);CHKERRQ(ierr);            /*   R <- B-R == B-A*X  */
40   } else {
41     ierr = VecCopy(B,R);CHKERRQ(ierr);                  /*   R <- B (X is 0)    */
42   }
43   ierr = KSP_PCApply(ksp,R,P);CHKERRQ(ierr);     /*   P   <- B*R         */
44   ierr = KSP_MatMult(ksp,Amat,P,AP);CHKERRQ(ierr);      /*   AP  <- A*P         */
45   ierr = VecCopy(P,RT);CHKERRQ(ierr);                   /*   RT  <- P           */
46   ierr = VecCopy(AP,ART);CHKERRQ(ierr);                 /*   ART <- AP          */
47   ierr = VecDotBegin(RT,ART,&btop);CHKERRQ(ierr);          /*   (RT,ART)           */
48 
49   if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
50     ierr = VecNormBegin(RT,NORM_2,&dp);CHKERRQ(ierr);        /*   dp <- RT'*RT       */
51     ierr = VecDotEnd   (RT,ART,&btop);CHKERRQ(ierr);           /*   (RT,ART)           */
52     ierr = VecNormEnd  (RT,NORM_2,&dp);CHKERRQ(ierr);        /*   dp <- RT'*RT       */
53     KSPCheckNorm(ksp,dp);
54   } else if (ksp->normtype == KSP_NORM_NONE) {
55       dp   = 0.0; /* meaningless value that is passed to monitor and convergence test */
56   } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
57     ierr = VecNormBegin(R,NORM_2,&dp);CHKERRQ(ierr);         /*   dp <- R'*R         */
58     ierr = VecDotEnd   (RT,ART,&btop);CHKERRQ(ierr);          /*   (RT,ART)           */
59     ierr = VecNormEnd  (R,NORM_2,&dp);CHKERRQ(ierr);        /*   dp <- RT'*RT       */
60     KSPCheckNorm(ksp,dp);
61   } else if (ksp->normtype == KSP_NORM_NATURAL) {
62     ierr = VecDotEnd   (RT,ART,&btop);CHKERRQ(ierr);           /*   (RT,ART)           */
63     dp   = PetscSqrtReal(PetscAbsScalar(btop));                  /* dp = sqrt(R,AR)      */
64   } else SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"KSPNormType of %d not supported",(int)ksp->normtype);
65   if (PetscAbsScalar(btop) < 0.0) {
66     ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
67     ierr        = PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");CHKERRQ(ierr);
68     PetscFunctionReturn(0);
69   }
70 
71   ksp->its   = 0;
72   ierr       = KSPMonitor(ksp,0,dp);CHKERRQ(ierr);
73   ierr       = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr);
74   ksp->rnorm = dp;
75   ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr);
76   ierr = KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr);
77   ierr = (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);
78   if (ksp->reason) PetscFunctionReturn(0);
79 
80   i = 0;
81   do {
82     ierr = KSP_PCApply(ksp,AP,Q);CHKERRQ(ierr);  /*   Q <- B* AP          */
83 
84     ierr = VecDot(AP,Q,&apq);CHKERRQ(ierr);
85     KSPCheckDot(ksp,apq);
86     if (PetscRealPart(apq) <= 0.0) {
87       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
88       ierr        = PetscInfo(ksp,"KSPSolve_CR:diverging due to indefinite or negative definite PC\n");CHKERRQ(ierr);
89       break;
90     }
91     ai = btop/apq;                                      /* ai = (RT,ART)/(AP,Q)  */
92 
93     ierr = VecAXPY(X,ai,P);CHKERRQ(ierr);              /*   X   <- X + ai*P     */
94     ierr = VecAXPY(RT,-ai,Q);CHKERRQ(ierr);             /*   RT  <- RT - ai*Q    */
95     ierr = KSP_MatMult(ksp,Amat,RT,ART);CHKERRQ(ierr);  /*   ART <-   A*RT       */
96     bbot = btop;
97     ierr = VecDotBegin(RT,ART,&btop);CHKERRQ(ierr);
98 
99     if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
100       ierr = VecNormBegin(RT,NORM_2,&dp);CHKERRQ(ierr);      /*   dp <- || RT ||      */
101       ierr = VecDotEnd   (RT,ART,&btop);CHKERRQ(ierr);
102       ierr = VecNormEnd  (RT,NORM_2,&dp);CHKERRQ(ierr);      /*   dp <- || RT ||      */
103       KSPCheckNorm(ksp,dp);
104     } else if (ksp->normtype == KSP_NORM_NATURAL) {
105       ierr = VecDotEnd(RT,ART,&btop);CHKERRQ(ierr);
106       dp   = PetscSqrtReal(PetscAbsScalar(btop));                /* dp = sqrt(R,AR)       */
107     } else if (ksp->normtype == KSP_NORM_NONE) {
108       ierr = VecDotEnd(RT,ART,&btop);CHKERRQ(ierr);
109       dp   = 0.0; /* meaningless value that is passed to monitor and convergence test */
110     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
111       ierr = VecAXPY(R,ai,AP);CHKERRQ(ierr);           /*   R   <- R - ai*AP    */
112       ierr = VecNormBegin(R,NORM_2,&dp);CHKERRQ(ierr);       /*   dp <- R'*R          */
113       ierr = VecDotEnd   (RT,ART,&btop);CHKERRQ(ierr);
114       ierr = VecNormEnd  (R,NORM_2,&dp);CHKERRQ(ierr);       /*   dp <- R'*R          */
115       KSPCheckNorm(ksp,dp);
116     } else SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"KSPNormType of %d not supported",(int)ksp->normtype);
117     if (PetscAbsScalar(btop) < 0.0) {
118       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
119       ierr        = PetscInfo(ksp,"diverging due to indefinite or negative definite PC\n");CHKERRQ(ierr);
120       break;
121     }
122 
123     ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr);
124     ksp->its++;
125     ksp->rnorm = dp;
126     ierr       = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr);
127 
128     ierr = KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr);
129     ierr = KSPMonitor(ksp,i+1,dp);CHKERRQ(ierr);
130     ierr = (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);
131     if (ksp->reason) break;
132 
133     bi   = btop/bbot;
134     ierr = VecAYPX(P,bi,RT);CHKERRQ(ierr);              /*   P <- RT + Bi P     */
135     ierr = VecAYPX(AP,bi,ART);CHKERRQ(ierr);            /*   AP <- ART + Bi AP  */
136     i++;
137   } while (i<ksp->max_it);
138   if (i >= ksp->max_it) ksp->reason =  KSP_DIVERGED_ITS;
139   PetscFunctionReturn(0);
140 }
141 
142 
143 /*MC
144      KSPCR - This code implements the (preconditioned) conjugate residuals method
145 
146    Options Database Keys:
147 .   see KSPSolve()
148 
149    Level: beginner
150 
151    Notes:
152     The operator and the preconditioner must be symmetric for this method. The
153           preconditioner must be POSITIVE-DEFINITE and the operator POSITIVE-SEMIDEFINITE.
154           Support only for left preconditioning.
155 
156    References:
157 .   1. - Magnus R. Hestenes and Eduard Stiefel, Methods of Conjugate Gradients for Solving Linear Systems,
158    Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379
159 
160 .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPCG
161 M*/
162 PETSC_EXTERN PetscErrorCode KSPCreate_CR(KSP ksp)
163 {
164   PetscErrorCode ierr;
165 
166   PetscFunctionBegin;
167   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);CHKERRQ(ierr);
168   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);CHKERRQ(ierr);
169   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);CHKERRQ(ierr);
170   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);CHKERRQ(ierr);
171 
172   ksp->ops->setup          = KSPSetUp_CR;
173   ksp->ops->solve          = KSPSolve_CR;
174   ksp->ops->destroy        = KSPDestroyDefault;
175   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
176   ksp->ops->buildresidual  = KSPBuildResidualDefault;
177   ksp->ops->setfromoptions = 0;
178   ksp->ops->view           = 0;
179   PetscFunctionReturn(0);
180 }
181