Actual source code: cgs.c
petsc-3.14.3 2021-01-09
2: /*
4: Note that for the complex numbers version, the VecDot() arguments
5: within the code MUST remain in the order given for correct computation
6: of inner products.
7: */
8: #include <petsc/private/kspimpl.h>
10: static PetscErrorCode KSPSetUp_CGS(KSP ksp)
11: {
15: KSPSetWorkVecs(ksp,7);
16: return(0);
17: }
19: static PetscErrorCode KSPSolve_CGS(KSP ksp)
20: {
22: PetscInt i;
23: PetscScalar rho,rhoold,a,s,b;
24: Vec X,B,V,P,R,RP,T,Q,U,AUQ;
25: PetscReal dp = 0.0;
26: PetscBool diagonalscale;
29: /* not sure what residual norm it does use, should use for right preconditioning */
31: PCGetDiagonalScale(ksp->pc,&diagonalscale);
32: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
34: X = ksp->vec_sol;
35: B = ksp->vec_rhs;
36: R = ksp->work[0];
37: RP = ksp->work[1];
38: V = ksp->work[2];
39: T = ksp->work[3];
40: Q = ksp->work[4];
41: P = ksp->work[5];
42: U = ksp->work[6];
43: AUQ = V;
45: /* Compute initial preconditioned residual */
46: KSPInitialResidual(ksp,X,V,T,R,B);
48: /* Test for nothing to do */
49: if (ksp->normtype != KSP_NORM_NONE) {
50: VecNorm(R,NORM_2,&dp);
51: KSPCheckNorm(ksp,dp);
52: if (ksp->normtype == KSP_NORM_NATURAL) dp *= dp;
53: } else dp = 0.0;
55: PetscObjectSAWsTakeAccess((PetscObject)ksp);
56: ksp->its = 0;
57: ksp->rnorm = dp;
58: PetscObjectSAWsGrantAccess((PetscObject)ksp);
59: KSPLogResidualHistory(ksp,dp);
60: KSPMonitor(ksp,0,dp);
61: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
62: if (ksp->reason) return(0);
64: /* Make the initial Rp == R */
65: VecCopy(R,RP);
66: /* added for Fidap */
67: /* Penalize Startup - Isaac Hasbani Trick for CGS
68: Since most initial conditions result in a mostly 0 residual,
69: we change all the 0 values in the vector RP to the maximum.
70: */
71: if (ksp->normtype == KSP_NORM_NATURAL) {
72: PetscReal vr0max;
73: PetscScalar *tmp_RP=NULL;
74: PetscInt numnp =0, *max_pos=NULL;
75: VecMax(RP, max_pos, &vr0max);
76: VecGetArray(RP, &tmp_RP);
77: VecGetLocalSize(RP, &numnp);
78: for (i=0; i<numnp; i++) {
79: if (tmp_RP[i] == 0.0) tmp_RP[i] = vr0max;
80: }
81: VecRestoreArray(RP, &tmp_RP);
82: }
83: /* end of addition for Fidap */
85: /* Set the initial conditions */
86: VecDot(R,RP,&rhoold); /* rhoold = (r,rp) */
87: VecCopy(R,U);
88: VecCopy(R,P);
89: KSP_PCApplyBAorAB(ksp,P,V,T);
91: i = 0;
92: do {
94: VecDot(V,RP,&s); /* s <- (v,rp) */
95: KSPCheckDot(ksp,s);
96: a = rhoold / s; /* a <- rho / s */
97: VecWAXPY(Q,-a,V,U); /* q <- u - a v */
98: VecWAXPY(T,1.0,U,Q); /* t <- u + q */
99: VecAXPY(X,a,T); /* x <- x + a (u + q) */
100: KSP_PCApplyBAorAB(ksp,T,AUQ,U);
101: VecAXPY(R,-a,AUQ); /* r <- r - a K (u + q) */
102: VecDot(R,RP,&rho); /* rho <- (r,rp) */
103: KSPCheckDot(ksp,rho);
104: if (ksp->normtype == KSP_NORM_NATURAL) {
105: dp = PetscAbsScalar(rho);
106: } else if (ksp->normtype != KSP_NORM_NONE) {
107: VecNorm(R,NORM_2,&dp);
108: KSPCheckNorm(ksp,dp);
109: } else dp = 0.0;
111: PetscObjectSAWsTakeAccess((PetscObject)ksp);
112: ksp->its++;
113: ksp->rnorm = dp;
114: PetscObjectSAWsGrantAccess((PetscObject)ksp);
115: KSPLogResidualHistory(ksp,dp);
116: KSPMonitor(ksp,i+1,dp);
117: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
118: if (ksp->reason) break;
120: b = rho / rhoold; /* b <- rho / rhoold */
121: VecWAXPY(U,b,Q,R); /* u <- r + b q */
122: VecAXPY(Q,b,P);
123: VecWAXPY(P,b,Q,U); /* p <- u + b(q + b p) */
124: KSP_PCApplyBAorAB(ksp,P,V,Q); /* v <- K p */
125: rhoold = rho;
126: i++;
127: } while (i<ksp->max_it);
128: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
130: KSPUnwindPreconditioner(ksp,X,T);
131: return(0);
132: }
134: /*MC
135: KSPCGS - This code implements the CGS (Conjugate Gradient Squared) method.
137: Options Database Keys:
138: . see KSPSolve()
140: Level: beginner
142: References:
143: . 1. - Sonneveld, 1989.
145: Notes:
146: Does not require a symmetric matrix. Does not apply transpose of the matrix.
147: Supports left and right preconditioning, but not symmetric.
149: Developer Notes:
150: Has this weird support for doing the convergence test with the natural norm, I assume this works only with
151: no preconditioning and symmetric positive definite operator.
153: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBCGS, KSPSetPCSide()
154: M*/
155: PETSC_EXTERN PetscErrorCode KSPCreate_CGS(KSP ksp)
156: {
160: ksp->data = (void*)0;
162: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
163: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
164: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
165: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_RIGHT,2);
166: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
167: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
169: ksp->ops->setup = KSPSetUp_CGS;
170: ksp->ops->solve = KSPSolve_CGS;
171: ksp->ops->destroy = KSPDestroyDefault;
172: ksp->ops->buildsolution = KSPBuildSolutionDefault;
173: ksp->ops->buildresidual = KSPBuildResidualDefault;
174: ksp->ops->setfromoptions = NULL;
175: ksp->ops->view = NULL;
176: return(0);
177: }