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