1 2 /* 3 This file implements the conjugate gradient method in PETSc as part of 4 KSP. You can use this as a starting point for implementing your own 5 Krylov method that is not provided with PETSc. 6 7 The following basic routines are required for each Krylov method. 8 KSPCreate_XXX() - Creates the Krylov context 9 KSPSetFromOptions_XXX() - Sets runtime options 10 KSPSolve_XXX() - Runs the Krylov method 11 KSPDestroy_XXX() - Destroys the Krylov context, freeing all 12 memory it needed 13 Here the "_XXX" denotes a particular implementation, in this case 14 we use _CG (e.g. KSPCreate_CG, KSPDestroy_CG). These routines are 15 are actually called via the common user interface routines 16 KSPSetType(), KSPSetFromOptions(), KSPSolve(), and KSPDestroy() so the 17 application code interface remains identical for all preconditioners. 18 19 Other basic routines for the KSP objects include 20 KSPSetUp_XXX() 21 KSPView_XXX() - Prints details of solver being used. 22 23 Detailed Notes: 24 By default, this code implements the CG (Conjugate Gradient) method, 25 which is valid for real symmetric (and complex Hermitian) positive 26 definite matrices. Note that for the complex Hermitian case, the 27 VecDot() arguments within the code MUST remain in the order given 28 for correct computation of inner products. 29 30 Reference: Hestenes and Steifel, 1952. 31 32 By switching to the indefinite vector inner product, VecTDot(), the 33 same code is used for the complex symmetric case as well. The user 34 must call KSPCGSetType(ksp,KSP_CG_SYMMETRIC) or use the option 35 -ksp_cg_type symmetric to invoke this variant for the complex case. 36 Note, however, that the complex symmetric code is NOT valid for 37 all such matrices ... and thus we don't recommend using this method. 38 */ 39 /* 40 cgimpl.h defines the simple data structured used to store information 41 related to the type of matrix (e.g. complex symmetric) being solved and 42 data used during the optional Lanczo process used to compute eigenvalues 43 */ 44 #include <../src/ksp/ksp/impls/cg/cgimpl.h> /*I "petscksp.h" I*/ 45 extern PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP,PetscReal*,PetscReal*); 46 extern PetscErrorCode KSPComputeEigenvalues_CG(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt*); 47 48 /* 49 KSPSetUp_CG - Sets up the workspace needed by the CG method. 50 51 This is called once, usually automatically by KSPSolve() or KSPSetUp() 52 but can be called directly by KSPSetUp() 53 */ 54 static PetscErrorCode KSPSetUp_CG(KSP ksp) 55 { 56 KSP_CG *cgP = (KSP_CG*)ksp->data; 57 PetscErrorCode ierr; 58 PetscInt maxit = ksp->max_it,nwork = 3; 59 60 PetscFunctionBegin; 61 /* get work vectors needed by CG */ 62 if (cgP->singlereduction) nwork += 2; 63 ierr = KSPSetWorkVecs(ksp,nwork);CHKERRQ(ierr); 64 65 /* 66 If user requested computations of eigenvalues then allocate work 67 work space needed 68 */ 69 if (ksp->calc_sings) { 70 /* get space to store tridiagonal matrix for Lanczos */ 71 ierr = PetscMalloc4(maxit+1,&cgP->e,maxit+1,&cgP->d,maxit+1,&cgP->ee,maxit+1,&cgP->dd);CHKERRQ(ierr); 72 ierr = PetscLogObjectMemory((PetscObject)ksp,2*(maxit+1)*(sizeof(PetscScalar)+sizeof(PetscReal)));CHKERRQ(ierr); 73 74 ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG; 75 ksp->ops->computeeigenvalues = KSPComputeEigenvalues_CG; 76 } 77 PetscFunctionReturn(0); 78 } 79 80 /* 81 A macro used in the following KSPSolve_CG and KSPSolve_CG_SingleReduction routines 82 */ 83 #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a)) 84 85 /* 86 KSPSolve_CG - This routine actually applies the conjugate gradient method 87 88 Note : this routine can be replaced with another one (see below) which implements 89 another variant of CG. 90 91 Input Parameter: 92 . ksp - the Krylov space object that was set to use conjugate gradient, by, for 93 example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG); 94 */ 95 static PetscErrorCode KSPSolve_CG(KSP ksp) 96 { 97 PetscErrorCode ierr; 98 PetscInt i,stored_max_it,eigs; 99 PetscScalar dpi = 0.0,a = 1.0,beta,betaold = 1.0,b = 0,*e = NULL,*d = NULL,dpiold; 100 PetscReal dp = 0.0; 101 Vec X,B,Z,R,P,W; 102 KSP_CG *cg; 103 Mat Amat,Pmat; 104 PetscBool diagonalscale; 105 106 PetscFunctionBegin; 107 ierr = PCGetDiagonalScale(ksp->pc,&diagonalscale);CHKERRQ(ierr); 108 if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name); 109 110 cg = (KSP_CG*)ksp->data; 111 eigs = ksp->calc_sings; 112 stored_max_it = ksp->max_it; 113 X = ksp->vec_sol; 114 B = ksp->vec_rhs; 115 R = ksp->work[0]; 116 Z = ksp->work[1]; 117 P = ksp->work[2]; 118 W = Z; 119 120 if (eigs) {e = cg->e; d = cg->d; e[0] = 0.0; } 121 ierr = PCGetOperators(ksp->pc,&Amat,&Pmat);CHKERRQ(ierr); 122 123 ksp->its = 0; 124 if (!ksp->guess_zero) { 125 ierr = KSP_MatMult(ksp,Amat,X,R);CHKERRQ(ierr); /* r <- b - Ax */ 126 ierr = VecAYPX(R,-1.0,B);CHKERRQ(ierr); 127 } else { 128 ierr = VecCopy(B,R);CHKERRQ(ierr); /* r <- b (x is 0) */ 129 } 130 131 switch (ksp->normtype) { 132 case KSP_NORM_PRECONDITIONED: 133 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 134 ierr = VecNorm(Z,NORM_2,&dp);CHKERRQ(ierr); /* dp <- z'*z = e'*A'*B'*B*A*e */ 135 KSPCheckNorm(ksp,dp); 136 break; 137 case KSP_NORM_UNPRECONDITIONED: 138 ierr = VecNorm(R,NORM_2,&dp);CHKERRQ(ierr); /* dp <- r'*r = e'*A'*A*e */ 139 KSPCheckNorm(ksp,dp); 140 break; 141 case KSP_NORM_NATURAL: 142 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 143 ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr); /* beta <- z'*r */ 144 KSPCheckDot(ksp,beta); 145 dp = PetscSqrtReal(PetscAbsScalar(beta)); /* dp <- r'*z = r'*B*r = e'*A'*B*A*e */ 146 break; 147 case KSP_NORM_NONE: 148 dp = 0.0; 149 break; 150 default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]); 151 } 152 ierr = KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr); 153 ierr = KSPMonitor(ksp,0,dp);CHKERRQ(ierr); 154 ksp->rnorm = dp; 155 156 ierr = (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); /* test for convergence */ 157 if (ksp->reason) PetscFunctionReturn(0); 158 159 if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) { 160 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 161 } 162 if (ksp->normtype != KSP_NORM_NATURAL) { 163 ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr); /* beta <- z'*r */ 164 KSPCheckDot(ksp,beta); 165 } 166 167 i = 0; 168 do { 169 ksp->its = i+1; 170 if (beta == 0.0) { 171 ksp->reason = KSP_CONVERGED_ATOL; 172 ierr = PetscInfo(ksp,"converged due to beta = 0\n");CHKERRQ(ierr); 173 break; 174 #if !defined(PETSC_USE_COMPLEX) 175 } else if ((i > 0) && (beta*betaold < 0.0)) { 176 if (ksp->errorifnotconverged) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Diverged due to indefinite preconditioner, beta %g, betaold %g",(double)beta,(double)betaold); 177 ksp->reason = KSP_DIVERGED_INDEFINITE_PC; 178 ierr = PetscInfo(ksp,"diverging due to indefinite preconditioner\n");CHKERRQ(ierr); 179 break; 180 #endif 181 } 182 if (!i) { 183 ierr = VecCopy(Z,P);CHKERRQ(ierr); /* p <- z */ 184 b = 0.0; 185 } else { 186 b = beta/betaold; 187 if (eigs) { 188 if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues"); 189 e[i] = PetscSqrtReal(PetscAbsScalar(b))/a; 190 } 191 ierr = VecAYPX(P,b,Z);CHKERRQ(ierr); /* p <- z + b* p */ 192 } 193 dpiold = dpi; 194 ierr = KSP_MatMult(ksp,Amat,P,W);CHKERRQ(ierr); /* w <- Ap */ 195 ierr = VecXDot(P,W,&dpi);CHKERRQ(ierr); /* dpi <- p'w */ 196 KSPCheckDot(ksp,dpi); 197 betaold = beta; 198 199 if ((dpi == 0.0) || ((i > 0) && ((PetscSign(PetscRealPart(dpi))*PetscSign(PetscRealPart(dpiold))) < 0.0))) { 200 if (ksp->errorifnotconverged) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Diverged due to indefinite matrix, dpi %g, dpiold %g",(double)PetscRealPart(dpi),(double)PetscRealPart(dpiold)); 201 ksp->reason = KSP_DIVERGED_INDEFINITE_MAT; 202 ierr = PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");CHKERRQ(ierr); 203 break; 204 } 205 a = beta/dpi; /* a = beta/p'w */ 206 if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b))*e[i] + 1.0/a; 207 ierr = VecAXPY(X,a,P);CHKERRQ(ierr); /* x <- x + ap */ 208 ierr = VecAXPY(R,-a,W);CHKERRQ(ierr); /* r <- r - aw */ 209 if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i+2) { 210 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 211 ierr = VecNorm(Z,NORM_2,&dp);CHKERRQ(ierr); /* dp <- z'*z */ 212 KSPCheckNorm(ksp,dp); 213 } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i+2) { 214 ierr = VecNorm(R,NORM_2,&dp);CHKERRQ(ierr); /* dp <- r'*r */ 215 KSPCheckNorm(ksp,dp); 216 } else if (ksp->normtype == KSP_NORM_NATURAL) { 217 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 218 ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr); /* beta <- r'*z */ 219 KSPCheckDot(ksp,beta); 220 dp = PetscSqrtReal(PetscAbsScalar(beta)); 221 } else { 222 dp = 0.0; 223 } 224 ksp->rnorm = dp; 225 ierr = KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr); 226 if (eigs) cg->ned = ksp->its; 227 ierr = KSPMonitor(ksp,i+1,dp);CHKERRQ(ierr); 228 ierr = (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); 229 if (ksp->reason) break; 230 231 if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i+2)) { 232 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 233 } 234 if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i+2)) { 235 ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr); /* beta <- z'*r */ 236 KSPCheckDot(ksp,beta); 237 } 238 239 i++; 240 } while (i<ksp->max_it); 241 if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS; 242 PetscFunctionReturn(0); 243 } 244 245 /* 246 KSPSolve_CG_SingleReduction 247 248 This variant of CG is identical in exact arithmetic to the standard algorithm, 249 but is rearranged to use only a single reduction stage per iteration, using additional 250 intermediate vectors. 251 252 See KSPCGUseSingleReduction_CG() 253 254 */ 255 static PetscErrorCode KSPSolve_CG_SingleReduction(KSP ksp) 256 { 257 PetscErrorCode ierr; 258 PetscInt i,stored_max_it,eigs; 259 PetscScalar dpi = 0.0,a = 1.0,beta,betaold = 1.0,b = 0,*e = NULL,*d = NULL,delta,dpiold,tmp[2]; 260 PetscReal dp = 0.0; 261 Vec X,B,Z,R,P,S,W,tmpvecs[2]; 262 KSP_CG *cg; 263 Mat Amat,Pmat; 264 PetscBool diagonalscale; 265 266 PetscFunctionBegin; 267 ierr = PCGetDiagonalScale(ksp->pc,&diagonalscale);CHKERRQ(ierr); 268 if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name); 269 270 cg = (KSP_CG*)ksp->data; 271 eigs = ksp->calc_sings; 272 stored_max_it = ksp->max_it; 273 X = ksp->vec_sol; 274 B = ksp->vec_rhs; 275 R = ksp->work[0]; 276 Z = ksp->work[1]; 277 P = ksp->work[2]; 278 S = ksp->work[3]; 279 W = ksp->work[4]; 280 281 if (eigs) {e = cg->e; d = cg->d; e[0] = 0.0; } 282 ierr = PCGetOperators(ksp->pc,&Amat,&Pmat);CHKERRQ(ierr); 283 284 ksp->its = 0; 285 if (!ksp->guess_zero) { 286 ierr = KSP_MatMult(ksp,Amat,X,R);CHKERRQ(ierr); /* r <- b - Ax */ 287 ierr = VecAYPX(R,-1.0,B);CHKERRQ(ierr); 288 } else { 289 ierr = VecCopy(B,R);CHKERRQ(ierr); /* r <- b (x is 0) */ 290 } 291 292 switch (ksp->normtype) { 293 case KSP_NORM_PRECONDITIONED: 294 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 295 ierr = VecNorm(Z,NORM_2,&dp);CHKERRQ(ierr); /* dp <- z'*z = e'*A'*B'*B*A'*e' */ 296 KSPCheckNorm(ksp,dp); 297 break; 298 case KSP_NORM_UNPRECONDITIONED: 299 ierr = VecNorm(R,NORM_2,&dp);CHKERRQ(ierr); /* dp <- r'*r = e'*A'*A*e */ 300 KSPCheckNorm(ksp,dp); 301 break; 302 case KSP_NORM_NATURAL: 303 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 304 ierr = KSP_MatMult(ksp,Amat,Z,S);CHKERRQ(ierr); 305 ierr = VecXDot(Z,S,&delta);CHKERRQ(ierr); /* delta <- z'*A*z = r'*B*A*B*r */ 306 ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr); /* beta <- z'*r */ 307 KSPCheckDot(ksp,beta); 308 dp = PetscSqrtReal(PetscAbsScalar(beta)); /* dp <- r'*z = r'*B*r = e'*A'*B*A*e */ 309 break; 310 case KSP_NORM_NONE: 311 dp = 0.0; 312 break; 313 default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]); 314 } 315 ierr = KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr); 316 ierr = KSPMonitor(ksp,0,dp);CHKERRQ(ierr); 317 ksp->rnorm = dp; 318 319 ierr = (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); /* test for convergence */ 320 if (ksp->reason) PetscFunctionReturn(0); 321 322 if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) { 323 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 324 } 325 if (ksp->normtype != KSP_NORM_NATURAL) { 326 ierr = KSP_MatMult(ksp,Amat,Z,S);CHKERRQ(ierr); 327 ierr = VecXDot(Z,S,&delta);CHKERRQ(ierr); /* delta <- z'*A*z = r'*B*A*B*r */ 328 ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr); /* beta <- z'*r */ 329 KSPCheckDot(ksp,beta); 330 } 331 332 i = 0; 333 do { 334 ksp->its = i+1; 335 if (beta == 0.0) { 336 ksp->reason = KSP_CONVERGED_ATOL; 337 ierr = PetscInfo(ksp,"converged due to beta = 0\n");CHKERRQ(ierr); 338 break; 339 #if !defined(PETSC_USE_COMPLEX) 340 } else if ((i > 0) && (beta*betaold < 0.0)) { 341 if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Diverged due to indefinite preconditioner"); 342 ksp->reason = KSP_DIVERGED_INDEFINITE_PC; 343 ierr = PetscInfo(ksp,"diverging due to indefinite preconditioner\n");CHKERRQ(ierr); 344 break; 345 #endif 346 } 347 if (!i) { 348 ierr = VecCopy(Z,P);CHKERRQ(ierr); /* p <- z */ 349 b = 0.0; 350 } else { 351 b = beta/betaold; 352 if (eigs) { 353 if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues"); 354 e[i] = PetscSqrtReal(PetscAbsScalar(b))/a; 355 } 356 ierr = VecAYPX(P,b,Z);CHKERRQ(ierr); /* p <- z + b* p */ 357 } 358 dpiold = dpi; 359 if (!i) { 360 ierr = KSP_MatMult(ksp,Amat,P,W);CHKERRQ(ierr); /* w <- Ap */ 361 ierr = VecXDot(P,W,&dpi);CHKERRQ(ierr); /* dpi <- p'w */ 362 } else { 363 ierr = VecAYPX(W,beta/betaold,S);CHKERRQ(ierr); /* w <- Ap */ 364 dpi = delta - beta*beta*dpiold/(betaold*betaold); /* dpi <- p'w */ 365 } 366 betaold = beta; 367 KSPCheckDot(ksp,beta); 368 369 if ((dpi == 0.0) || ((i > 0) && (PetscRealPart(dpi*dpiold) <= 0.0))) { 370 if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Diverged due to indefinite matrix"); 371 ksp->reason = KSP_DIVERGED_INDEFINITE_MAT; 372 ierr = PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");CHKERRQ(ierr); 373 break; 374 } 375 a = beta/dpi; /* a = beta/p'w */ 376 if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b))*e[i] + 1.0/a; 377 ierr = VecAXPY(X,a,P);CHKERRQ(ierr); /* x <- x + ap */ 378 ierr = VecAXPY(R,-a,W);CHKERRQ(ierr); /* r <- r - aw */ 379 if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i+2) { 380 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 381 ierr = KSP_MatMult(ksp,Amat,Z,S);CHKERRQ(ierr); 382 ierr = VecNorm(Z,NORM_2,&dp);CHKERRQ(ierr); /* dp <- z'*z */ 383 KSPCheckNorm(ksp,dp); 384 } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i+2) { 385 ierr = VecNorm(R,NORM_2,&dp);CHKERRQ(ierr); /* dp <- r'*r */ 386 KSPCheckNorm(ksp,dp); 387 } else if (ksp->normtype == KSP_NORM_NATURAL) { 388 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 389 tmpvecs[0] = S; tmpvecs[1] = R; 390 ierr = KSP_MatMult(ksp,Amat,Z,S);CHKERRQ(ierr); 391 ierr = VecMDot(Z,2,tmpvecs,tmp);CHKERRQ(ierr); /* delta <- z'*A*z = r'*B*A*B*r */ 392 delta = tmp[0]; beta = tmp[1]; /* beta <- z'*r */ 393 KSPCheckDot(ksp,beta); 394 dp = PetscSqrtReal(PetscAbsScalar(beta)); /* dp <- r'*z = r'*B*r = e'*A'*B*A*e */ 395 } else { 396 dp = 0.0; 397 } 398 ksp->rnorm = dp; 399 ierr = KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr); 400 if (eigs) cg->ned = ksp->its; 401 ierr = KSPMonitor(ksp,i+1,dp);CHKERRQ(ierr); 402 ierr = (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); 403 if (ksp->reason) break; 404 405 if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i+2)) { 406 ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ 407 ierr = KSP_MatMult(ksp,Amat,Z,S);CHKERRQ(ierr); 408 } 409 if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i+2)) { 410 tmpvecs[0] = S; tmpvecs[1] = R; 411 ierr = VecMDot(Z,2,tmpvecs,tmp);CHKERRQ(ierr); 412 delta = tmp[0]; beta = tmp[1]; /* delta <- z'*A*z = r'*B'*A*B*r */ 413 KSPCheckDot(ksp,beta); /* beta <- z'*r */ 414 } 415 416 i++; 417 } while (i<ksp->max_it); 418 if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS; 419 PetscFunctionReturn(0); 420 } 421 422 /* 423 KSPDestroy_CG - Frees resources allocated in KSPSetup_CG and clears function 424 compositions from KSPCreate_CG. If adding your own KSP implementation, 425 you must be sure to free all allocated resources here to prevent 426 leaks. 427 */ 428 PetscErrorCode KSPDestroy_CG(KSP ksp) 429 { 430 KSP_CG *cg = (KSP_CG*)ksp->data; 431 PetscErrorCode ierr; 432 433 PetscFunctionBegin; 434 /* free space used for singular value calculations */ 435 if (ksp->calc_sings) { 436 ierr = PetscFree4(cg->e,cg->d,cg->ee,cg->dd);CHKERRQ(ierr); 437 } 438 ierr = KSPDestroyDefault(ksp);CHKERRQ(ierr); 439 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",NULL);CHKERRQ(ierr); 440 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPCGUseSingleReduction_C",NULL);CHKERRQ(ierr); 441 PetscFunctionReturn(0); 442 } 443 444 /* 445 KSPView_CG - Prints information about the current Krylov method being used. 446 If your Krylov method has special options or flags that information 447 should be printed here. 448 */ 449 PetscErrorCode KSPView_CG(KSP ksp,PetscViewer viewer) 450 { 451 KSP_CG *cg = (KSP_CG*)ksp->data; 452 PetscErrorCode ierr; 453 PetscBool iascii; 454 455 PetscFunctionBegin; 456 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 457 if (iascii) { 458 #if defined(PETSC_USE_COMPLEX) 459 ierr = PetscViewerASCIIPrintf(viewer," variant %s\n",KSPCGTypes[cg->type]);CHKERRQ(ierr); 460 #endif 461 if (cg->singlereduction) { 462 ierr = PetscViewerASCIIPrintf(viewer," using single-reduction variant\n");CHKERRQ(ierr); 463 } 464 } 465 PetscFunctionReturn(0); 466 } 467 468 /* 469 KSPSetFromOptions_CG - Checks the options database for options related to the 470 conjugate gradient method. 471 */ 472 PetscErrorCode KSPSetFromOptions_CG(PetscOptionItems *PetscOptionsObject,KSP ksp) 473 { 474 PetscErrorCode ierr; 475 KSP_CG *cg = (KSP_CG*)ksp->data; 476 PetscBool flg; 477 478 PetscFunctionBegin; 479 ierr = PetscOptionsHead(PetscOptionsObject,"KSP CG and CGNE options");CHKERRQ(ierr); 480 #if defined(PETSC_USE_COMPLEX) 481 ierr = PetscOptionsEnum("-ksp_cg_type","Matrix is Hermitian or complex symmetric","KSPCGSetType",KSPCGTypes,(PetscEnum)cg->type, 482 (PetscEnum*)&cg->type,NULL);CHKERRQ(ierr); 483 #endif 484 ierr = PetscOptionsBool("-ksp_cg_single_reduction","Merge inner products into single MPIU_Allreduce()","KSPCGUseSingleReduction",cg->singlereduction,&cg->singlereduction,&flg);CHKERRQ(ierr); 485 if (flg) { 486 ierr = KSPCGUseSingleReduction(ksp,cg->singlereduction);CHKERRQ(ierr); 487 } 488 ierr = PetscOptionsTail();CHKERRQ(ierr); 489 PetscFunctionReturn(0); 490 } 491 492 /* 493 KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method. 494 This routine is registered below in KSPCreate_CG() and called from the 495 routine KSPCGSetType() (see the file cgtype.c). 496 */ 497 PetscErrorCode KSPCGSetType_CG(KSP ksp,KSPCGType type) 498 { 499 KSP_CG *cg = (KSP_CG*)ksp->data; 500 501 PetscFunctionBegin; 502 cg->type = type; 503 PetscFunctionReturn(0); 504 } 505 506 /* 507 KSPCGUseSingleReduction_CG 508 509 This routine sets a flag to use a variant of CG. Note that (in somewhat 510 atypical fashion) it also swaps out the routine called when KSPSolve() 511 is invoked. 512 */ 513 static PetscErrorCode KSPCGUseSingleReduction_CG(KSP ksp,PetscBool flg) 514 { 515 KSP_CG *cg = (KSP_CG*)ksp->data; 516 517 PetscFunctionBegin; 518 cg->singlereduction = flg; 519 if (cg->singlereduction) { 520 ksp->ops->solve = KSPSolve_CG_SingleReduction; 521 } else { 522 ksp->ops->solve = KSPSolve_CG; 523 } 524 PetscFunctionReturn(0); 525 } 526 527 /* 528 KSPCreate_CG - Creates the data structure for the Krylov method CG and sets the 529 function pointers for all the routines it needs to call (KSPSolve_CG() etc) 530 531 It must be labeled as PETSC_EXTERN to be dynamically linkable in C++ 532 */ 533 /*MC 534 KSPCG - The Preconditioned Conjugate Gradient (PCG) iterative method 535 536 Options Database Keys: 537 + -ksp_cg_type Hermitian - (for complex matrices only) indicates the matrix is Hermitian, see KSPCGSetType() 538 . -ksp_cg_type symmetric - (for complex matrices only) indicates the matrix is symmetric 539 - -ksp_cg_single_reduction - performs both inner products needed in the algorithm with a single MPIU_Allreduce() call, see KSPCGUseSingleReduction() 540 541 Level: beginner 542 543 Notes: 544 The PCG method requires both the matrix and preconditioner to be symmetric positive (or negative) (semi) definite. 545 546 Only left preconditioning is supported; there are several ways to motivate preconditioned CG, but they all produce the same algorithm. 547 One can interpret preconditioning A with B to mean any of the following\: 548 .n (1) Solve a left-preconditioned system BAx = Bb, using inv(B) to define an inner product in the algorithm. 549 .n (2) Solve a right-preconditioned system ABy = b, x = By, using B to define an inner product in the algorithm. 550 .n (3) Solve a symmetrically-preconditioned system, E^TAEy = E^Tb, x = Ey, where B = EE^T. 551 .n (4) Solve Ax=b with CG, but use the inner product defined by B to define the method [2]. 552 .n In all cases, the resulting algorithm only requires application of B to vectors. 553 554 For complex numbers there are two different CG methods, one for Hermitian symmetric matrices and one for non-Hermitian symmetric matrices. Use 555 KSPCGSetType() to indicate which type you are using. 556 557 Developer Notes: 558 KSPSolve_CG() should actually query the matrix to determine if it is Hermitian symmetric or not and NOT require the user to 559 indicate it to the KSP object. 560 561 References: 562 + 1. - Magnus R. Hestenes and Eduard Stiefel, Methods of Conjugate Gradients for Solving Linear Systems, 563 Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379 564 - 2. - Josef Malek and Zdenek Strakos, Preconditioning and the Conjugate Gradient Method in the Context of Solving PDEs, 565 SIAM, 2014. 566 567 .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, 568 KSPCGSetType(), KSPCGUseSingleReduction(), KSPPIPECG, KSPGROPPCG 569 570 M*/ 571 PETSC_EXTERN PetscErrorCode KSPCreate_CG(KSP ksp) 572 { 573 PetscErrorCode ierr; 574 KSP_CG *cg; 575 576 PetscFunctionBegin; 577 ierr = PetscNewLog(ksp,&cg);CHKERRQ(ierr); 578 #if !defined(PETSC_USE_COMPLEX) 579 cg->type = KSP_CG_SYMMETRIC; 580 #else 581 cg->type = KSP_CG_HERMITIAN; 582 #endif 583 ksp->data = (void*)cg; 584 585 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);CHKERRQ(ierr); 586 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);CHKERRQ(ierr); 587 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);CHKERRQ(ierr); 588 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);CHKERRQ(ierr); 589 590 /* 591 Sets the functions that are associated with this data structure 592 (in C++ this is the same as defining virtual functions) 593 */ 594 ksp->ops->setup = KSPSetUp_CG; 595 ksp->ops->solve = KSPSolve_CG; 596 ksp->ops->destroy = KSPDestroy_CG; 597 ksp->ops->view = KSPView_CG; 598 ksp->ops->setfromoptions = KSPSetFromOptions_CG; 599 ksp->ops->buildsolution = KSPBuildSolutionDefault; 600 ksp->ops->buildresidual = KSPBuildResidualDefault; 601 602 /* 603 Attach the function KSPCGSetType_CG() to this object. The routine 604 KSPCGSetType() checks for this attached function and calls it if it finds 605 it. (Sort of like a dynamic member function that can be added at run time 606 */ 607 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",KSPCGSetType_CG);CHKERRQ(ierr); 608 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPCGUseSingleReduction_C",KSPCGUseSingleReduction_CG);CHKERRQ(ierr); 609 PetscFunctionReturn(0); 610 } 611