Actual source code: cr.c
petsc-3.14.2 2020-12-03
2: #include <petsc/private/kspimpl.h>
4: static PetscErrorCode KSPSetUp_CR(KSP ksp)
5: {
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: KSPSetWorkVecs(ksp,6);
12: return(0);
13: }
15: static PetscErrorCode KSPSolve_CR(KSP ksp)
16: {
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;
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];
35: /* R is the true residual norm, RT is the preconditioned residual norm */
36: PCGetOperators(ksp->pc,&Amat,&Pmat);
37: if (!ksp->guess_zero) {
38: KSP_MatMult(ksp,Amat,X,R); /* R <- A*X */
39: VecAYPX(R,-1.0,B); /* R <- B-R == B-A*X */
40: } else {
41: VecCopy(B,R); /* R <- B (X is 0) */
42: }
43: KSP_PCApply(ksp,R,P); /* P <- B*R */
44: KSP_MatMult(ksp,Amat,P,AP); /* AP <- A*P */
45: VecCopy(P,RT); /* RT <- P */
46: VecCopy(AP,ART); /* ART <- AP */
47: VecDotBegin(RT,ART,&btop); /* (RT,ART) */
49: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
50: VecNormBegin(RT,NORM_2,&dp); /* dp <- RT'*RT */
51: VecDotEnd (RT,ART,&btop); /* (RT,ART) */
52: VecNormEnd (RT,NORM_2,&dp); /* 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: VecDotEnd (RT,ART,&btop); /* (RT,ART) */
57: } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
58: VecNormBegin(R,NORM_2,&dp); /* dp <- R'*R */
59: VecDotEnd (RT,ART,&btop); /* (RT,ART) */
60: VecNormEnd (R,NORM_2,&dp); /* dp <- RT'*RT */
61: KSPCheckNorm(ksp,dp);
62: } else if (ksp->normtype == KSP_NORM_NATURAL) {
63: VecDotEnd (RT,ART,&btop); /* (RT,ART) */
64: dp = PetscSqrtReal(PetscAbsScalar(btop)); /* dp = sqrt(R,AR) */
65: } else SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"KSPNormType of %d not supported",(int)ksp->normtype);
66: if (PetscAbsScalar(btop) < 0.0) {
67: ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
68: PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");
69: return(0);
70: }
72: ksp->its = 0;
73: KSPMonitor(ksp,0,dp);
74: PetscObjectSAWsTakeAccess((PetscObject)ksp);
75: ksp->rnorm = dp;
76: PetscObjectSAWsGrantAccess((PetscObject)ksp);
77: KSPLogResidualHistory(ksp,dp);
78: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
79: if (ksp->reason) return(0);
81: i = 0;
82: do {
83: KSP_PCApply(ksp,AP,Q); /* Q <- B* AP */
85: VecDot(AP,Q,&apq);
86: KSPCheckDot(ksp,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: KSPCheckNorm(ksp,dp);
105: } else if (ksp->normtype == KSP_NORM_NATURAL) {
106: VecDotEnd(RT,ART,&btop);
107: dp = PetscSqrtReal(PetscAbsScalar(btop)); /* dp = sqrt(R,AR) */
108: } else if (ksp->normtype == KSP_NORM_NONE) {
109: VecDotEnd(RT,ART,&btop);
110: dp = 0.0; /* meaningless value that is passed to monitor and convergence test */
111: } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
112: VecAXPY(R,ai,AP); /* R <- R - ai*AP */
113: VecNormBegin(R,NORM_2,&dp); /* dp <- R'*R */
114: VecDotEnd (RT,ART,&btop);
115: VecNormEnd (R,NORM_2,&dp); /* dp <- R'*R */
116: KSPCheckNorm(ksp,dp);
117: } else SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"KSPNormType of %d not supported",(int)ksp->normtype);
118: if (PetscAbsScalar(btop) < 0.0) {
119: ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
120: PetscInfo(ksp,"diverging due to indefinite or negative definite PC\n");
121: break;
122: }
124: PetscObjectSAWsTakeAccess((PetscObject)ksp);
125: ksp->its++;
126: ksp->rnorm = dp;
127: PetscObjectSAWsGrantAccess((PetscObject)ksp);
129: KSPLogResidualHistory(ksp,dp);
130: KSPMonitor(ksp,i+1,dp);
131: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
132: if (ksp->reason) break;
134: bi = btop/bbot;
135: VecAYPX(P,bi,RT); /* P <- RT + Bi P */
136: VecAYPX(AP,bi,ART); /* AP <- ART + Bi AP */
137: i++;
138: } while (i<ksp->max_it);
139: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
140: return(0);
141: }
144: /*MC
145: KSPCR - This code implements the (preconditioned) conjugate residuals method
147: Options Database Keys:
148: . see KSPSolve()
150: Level: beginner
152: Notes:
153: The operator and the preconditioner must be symmetric for this method. The
154: preconditioner must be POSITIVE-DEFINITE and the operator POSITIVE-SEMIDEFINITE.
155: Support only for left preconditioning.
157: References:
158: . 1. - Magnus R. Hestenes and Eduard Stiefel, Methods of Conjugate Gradients for Solving Linear Systems,
159: Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379
161: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPCG
162: M*/
163: PETSC_EXTERN PetscErrorCode KSPCreate_CR(KSP ksp)
164: {
168: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
169: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
170: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
171: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
173: ksp->ops->setup = KSPSetUp_CR;
174: ksp->ops->solve = KSPSolve_CR;
175: ksp->ops->destroy = KSPDestroyDefault;
176: ksp->ops->buildsolution = KSPBuildSolutionDefault;
177: ksp->ops->buildresidual = KSPBuildResidualDefault;
178: ksp->ops->setfromoptions = NULL;
179: ksp->ops->view = NULL;
180: return(0);
181: }