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