1 /* 2 This file implements the FCG (Flexible Conjugate Gradient) method 3 4 */ 5 6 #include <../src/ksp/ksp/impls/fcg/fcgimpl.h> /*I "petscksp.h" I*/ 7 extern PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP,PetscReal*,PetscReal*); 8 extern PetscErrorCode KSPComputeEigenvalues_CG(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt*); 9 10 #define KSPFCG_DEFAULT_MMAX 30 /* maximum number of search directions to keep */ 11 #define KSPFCG_DEFAULT_NPREALLOC 10 /* number of search directions to preallocate */ 12 #define KSPFCG_DEFAULT_VECB 5 /* number of search directions to allocate each time new direction vectors are needed */ 13 #define KSPFCG_DEFAULT_TRUNCSTRAT KSP_FCD_TRUNC_TYPE_NOTAY 14 15 static PetscErrorCode KSPAllocateVectors_FCG(KSP ksp, PetscInt nvecsneeded, PetscInt chunksize) 16 { 17 PetscInt i; 18 KSP_FCG *fcg = (KSP_FCG*)ksp->data; 19 PetscInt nnewvecs, nvecsprev; 20 21 PetscFunctionBegin; 22 /* Allocate enough new vectors to add chunksize new vectors, reach nvecsneedtotal, or to reach mmax+1, whichever is smallest */ 23 if (fcg->nvecs < PetscMin(fcg->mmax+1,nvecsneeded)) { 24 nvecsprev = fcg->nvecs; 25 nnewvecs = PetscMin(PetscMax(nvecsneeded-fcg->nvecs,chunksize),fcg->mmax+1-fcg->nvecs); 26 PetscCall(KSPCreateVecs(ksp,nnewvecs,&fcg->pCvecs[fcg->nchunks],0,NULL)); 27 PetscCall(PetscLogObjectParents((PetscObject)ksp,nnewvecs,fcg->pCvecs[fcg->nchunks])); 28 PetscCall(KSPCreateVecs(ksp,nnewvecs,&fcg->pPvecs[fcg->nchunks],0,NULL)); 29 PetscCall(PetscLogObjectParents((PetscObject)ksp,nnewvecs,fcg->pPvecs[fcg->nchunks])); 30 fcg->nvecs += nnewvecs; 31 for (i=0;i<nnewvecs;++i) { 32 fcg->Cvecs[nvecsprev + i] = fcg->pCvecs[fcg->nchunks][i]; 33 fcg->Pvecs[nvecsprev + i] = fcg->pPvecs[fcg->nchunks][i]; 34 } 35 fcg->chunksizes[fcg->nchunks] = nnewvecs; 36 ++fcg->nchunks; 37 } 38 PetscFunctionReturn(0); 39 } 40 41 static PetscErrorCode KSPSetUp_FCG(KSP ksp) 42 { 43 KSP_FCG *fcg = (KSP_FCG*)ksp->data; 44 PetscInt maxit = ksp->max_it; 45 const PetscInt nworkstd = 2; 46 47 PetscFunctionBegin; 48 49 /* Allocate "standard" work vectors (not including the basis and transformed basis vectors) */ 50 PetscCall(KSPSetWorkVecs(ksp,nworkstd)); 51 52 /* Allocated space for pointers to additional work vectors 53 note that mmax is the number of previous directions, so we add 1 for the current direction, 54 and an extra 1 for the prealloc (which might be empty) */ 55 PetscCall(PetscMalloc5(fcg->mmax+1,&fcg->Pvecs,fcg->mmax+1,&fcg->Cvecs,fcg->mmax+1,&fcg->pPvecs,fcg->mmax+1,&fcg->pCvecs,fcg->mmax+2,&fcg->chunksizes)); 56 PetscCall(PetscLogObjectMemory((PetscObject)ksp,2*(fcg->mmax+1)*sizeof(Vec*) + 2*(fcg->mmax + 1)*sizeof(Vec**) + (fcg->mmax + 2)*sizeof(PetscInt))); 57 58 /* If the requested number of preallocated vectors is greater than mmax reduce nprealloc */ 59 if (fcg->nprealloc > fcg->mmax+1) { 60 PetscCall(PetscInfo(NULL,"Requested nprealloc=%" PetscInt_FMT " is greater than m_max+1=%" PetscInt_FMT ". Resetting nprealloc = m_max+1.\n",fcg->nprealloc, fcg->mmax+1)); 61 } 62 63 /* Preallocate additional work vectors */ 64 PetscCall(KSPAllocateVectors_FCG(ksp,fcg->nprealloc,fcg->nprealloc)); 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 PetscCall(PetscMalloc4(maxit,&fcg->e,maxit,&fcg->d,maxit,&fcg->ee,maxit,&fcg->dd)); 72 PetscCall(PetscLogObjectMemory((PetscObject)ksp,2*(maxit+1)*(sizeof(PetscScalar)+sizeof(PetscReal)))); 73 74 ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG; 75 ksp->ops->computeeigenvalues = KSPComputeEigenvalues_CG; 76 } 77 PetscFunctionReturn(0); 78 } 79 80 static PetscErrorCode KSPSolve_FCG(KSP ksp) 81 { 82 PetscInt i,k,idx,mi; 83 KSP_FCG *fcg = (KSP_FCG*)ksp->data; 84 PetscScalar alpha=0.0,beta = 0.0,dpi,s; 85 PetscReal dp=0.0; 86 Vec B,R,Z,X,Pcurr,Ccurr; 87 Mat Amat,Pmat; 88 PetscInt eigs = ksp->calc_sings; /* Variables for eigen estimation - START*/ 89 PetscInt stored_max_it = ksp->max_it; 90 PetscScalar alphaold = 0,betaold = 1.0,*e = NULL,*d = NULL;/* Variables for eigen estimation - FINISH */ 91 92 PetscFunctionBegin; 93 94 #define VecXDot(x,y,a) (((fcg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a)) 95 #define VecXMDot(a,b,c,d) (((fcg->type) == (KSP_CG_HERMITIAN)) ? VecMDot(a,b,c,d) : VecMTDot(a,b,c,d)) 96 97 X = ksp->vec_sol; 98 B = ksp->vec_rhs; 99 R = ksp->work[0]; 100 Z = ksp->work[1]; 101 102 PetscCall(PCGetOperators(ksp->pc,&Amat,&Pmat)); 103 if (eigs) {e = fcg->e; d = fcg->d; e[0] = 0.0; } 104 /* Compute initial residual needed for convergence check*/ 105 ksp->its = 0; 106 if (!ksp->guess_zero) { 107 PetscCall(KSP_MatMult(ksp,Amat,X,R)); 108 PetscCall(VecAYPX(R,-1.0,B)); /* r <- b - Ax */ 109 } else { 110 PetscCall(VecCopy(B,R)); /* r <- b (x is 0) */ 111 } 112 switch (ksp->normtype) { 113 case KSP_NORM_PRECONDITIONED: 114 PetscCall(KSP_PCApply(ksp,R,Z)); /* z <- Br */ 115 PetscCall(VecNorm(Z,NORM_2,&dp)); /* dp <- dqrt(z'*z) = sqrt(e'*A'*B'*B*A*e) */ 116 KSPCheckNorm(ksp,dp); 117 break; 118 case KSP_NORM_UNPRECONDITIONED: 119 PetscCall(VecNorm(R,NORM_2,&dp)); /* dp <- sqrt(r'*r) = sqrt(e'*A'*A*e) */ 120 KSPCheckNorm(ksp,dp); 121 break; 122 case KSP_NORM_NATURAL: 123 PetscCall(KSP_PCApply(ksp,R,Z)); /* z <- Br */ 124 PetscCall(VecXDot(R,Z,&s)); 125 KSPCheckDot(ksp,s); 126 dp = PetscSqrtReal(PetscAbsScalar(s)); /* dp <- sqrt(r'*z) = sqrt(e'*A'*B*A*e) */ 127 break; 128 case KSP_NORM_NONE: 129 dp = 0.0; 130 break; 131 default: SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]); 132 } 133 134 /* Initial Convergence Check */ 135 PetscCall(KSPLogResidualHistory(ksp,dp)); 136 PetscCall(KSPMonitor(ksp,0,dp)); 137 ksp->rnorm = dp; 138 if (ksp->normtype == KSP_NORM_NONE) { 139 PetscCall(KSPConvergedSkip(ksp,0,dp,&ksp->reason,ksp->cnvP)); 140 } else { 141 PetscCall((*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP)); 142 } 143 if (ksp->reason) PetscFunctionReturn(0); 144 145 /* Apply PC if not already done for convergence check */ 146 if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->normtype == KSP_NORM_NONE) { 147 PetscCall(KSP_PCApply(ksp,R,Z)); /* z <- Br */ 148 } 149 150 i = 0; 151 do { 152 ksp->its = i+1; 153 154 /* If needbe, allocate a new chunk of vectors in P and C */ 155 PetscCall(KSPAllocateVectors_FCG(ksp,i+1,fcg->vecb)); 156 157 /* Note that we wrap around and start clobbering old vectors */ 158 idx = i % (fcg->mmax+1); 159 Pcurr = fcg->Pvecs[idx]; 160 Ccurr = fcg->Cvecs[idx]; 161 162 /* number of old directions to orthogonalize against */ 163 switch(fcg->truncstrat) { 164 case KSP_FCD_TRUNC_TYPE_STANDARD: 165 mi = fcg->mmax; 166 break; 167 case KSP_FCD_TRUNC_TYPE_NOTAY: 168 mi = ((i-1) % fcg->mmax)+1; 169 break; 170 default: 171 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unrecognized Truncation Strategy"); 172 } 173 174 /* Compute a new column of P (Currently does not support modified G-S or iterative refinement)*/ 175 PetscCall(VecCopy(Z,Pcurr)); 176 177 { 178 PetscInt l,ndots; 179 180 l = PetscMax(0,i-mi); 181 ndots = i-l; 182 if (ndots) { 183 PetscInt j; 184 Vec *Pold, *Cold; 185 PetscScalar *dots; 186 187 PetscCall(PetscMalloc3(ndots,&dots,ndots,&Cold,ndots,&Pold)); 188 for (k=l,j=0;j<ndots;++k,++j) { 189 idx = k % (fcg->mmax+1); 190 Cold[j] = fcg->Cvecs[idx]; 191 Pold[j] = fcg->Pvecs[idx]; 192 } 193 PetscCall(VecXMDot(Z,ndots,Cold,dots)); 194 for (k=0;k<ndots;++k) { 195 dots[k] = -dots[k]; 196 } 197 PetscCall(VecMAXPY(Pcurr,ndots,dots,Pold)); 198 PetscCall(PetscFree3(dots,Cold,Pold)); 199 } 200 } 201 202 /* Update X and R */ 203 betaold = beta; 204 PetscCall(VecXDot(Pcurr,R,&beta)); /* beta <- pi'*r */ 205 KSPCheckDot(ksp,beta); 206 PetscCall(KSP_MatMult(ksp,Amat,Pcurr,Ccurr)); /* w <- A*pi (stored in ci) */ 207 PetscCall(VecXDot(Pcurr,Ccurr,&dpi)); /* dpi <- pi'*w */ 208 alphaold = alpha; 209 alpha = beta / dpi; /* alpha <- beta/dpi */ 210 PetscCall(VecAXPY(X,alpha,Pcurr)); /* x <- x + alpha * pi */ 211 PetscCall(VecAXPY(R,-alpha,Ccurr)); /* r <- r - alpha * wi */ 212 213 /* Compute norm for convergence check */ 214 switch (ksp->normtype) { 215 case KSP_NORM_PRECONDITIONED: 216 PetscCall(KSP_PCApply(ksp,R,Z)); /* z <- Br */ 217 PetscCall(VecNorm(Z,NORM_2,&dp)); /* dp <- sqrt(z'*z) = sqrt(e'*A'*B'*B*A*e) */ 218 KSPCheckNorm(ksp,dp); 219 break; 220 case KSP_NORM_UNPRECONDITIONED: 221 PetscCall(VecNorm(R,NORM_2,&dp)); /* dp <- sqrt(r'*r) = sqrt(e'*A'*A*e) */ 222 KSPCheckNorm(ksp,dp); 223 break; 224 case KSP_NORM_NATURAL: 225 PetscCall(KSP_PCApply(ksp,R,Z)); /* z <- Br */ 226 PetscCall(VecXDot(R,Z,&s)); 227 KSPCheckDot(ksp,s); 228 dp = PetscSqrtReal(PetscAbsScalar(s)); /* dp <- sqrt(r'*z) = sqrt(e'*A'*B*A*e) */ 229 break; 230 case KSP_NORM_NONE: 231 dp = 0.0; 232 break; 233 default: SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]); 234 } 235 236 /* Check for convergence */ 237 ksp->rnorm = dp; 238 PetscCall(KSPLogResidualHistory(ksp,dp)); 239 PetscCall(KSPMonitor(ksp,i+1,dp)); 240 PetscCall((*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP)); 241 if (ksp->reason) break; 242 243 /* Apply PC if not already done for convergence check */ 244 if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->normtype == KSP_NORM_NONE) { 245 PetscCall(KSP_PCApply(ksp,R,Z)); /* z <- Br */ 246 } 247 248 /* Compute current C (which is W/dpi) */ 249 PetscCall(VecScale(Ccurr,1.0/dpi)); /* w <- ci/dpi */ 250 251 if (eigs) { 252 if (i > 0) { 253 PetscCheck(ksp->max_it == stored_max_it,PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues"); 254 e[i] = PetscSqrtReal(PetscAbsScalar(beta/betaold))/alphaold; 255 d[i] = PetscSqrtReal(PetscAbsScalar(beta/betaold))*e[i] + 1.0/alpha; 256 } else { 257 d[i] = PetscSqrtReal(PetscAbsScalar(beta))*e[i] + 1.0/alpha; 258 } 259 fcg->ned = ksp->its-1; 260 } 261 ++i; 262 } while (i<ksp->max_it); 263 if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS; 264 PetscFunctionReturn(0); 265 } 266 267 static PetscErrorCode KSPDestroy_FCG(KSP ksp) 268 { 269 PetscInt i; 270 KSP_FCG *fcg = (KSP_FCG*)ksp->data; 271 272 PetscFunctionBegin; 273 274 /* Destroy "standard" work vecs */ 275 VecDestroyVecs(ksp->nwork,&ksp->work); 276 277 /* Destroy P and C vectors and the arrays that manage pointers to them */ 278 if (fcg->nvecs) { 279 for (i=0;i<fcg->nchunks;++i) { 280 PetscCall(VecDestroyVecs(fcg->chunksizes[i],&fcg->pPvecs[i])); 281 PetscCall(VecDestroyVecs(fcg->chunksizes[i],&fcg->pCvecs[i])); 282 } 283 } 284 PetscCall(PetscFree5(fcg->Pvecs,fcg->Cvecs,fcg->pPvecs,fcg->pCvecs,fcg->chunksizes)); 285 /* free space used for singular value calculations */ 286 if (ksp->calc_sings) PetscCall(PetscFree4(fcg->e,fcg->d,fcg->ee,fcg->dd)); 287 PetscCall(KSPDestroyDefault(ksp)); 288 PetscFunctionReturn(0); 289 } 290 291 static PetscErrorCode KSPView_FCG(KSP ksp,PetscViewer viewer) 292 { 293 KSP_FCG *fcg = (KSP_FCG*)ksp->data; 294 PetscBool iascii,isstring; 295 const char *truncstr; 296 297 PetscFunctionBegin; 298 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 299 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring)); 300 301 if (fcg->truncstrat == KSP_FCD_TRUNC_TYPE_STANDARD) truncstr = "Using standard truncation strategy"; 302 else if (fcg->truncstrat == KSP_FCD_TRUNC_TYPE_NOTAY) truncstr = "Using Notay's truncation strategy"; 303 else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Undefined FCG truncation strategy"); 304 305 if (iascii) { 306 PetscCall(PetscViewerASCIIPrintf(viewer," m_max=%" PetscInt_FMT "\n",fcg->mmax)); 307 PetscCall(PetscViewerASCIIPrintf(viewer," preallocated %" PetscInt_FMT " directions\n",PetscMin(fcg->nprealloc,fcg->mmax+1))); 308 PetscCall(PetscViewerASCIIPrintf(viewer," %s\n",truncstr)); 309 } else if (isstring) { 310 PetscCall(PetscViewerStringSPrintf(viewer,"m_max %" PetscInt_FMT " nprealloc %" PetscInt_FMT " %s",fcg->mmax,fcg->nprealloc,truncstr)); 311 } 312 PetscFunctionReturn(0); 313 } 314 315 /*@ 316 KSPFCGSetMmax - set the maximum number of previous directions FCG will store for orthogonalization 317 318 Note: mmax + 1 directions are stored (mmax previous ones along with a current one) 319 and whether all are used in each iteration also depends on the truncation strategy 320 (see KSPFCGSetTruncationType()) 321 322 Logically Collective on ksp 323 324 Input Parameters: 325 + ksp - the Krylov space context 326 - mmax - the maximum number of previous directions to orthogonalize againt 327 328 Level: intermediate 329 330 .seealso: `KSPFCG`, `KSPFCGGetTruncationType()`, `KSPFCGGetNprealloc()` 331 @*/ 332 PetscErrorCode KSPFCGSetMmax(KSP ksp,PetscInt mmax) 333 { 334 KSP_FCG *fcg = (KSP_FCG*)ksp->data; 335 336 PetscFunctionBegin; 337 PetscValidHeaderSpecific(ksp,KSP_CLASSID,1); 338 PetscValidLogicalCollectiveInt(ksp,mmax,2); 339 fcg->mmax = mmax; 340 PetscFunctionReturn(0); 341 } 342 343 /*@ 344 KSPFCGGetMmax - get the maximum number of previous directions FCG will store 345 346 Note: FCG stores mmax+1 directions at most (mmax previous ones, and one current one) 347 348 Not Collective 349 350 Input Parameter: 351 . ksp - the Krylov space context 352 353 Output Parameter: 354 . mmax - the maximum number of previous directions allowed for orthogonalization 355 356 Level: intermediate 357 358 .seealso: `KSPFCG`, `KSPFCGGetTruncationType()`, `KSPFCGGetNprealloc()`, `KSPFCGSetMmax()` 359 @*/ 360 361 PetscErrorCode KSPFCGGetMmax(KSP ksp,PetscInt *mmax) 362 { 363 KSP_FCG *fcg=(KSP_FCG*)ksp->data; 364 365 PetscFunctionBegin; 366 PetscValidHeaderSpecific(ksp,KSP_CLASSID,1); 367 *mmax = fcg->mmax; 368 PetscFunctionReturn(0); 369 } 370 371 /*@ 372 KSPFCGSetNprealloc - set the number of directions to preallocate with FCG 373 374 Logically Collective on ksp 375 376 Input Parameters: 377 + ksp - the Krylov space context 378 - nprealloc - the number of vectors to preallocate 379 380 Level: advanced 381 382 Options Database: 383 . -ksp_fcg_nprealloc <N> - number of directions to preallocate 384 385 .seealso: `KSPFCG`, `KSPFCGGetTruncationType()`, `KSPFCGGetNprealloc()` 386 @*/ 387 PetscErrorCode KSPFCGSetNprealloc(KSP ksp,PetscInt nprealloc) 388 { 389 KSP_FCG *fcg=(KSP_FCG*)ksp->data; 390 391 PetscFunctionBegin; 392 PetscValidHeaderSpecific(ksp,KSP_CLASSID,1); 393 PetscValidLogicalCollectiveInt(ksp,nprealloc,2); 394 PetscCheck(nprealloc <= fcg->mmax+1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Cannot preallocate more than m_max+1 vectors"); 395 fcg->nprealloc = nprealloc; 396 PetscFunctionReturn(0); 397 } 398 399 /*@ 400 KSPFCGGetNprealloc - get the number of directions preallocate by FCG 401 402 Not Collective 403 404 Input Parameter: 405 . ksp - the Krylov space context 406 407 Output Parameter: 408 . nprealloc - the number of directions preallocated 409 410 Level: advanced 411 412 .seealso: `KSPFCG`, `KSPFCGGetTruncationType()`, `KSPFCGSetNprealloc()` 413 @*/ 414 PetscErrorCode KSPFCGGetNprealloc(KSP ksp,PetscInt *nprealloc) 415 { 416 KSP_FCG *fcg=(KSP_FCG*)ksp->data; 417 418 PetscFunctionBegin; 419 PetscValidHeaderSpecific(ksp,KSP_CLASSID,1); 420 *nprealloc = fcg->nprealloc; 421 PetscFunctionReturn(0); 422 } 423 424 /*@ 425 KSPFCGSetTruncationType - specify how many of its stored previous directions FCG uses during orthoganalization 426 427 Logically Collective on ksp 428 429 KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions 430 KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1,.. 431 432 Input Parameters: 433 + ksp - the Krylov space context 434 - truncstrat - the choice of strategy 435 436 Level: intermediate 437 438 Options Database: 439 . -ksp_fcg_truncation_type <standard, notay> - specify how many of its stored previous directions FCG uses during orthoganalization 440 441 .seealso: `KSPFCDTruncationType`, `KSPFCGGetTruncationType` 442 @*/ 443 PetscErrorCode KSPFCGSetTruncationType(KSP ksp,KSPFCDTruncationType truncstrat) 444 { 445 KSP_FCG *fcg=(KSP_FCG*)ksp->data; 446 447 PetscFunctionBegin; 448 PetscValidHeaderSpecific(ksp,KSP_CLASSID,1); 449 PetscValidLogicalCollectiveEnum(ksp,truncstrat,2); 450 fcg->truncstrat=truncstrat; 451 PetscFunctionReturn(0); 452 } 453 454 /*@ 455 KSPFCGGetTruncationType - get the truncation strategy employed by FCG 456 457 Not Collective 458 459 Input Parameter: 460 . ksp - the Krylov space context 461 462 Output Parameter: 463 . truncstrat - the strategy type 464 465 Level: intermediate 466 467 .seealso: `KSPFCG`, `KSPFCGSetTruncationType`, `KSPFCDTruncationType` 468 @*/ 469 PetscErrorCode KSPFCGGetTruncationType(KSP ksp,KSPFCDTruncationType *truncstrat) 470 { 471 KSP_FCG *fcg=(KSP_FCG*)ksp->data; 472 473 PetscFunctionBegin; 474 PetscValidHeaderSpecific(ksp,KSP_CLASSID,1); 475 *truncstrat=fcg->truncstrat; 476 PetscFunctionReturn(0); 477 } 478 479 static PetscErrorCode KSPSetFromOptions_FCG(PetscOptionItems *PetscOptionsObject,KSP ksp) 480 { 481 KSP_FCG *fcg=(KSP_FCG*)ksp->data; 482 PetscInt mmax,nprealloc; 483 PetscBool flg; 484 485 PetscFunctionBegin; 486 PetscOptionsHeadBegin(PetscOptionsObject,"KSP FCG Options"); 487 PetscCall(PetscOptionsInt("-ksp_fcg_mmax","Maximum number of search directions to store","KSPFCGSetMmax",fcg->mmax,&mmax,&flg)); 488 if (flg) PetscCall(KSPFCGSetMmax(ksp,mmax)); 489 PetscCall(PetscOptionsInt("-ksp_fcg_nprealloc","Number of directions to preallocate","KSPFCGSetNprealloc",fcg->nprealloc,&nprealloc,&flg)); 490 if (flg) PetscCall(KSPFCGSetNprealloc(ksp,nprealloc)); 491 PetscCall(PetscOptionsEnum("-ksp_fcg_truncation_type","Truncation approach for directions","KSPFCGSetTruncationType",KSPFCDTruncationTypes,(PetscEnum)fcg->truncstrat,(PetscEnum*)&fcg->truncstrat,NULL)); 492 PetscOptionsHeadEnd(); 493 PetscFunctionReturn(0); 494 } 495 496 /*MC 497 KSPFCG - Implements the Flexible Conjugate Gradient method (FCG) 498 499 Options Database Keys: 500 + -ksp_fcg_mmax <N> - maximum number of search directions 501 . -ksp_fcg_nprealloc <N> - number of directions to preallocate 502 - -ksp_fcg_truncation_type <standard,notay> - truncation approach for directions 503 504 Contributed by Patrick Sanan 505 506 Notes: 507 Supports left preconditioning only. 508 509 Level: beginner 510 511 References: 512 + * - Notay, Y."Flexible Conjugate Gradients", SIAM J. Sci. Comput. 22:4, 2000 513 - * - Axelsson, O. and Vassilevski, P. S. "A Black Box Generalized Conjugate Gradient Solver with Inner Iterations and Variable step Preconditioning", 514 SIAM J. Matrix Anal. Appl. 12:4, 1991 515 516 .seealso `:` `KSPGCR`, `KSPFGMRES`, `KSPCG`, `KSPFCGSetMmax()`, `KSPFCGGetMmax()`, `KSPFCGSetNprealloc()`, `KSPFCGGetNprealloc()`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()` 517 518 M*/ 519 PETSC_EXTERN PetscErrorCode KSPCreate_FCG(KSP ksp) 520 { 521 KSP_FCG *fcg; 522 523 PetscFunctionBegin; 524 PetscCall(PetscNewLog(ksp,&fcg)); 525 #if !defined(PETSC_USE_COMPLEX) 526 fcg->type = KSP_CG_SYMMETRIC; 527 #else 528 fcg->type = KSP_CG_HERMITIAN; 529 #endif 530 fcg->mmax = KSPFCG_DEFAULT_MMAX; 531 fcg->nprealloc = KSPFCG_DEFAULT_NPREALLOC; 532 fcg->nvecs = 0; 533 fcg->vecb = KSPFCG_DEFAULT_VECB; 534 fcg->nchunks = 0; 535 fcg->truncstrat = KSPFCG_DEFAULT_TRUNCSTRAT; 536 537 ksp->data = (void*)fcg; 538 539 PetscCall(KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2)); 540 PetscCall(KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,1)); 541 PetscCall(KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,1)); 542 PetscCall(KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1)); 543 544 ksp->ops->setup = KSPSetUp_FCG; 545 ksp->ops->solve = KSPSolve_FCG; 546 ksp->ops->destroy = KSPDestroy_FCG; 547 ksp->ops->view = KSPView_FCG; 548 ksp->ops->setfromoptions = KSPSetFromOptions_FCG; 549 ksp->ops->buildsolution = KSPBuildSolutionDefault; 550 ksp->ops->buildresidual = KSPBuildResidualDefault; 551 PetscFunctionReturn(0); 552 } 553