Actual source code: cr.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: #include <petsc/private/kspimpl.h>

  6: static PetscErrorCode KSPSetUp_CR(KSP ksp)
  7: {

 11:   if (ksp->pc_side == PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"no right preconditioning for KSPCR");
 12:   else if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"no symmetric preconditioning for KSPCR");
 13:   KSPSetWorkVecs(ksp,6);
 14:   return(0);
 15: }

 19: static PetscErrorCode  KSPSolve_CR(KSP ksp)
 20: {
 22:   PetscInt       i = 0;
 23:   PetscReal      dp;
 24:   PetscScalar    ai, bi;
 25:   PetscScalar    apq,btop, bbot;
 26:   Vec            X,B,R,RT,P,AP,ART,Q;
 27:   Mat            Amat, Pmat;

 30:   X   = ksp->vec_sol;
 31:   B   = ksp->vec_rhs;
 32:   R   = ksp->work[0];
 33:   RT  = ksp->work[1];
 34:   P   = ksp->work[2];
 35:   AP  = ksp->work[3];
 36:   ART = ksp->work[4];
 37:   Q   = ksp->work[5];

 39:   /* R is the true residual norm, RT is the preconditioned residual norm */
 40:   PCGetOperators(ksp->pc,&Amat,&Pmat);
 41:   if (!ksp->guess_zero) {
 42:     KSP_MatMult(ksp,Amat,X,R);     /*   R <- A*X           */
 43:     VecAYPX(R,-1.0,B);            /*   R <- B-R == B-A*X  */
 44:   } else {
 45:     VecCopy(B,R);                  /*   R <- B (X is 0)    */
 46:   }
 47:   KSP_PCApply(ksp,R,P);     /*   P   <- B*R         */
 48:   KSP_MatMult(ksp,Amat,P,AP);      /*   AP  <- A*P         */
 49:   VecCopy(P,RT);                   /*   RT  <- P           */
 50:   VecCopy(AP,ART);                 /*   ART <- AP          */
 51:   VecDotBegin(RT,ART,&btop);          /*   (RT,ART)           */

 53:   if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
 54:     VecNormBegin(RT,NORM_2,&dp);        /*   dp <- RT'*RT       */
 55:     VecDotEnd   (RT,ART,&btop);           /*   (RT,ART)           */
 56:     VecNormEnd  (RT,NORM_2,&dp);        /*   dp <- RT'*RT       */
 57:   } else if (ksp->normtype == KSP_NORM_NONE) {
 58:       dp   = 0.0; /* meaningless value that is passed to monitor and convergence test */
 59:   } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
 60:     VecNormBegin(R,NORM_2,&dp);         /*   dp <- R'*R         */
 61:     VecDotEnd   (RT,ART,&btop);          /*   (RT,ART)           */
 62:     VecNormEnd  (R,NORM_2,&dp);        /*   dp <- RT'*RT       */
 63:   } else if (ksp->normtype == KSP_NORM_NATURAL) {
 64:     VecDotEnd   (RT,ART,&btop);           /*   (RT,ART)           */
 65:     dp   = PetscSqrtReal(PetscAbsScalar(btop));                  /* dp = sqrt(R,AR)      */
 66:   } else SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"KSPNormType of %d not supported",(int)ksp->normtype);
 67:   if (PetscAbsScalar(btop) < 0.0) {
 68:     ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
 69:     PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");
 70:     return(0);
 71:   }

 73:   ksp->its   = 0;
 74:   KSPMonitor(ksp,0,dp);
 75:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
 76:   ksp->rnorm = dp;
 77:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
 78:   KSPLogResidualHistory(ksp,dp);
 79:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
 80:   if (ksp->reason) return(0);

 82:   i = 0;
 83:   do {
 84:     KSP_PCApply(ksp,AP,Q);  /*   Q <- B* AP          */

 86:     VecDot(AP,Q,&apq);
 87:     if (PetscRealPart(apq) <= 0.0) {
 88:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
 89:       PetscInfo(ksp,"KSPSolve_CR:diverging due to indefinite or negative definite PC\n");
 90:       break;
 91:     }
 92:     ai = btop/apq;                                      /* ai = (RT,ART)/(AP,Q)  */

 94:     VecAXPY(X,ai,P);              /*   X   <- X + ai*P     */
 95:     VecAXPY(RT,-ai,Q);             /*   RT  <- RT - ai*Q    */
 96:     KSP_MatMult(ksp,Amat,RT,ART);  /*   ART <-   A*RT       */
 97:     bbot = btop;
 98:     VecDotBegin(RT,ART,&btop);

100:     if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
101:       VecNormBegin(RT,NORM_2,&dp);      /*   dp <- || RT ||      */
102:       VecDotEnd   (RT,ART,&btop);
103:       VecNormEnd  (RT,NORM_2,&dp);      /*   dp <- || RT ||      */
104:     } else if (ksp->normtype == KSP_NORM_NATURAL) {
105:       VecDotEnd(RT,ART,&btop);
106:       dp   = PetscSqrtReal(PetscAbsScalar(btop));                /* dp = sqrt(R,AR)       */
107:     } else if (ksp->normtype == KSP_NORM_NONE) {
108:       VecDotEnd(RT,ART,&btop);
109:       dp   = 0.0; /* meaningless value that is passed to monitor and convergence test */
110:     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
111:       VecAXPY(R,ai,AP);           /*   R   <- R - ai*AP    */
112:       VecNormBegin(R,NORM_2,&dp);       /*   dp <- R'*R          */
113:       VecDotEnd   (RT,ART,&btop);
114:       VecNormEnd  (R,NORM_2,&dp);       /*   dp <- R'*R          */
115:     } else SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"KSPNormType of %d not supported",(int)ksp->normtype);
116:     if (PetscAbsScalar(btop) < 0.0) {
117:       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
118:       PetscInfo(ksp,"diverging due to indefinite or negative definite PC\n");
119:       break;
120:     }

122:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
123:     ksp->its++;
124:     ksp->rnorm = dp;
125:     PetscObjectSAWsGrantAccess((PetscObject)ksp);

127:     KSPLogResidualHistory(ksp,dp);
128:     KSPMonitor(ksp,i+1,dp);
129:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
130:     if (ksp->reason) break;

132:     bi   = btop/bbot;
133:     VecAYPX(P,bi,RT);              /*   P <- RT + Bi P     */
134:     VecAYPX(AP,bi,ART);            /*   AP <- ART + Bi AP  */
135:     i++;
136:   } while (i<ksp->max_it);
137:   if (i >= ksp->max_it) ksp->reason =  KSP_DIVERGED_ITS;
138:   return(0);
139: }


142: /*MC
143:      KSPCR - This code implements the (preconditioned) conjugate residuals method

145:    Options Database Keys:
146: .   see KSPSolve()

148:    Level: beginner

150:    Notes: The operator and the preconditioner must be symmetric for this method. The
151:           preconditioner must be POSITIVE-DEFINITE and the operator POSITIVE-SEMIDEFINITE.
152:           Support only for left preconditioning.

154:    References:
155: .   1. - Magnus R. Hestenes and Eduard Stiefel, Methods of Conjugate Gradients for Solving Linear Systems, 
156:    Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379

158: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPCG
159: M*/
162: PETSC_EXTERN PetscErrorCode KSPCreate_CR(KSP ksp)
163: {

167:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
168:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
169:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);

171:   ksp->ops->setup          = KSPSetUp_CR;
172:   ksp->ops->solve          = KSPSolve_CR;
173:   ksp->ops->destroy        = KSPDestroyDefault;
174:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
175:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
176:   ksp->ops->setfromoptions = 0;
177:   ksp->ops->view           = 0;
178:   return(0);
179: }