1 /* 2 Contributed by Patrick Sanan and Sascha M. Schnepp 3 */ 4 5 #include <../src/ksp/ksp/impls/gmres/pipefgmres/pipefgmresimpl.h> /*I "petscksp.h" I*/ 6 7 static PetscBool cited = PETSC_FALSE; 8 static const char citation[] = 9 "@article{SSM2016,\n" 10 " author = {P. Sanan and S.M. Schnepp and D.A. May},\n" 11 " title = {Pipelined, Flexible Krylov Subspace Methods},\n" 12 " journal = {SIAM Journal on Scientific Computing},\n" 13 " volume = {38},\n" 14 " number = {5},\n" 15 " pages = {C441-C470},\n" 16 " year = {2016},\n" 17 " doi = {10.1137/15M1049130},\n" 18 " URL = {http://dx.doi.org/10.1137/15M1049130},\n" 19 " eprint = {http://dx.doi.org/10.1137/15M1049130}\n" 20 "}\n"; 21 22 #define PIPEFGMRES_DELTA_DIRECTIONS 10 23 #define PIPEFGMRES_DEFAULT_MAXK 30 24 25 static PetscErrorCode KSPPIPEFGMRESGetNewVectors(KSP,PetscInt); 26 static PetscErrorCode KSPPIPEFGMRESUpdateHessenberg(KSP,PetscInt,PetscBool*,PetscReal*); 27 static PetscErrorCode KSPPIPEFGMRESBuildSoln(PetscScalar*,Vec,Vec,KSP,PetscInt); 28 extern PetscErrorCode KSPReset_PIPEFGMRES(KSP); 29 30 /* 31 32 KSPSetUp_PIPEFGMRES - Sets up the workspace needed by pipefgmres. 33 34 This is called once, usually automatically by KSPSolve() or KSPSetUp(), 35 but can be called directly by KSPSetUp(). 36 37 */ 38 static PetscErrorCode KSPSetUp_PIPEFGMRES(KSP ksp) 39 { 40 PetscErrorCode ierr; 41 PetscInt k; 42 KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data; 43 const PetscInt max_k = pipefgmres->max_k; 44 45 PetscFunctionBegin; 46 ierr = KSPSetUp_GMRES(ksp);CHKERRQ(ierr); 47 48 ierr = PetscMalloc1((VEC_OFFSET+max_k),&pipefgmres->prevecs);CHKERRQ(ierr); 49 ierr = PetscMalloc1((VEC_OFFSET+max_k),&pipefgmres->prevecs_user_work);CHKERRQ(ierr); 50 ierr = PetscLogObjectMemory((PetscObject)ksp,(VEC_OFFSET+max_k)*(2*sizeof(void*)));CHKERRQ(ierr); 51 52 ierr = KSPCreateVecs(ksp,pipefgmres->vv_allocated,&pipefgmres->prevecs_user_work[0],0,NULL);CHKERRQ(ierr); 53 ierr = PetscLogObjectParents(ksp,pipefgmres->vv_allocated,pipefgmres->prevecs_user_work[0]);CHKERRQ(ierr); 54 for (k=0; k < pipefgmres->vv_allocated; k++) { 55 pipefgmres->prevecs[k] = pipefgmres->prevecs_user_work[0][k]; 56 } 57 58 ierr = PetscMalloc1((VEC_OFFSET+max_k),&pipefgmres->zvecs);CHKERRQ(ierr); 59 ierr = PetscMalloc1((VEC_OFFSET+max_k),&pipefgmres->zvecs_user_work);CHKERRQ(ierr); 60 ierr = PetscLogObjectMemory((PetscObject)ksp,(VEC_OFFSET+max_k)*(2*sizeof(void*)));CHKERRQ(ierr); 61 62 ierr = PetscMalloc1((VEC_OFFSET+max_k),&pipefgmres->redux);CHKERRQ(ierr); 63 ierr = PetscLogObjectMemory((PetscObject)ksp,(VEC_OFFSET+max_k)*(sizeof(void*)));CHKERRQ(ierr); 64 65 ierr = KSPCreateVecs(ksp,pipefgmres->vv_allocated,&pipefgmres->zvecs_user_work[0],0,NULL);CHKERRQ(ierr); 66 ierr = PetscLogObjectParents(ksp,pipefgmres->vv_allocated,pipefgmres->zvecs_user_work[0]);CHKERRQ(ierr); 67 for (k=0; k < pipefgmres->vv_allocated; k++) { 68 pipefgmres->zvecs[k] = pipefgmres->zvecs_user_work[0][k]; 69 } 70 71 PetscFunctionReturn(0); 72 } 73 74 /* 75 76 KSPPIPEFGMRESCycle - Run pipefgmres, possibly with restart. Return residual 77 history if requested. 78 79 input parameters: 80 . pipefgmres - structure containing parameters and work areas 81 82 output parameters: 83 . itcount - number of iterations used. If null, ignored. 84 . converged - 0 if not converged 85 86 Notes: 87 On entry, the value in vector VEC_VV(0) should be 88 the initial residual. 89 90 91 */ 92 93 static PetscErrorCode KSPPIPEFGMRESCycle(PetscInt *itcount,KSP ksp) 94 { 95 KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)(ksp->data); 96 PetscReal res_norm; 97 PetscReal hapbnd,tt; 98 PetscScalar *hh,*hes,*lhh,shift = pipefgmres->shift; 99 PetscBool hapend = PETSC_FALSE; /* indicates happy breakdown ending */ 100 PetscErrorCode ierr; 101 PetscInt loc_it; /* local count of # of dir. in Krylov space */ 102 PetscInt max_k = pipefgmres->max_k; /* max # of directions Krylov space */ 103 PetscInt i,j,k; 104 Mat Amat,Pmat; 105 Vec Q,W; /* Pipelining vectors */ 106 Vec *redux = pipefgmres->redux; /* workspace for single reduction */ 107 108 PetscFunctionBegin; 109 if (itcount) *itcount = 0; 110 111 /* Assign simpler names to these vectors, allocated as pipelining workspace */ 112 Q = VEC_Q; 113 W = VEC_W; 114 115 /* Allocate memory for orthogonalization work (freed in the GMRES Destroy routine)*/ 116 /* Note that we add an extra value here to allow for a single reduction */ 117 if (!pipefgmres->orthogwork) { ierr = PetscMalloc1(pipefgmres->max_k + 2 ,&pipefgmres->orthogwork);CHKERRQ(ierr); 118 } 119 lhh = pipefgmres->orthogwork; 120 121 /* Number of pseudo iterations since last restart is the number 122 of prestart directions */ 123 loc_it = 0; 124 125 /* note: (pipefgmres->it) is always set one less than (loc_it) It is used in 126 KSPBUILDSolution_PIPEFGMRES, where it is passed to KSPPIPEFGMRESBuildSoln. 127 Note that when KSPPIPEFGMRESBuildSoln is called from this function, 128 (loc_it -1) is passed, so the two are equivalent */ 129 pipefgmres->it = (loc_it - 1); 130 131 /* initial residual is in VEC_VV(0) - compute its norm*/ 132 ierr = VecNorm(VEC_VV(0),NORM_2,&res_norm);CHKERRQ(ierr); 133 134 /* first entry in right-hand-side of hessenberg system is just 135 the initial residual norm */ 136 *RS(0) = res_norm; 137 138 ksp->rnorm = res_norm; 139 ierr = KSPLogResidualHistory(ksp,res_norm);CHKERRQ(ierr); 140 ierr = KSPMonitor(ksp,ksp->its,res_norm);CHKERRQ(ierr); 141 142 /* check for the convergence - maybe the current guess is good enough */ 143 ierr = (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); 144 if (ksp->reason) { 145 if (itcount) *itcount = 0; 146 PetscFunctionReturn(0); 147 } 148 149 /* scale VEC_VV (the initial residual) */ 150 ierr = VecScale(VEC_VV(0),1.0/res_norm);CHKERRQ(ierr); 151 152 /* Fill the pipeline */ 153 ierr = KSP_PCApply(ksp,VEC_VV(loc_it),PREVEC(loc_it));CHKERRQ(ierr); 154 ierr = PCGetOperators(ksp->pc,&Amat,&Pmat);CHKERRQ(ierr); 155 ierr = KSP_MatMult(ksp,Amat,PREVEC(loc_it),ZVEC(loc_it));CHKERRQ(ierr); 156 ierr = VecAXPY(ZVEC(loc_it),-shift,VEC_VV(loc_it));CHKERRQ(ierr); /* Note shift */ 157 158 /* MAIN ITERATION LOOP BEGINNING*/ 159 /* keep iterating until we have converged OR generated the max number 160 of directions OR reached the max number of iterations for the method */ 161 while (!ksp->reason && loc_it < max_k && ksp->its < ksp->max_it) { 162 if (loc_it) { 163 ierr = KSPLogResidualHistory(ksp,res_norm);CHKERRQ(ierr); 164 ierr = KSPMonitor(ksp,ksp->its,res_norm);CHKERRQ(ierr); 165 } 166 pipefgmres->it = (loc_it - 1); 167 168 /* see if more space is needed for work vectors */ 169 if (pipefgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) { 170 ierr = KSPPIPEFGMRESGetNewVectors(ksp,loc_it+1);CHKERRQ(ierr); 171 /* (loc_it+1) is passed in as number of the first vector that should 172 be allocated */ 173 } 174 175 /* Note that these inner products are with "Z" now, so 176 in particular, lhh[loc_it] is the 'barred' or 'shifted' value, 177 not the value from the equivalent FGMRES run (even in exact arithmetic) 178 That is, the H we need for the Arnoldi relation is different from the 179 coefficients we use in the orthogonalization process,because of the shift */ 180 181 /* Do some local twiddling to allow for a single reduction */ 182 for(i=0;i<loc_it+1;i++){ 183 redux[i] = VEC_VV(i); 184 } 185 redux[loc_it+1] = ZVEC(loc_it); 186 187 /* note the extra dot product which ends up in lh[loc_it+1], which computes ||z||^2 */ 188 ierr = VecMDotBegin(ZVEC(loc_it),loc_it+2,redux,lhh);CHKERRQ(ierr); 189 190 /* Start the split reduction (This actually calls the MPI_Iallreduce, otherwise, the reduction is simply delayed until the "end" call)*/ 191 ierr = PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)ZVEC(loc_it)));CHKERRQ(ierr); 192 193 /* The work to be overlapped with the inner products follows. 194 This is application of the preconditioner and the operator 195 to compute intermediate quantites which will be combined (locally) 196 with the results of the inner products. 197 */ 198 ierr = KSP_PCApply(ksp,ZVEC(loc_it),Q);CHKERRQ(ierr); 199 ierr = PCGetOperators(ksp->pc,&Amat,&Pmat);CHKERRQ(ierr); 200 ierr = KSP_MatMult(ksp,Amat,Q,W);CHKERRQ(ierr); 201 202 /* Compute inner products of the new direction with previous directions, 203 and the norm of the to-be-orthogonalized direction "Z". 204 This information is enough to build the required entries 205 of H. The inner product with VEC_VV(it_loc) is 206 *different* than in the standard FGMRES and need to be dealt with specially. 207 That is, for standard FGMRES the orthogonalization coefficients are the same 208 as the coefficients used in the Arnoldi relation to reconstruct, but here this 209 is not true (albeit only for the one entry of H which we "unshift" below. */ 210 211 /* Finish the dot product, retrieving the extra entry */ 212 ierr = VecMDotEnd(ZVEC(loc_it),loc_it+2,redux,lhh);CHKERRQ(ierr); 213 tt = PetscRealPart(lhh[loc_it+1]); 214 215 /* Hessenberg entries, and entries for (naive) classical Graham-Schmidt 216 Note that the Hessenberg entries require a shift, as these are for the 217 relation AU = VH, which is wrt unshifted basis vectors */ 218 hh = HH(0,loc_it); hes=HES(0,loc_it); 219 for (j=0; j<loc_it; j++) { 220 hh[j] = lhh[j]; 221 hes[j] = lhh[j]; 222 } 223 hh[loc_it] = lhh[loc_it] + shift; 224 hes[loc_it] = lhh[loc_it] + shift; 225 226 /* we delay applying the shift here */ 227 for (j=0; j<=loc_it; j++) { 228 lhh[j] = -lhh[j]; /* flip sign */ 229 } 230 231 /* Compute the norm of the un-normalized new direction using the rearranged formula 232 Note that these are shifted ("barred") quantities */ 233 for(k=0;k<=loc_it;k++) tt -= ((PetscReal)(PetscAbsScalar(lhh[k]) * PetscAbsScalar(lhh[k]))); 234 if (tt < 0.0) { 235 /* If we detect square root breakdown in the norm, we must restart the algorithm. 236 Here this means we simply break the current loop and reconstruct the solution 237 using the basis we have computed thus far. Note that by breaking immediately, 238 we do not update the iteration count, so computation done in this iteration 239 should be disregarded. 240 */ 241 ierr = PetscInfo1(ksp,"Restart due to square root breakdown at it = \n",ksp->its);CHKERRQ(ierr); 242 break; 243 } else { 244 tt = PetscSqrtReal(tt); 245 } 246 247 /* new entry in hessenburg is the 2-norm of our new direction */ 248 hh[loc_it+1] = tt; 249 hes[loc_it+1] = tt; 250 251 /* The recurred computation for the new direction 252 The division by tt is delayed to the happy breakdown check later 253 Note placement BEFORE the unshift 254 */ 255 ierr = VecCopy(ZVEC(loc_it),VEC_VV(loc_it+1));CHKERRQ(ierr); 256 ierr = VecMAXPY(VEC_VV(loc_it+1),loc_it+1,lhh,&VEC_VV(0));CHKERRQ(ierr); 257 /* (VEC_VV(loc_it+1) is not normalized yet) */ 258 259 /* The recurred computation for the preconditioned vector (u) */ 260 ierr = VecCopy(Q,PREVEC(loc_it+1));CHKERRQ(ierr); 261 ierr = VecMAXPY(PREVEC(loc_it+1),loc_it+1,lhh,&PREVEC(0));CHKERRQ(ierr); 262 ierr = VecScale(PREVEC(loc_it+1),1.0/tt);CHKERRQ(ierr); 263 264 /* Unshift an entry in the GS coefficients ("removing the bar") */ 265 lhh[loc_it] -= shift; 266 267 /* The recurred computation for z (Au) 268 Note placement AFTER the "unshift" */ 269 ierr = VecCopy(W,ZVEC(loc_it+1));CHKERRQ(ierr); 270 ierr = VecMAXPY(ZVEC(loc_it+1),loc_it+1,lhh,&ZVEC(0));CHKERRQ(ierr); 271 ierr = VecScale(ZVEC(loc_it+1),1.0/tt);CHKERRQ(ierr); 272 273 /* Happy Breakdown Check */ 274 hapbnd = PetscAbsScalar((tt) / *RS(loc_it)); 275 /* RS(loc_it) contains the res_norm from the last iteration */ 276 hapbnd = PetscMin(pipefgmres->haptol,hapbnd); 277 if (tt > hapbnd) { 278 /* scale new direction by its norm */ 279 ierr = VecScale(VEC_VV(loc_it+1),1.0/tt);CHKERRQ(ierr); 280 } else { 281 /* This happens when the solution is exactly reached. */ 282 /* So there is no new direction... */ 283 ierr = VecSet(VEC_TEMP,0.0);CHKERRQ(ierr); /* set VEC_TEMP to 0 */ 284 hapend = PETSC_TRUE; 285 } 286 /* note that for pipefgmres we could get HES(loc_it+1, loc_it) = 0 and the 287 current solution would not be exact if HES was singular. Note that 288 HH non-singular implies that HES is not singular, and HES is guaranteed 289 to be nonsingular when PREVECS are linearly independent and A is 290 nonsingular (in GMRES, the nonsingularity of A implies the nonsingularity 291 of HES). So we should really add a check to verify that HES is nonsingular.*/ 292 293 /* Note that to be thorough, in debug mode, one could call a LAPACK routine 294 here to check that the hessenberg matrix is indeed non-singular (since 295 FGMRES does not guarantee this) */ 296 297 /* Now apply rotations to new col of hessenberg (and right side of system), 298 calculate new rotation, and get new residual norm at the same time*/ 299 ierr = KSPPIPEFGMRESUpdateHessenberg(ksp,loc_it,&hapend,&res_norm);CHKERRQ(ierr); 300 if (ksp->reason) break; 301 302 loc_it++; 303 pipefgmres->it = (loc_it-1); /* Add this here in case it has converged */ 304 305 ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); 306 ksp->its++; 307 ksp->rnorm = res_norm; 308 ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); 309 310 ierr = (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); 311 312 /* Catch error in happy breakdown and signal convergence and break from loop */ 313 if (hapend) { 314 if (!ksp->reason) { 315 if (ksp->errorifnotconverged) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"You reached the happy break down, but convergence was not indicated. Residual norm = %g",(double)res_norm); 316 else { 317 ksp->reason = KSP_DIVERGED_BREAKDOWN; 318 break; 319 } 320 } 321 } 322 } 323 /* END OF ITERATION LOOP */ 324 ierr = KSPLogResidualHistory(ksp,res_norm);CHKERRQ(ierr); 325 326 /* 327 Monitor if we know that we will not return for a restart */ 328 if (loc_it && (ksp->reason || ksp->its >= ksp->max_it)) { 329 ierr = KSPMonitor(ksp,ksp->its,res_norm);CHKERRQ(ierr); 330 } 331 332 if (itcount) *itcount = loc_it; 333 334 /* 335 Down here we have to solve for the "best" coefficients of the Krylov 336 columns, add the solution values together, and possibly unwind the 337 preconditioning from the solution 338 */ 339 340 /* Form the solution (or the solution so far) */ 341 /* Note: must pass in (loc_it-1) for iteration count so that KSPPIPEGMRESIIBuildSoln 342 properly navigates */ 343 344 ierr = KSPPIPEFGMRESBuildSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);CHKERRQ(ierr); 345 346 PetscFunctionReturn(0); 347 } 348 349 /* 350 KSPSolve_PIPEFGMRES - This routine applies the PIPEFGMRES method. 351 352 353 Input Parameter: 354 . ksp - the Krylov space object that was set to use pipefgmres 355 356 Output Parameter: 357 . outits - number of iterations used 358 359 */ 360 static PetscErrorCode KSPSolve_PIPEFGMRES(KSP ksp) 361 { 362 PetscErrorCode ierr; 363 PetscInt its,itcount; 364 KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data; 365 PetscBool guess_zero = ksp->guess_zero; 366 367 PetscFunctionBegin; 368 369 /* We have not checked these routines for use with complex numbers. The inner products 370 are likely not defined correctly for that case */ 371 #if (defined(PETSC_USE_COMPLEX) && !defined(PETSC_SKIP_COMPLEX)) 372 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"PIPEFGMRES has not been implemented for use with complex scalars"); 373 #endif 374 375 ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 376 377 if (ksp->calc_sings && !pipefgmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called"); 378 ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); 379 ksp->its = 0; 380 ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); 381 382 itcount = 0; 383 ksp->reason = KSP_CONVERGED_ITERATING; 384 while (!ksp->reason) { 385 ierr = KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);CHKERRQ(ierr); 386 ierr = KSPPIPEFGMRESCycle(&its,ksp);CHKERRQ(ierr); 387 itcount += its; 388 if (itcount >= ksp->max_it) { 389 if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS; 390 break; 391 } 392 ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */ 393 } 394 ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */ 395 PetscFunctionReturn(0); 396 } 397 398 static PetscErrorCode KSPDestroy_PIPEFGMRES(KSP ksp) 399 { 400 PetscErrorCode ierr; 401 402 PetscFunctionBegin; 403 ierr = KSPReset_PIPEFGMRES(ksp);CHKERRQ(ierr); 404 ierr = KSPDestroy_GMRES(ksp);CHKERRQ(ierr); 405 PetscFunctionReturn(0); 406 } 407 408 /* 409 KSPPIPEFGMRESBuildSoln - create the solution from the starting vector and the 410 current iterates. 411 412 Input parameters: 413 nrs - work area of size it + 1. 414 vguess - index of initial guess 415 vdest - index of result. Note that vguess may == vdest (replace 416 guess with the solution). 417 it - HH upper triangular part is a block of size (it+1) x (it+1) 418 419 This is an internal routine that knows about the PIPEFGMRES internals. 420 */ 421 static PetscErrorCode KSPPIPEFGMRESBuildSoln(PetscScalar *nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it) 422 { 423 PetscScalar tt; 424 PetscErrorCode ierr; 425 PetscInt k,j; 426 KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)(ksp->data); 427 428 PetscFunctionBegin; 429 /* Solve for solution vector that minimizes the residual */ 430 431 if (it < 0) { /* no pipefgmres steps have been performed */ 432 ierr = VecCopy(vguess,vdest);CHKERRQ(ierr); /* VecCopy() is smart, exits immediately if vguess == vdest */ 433 PetscFunctionReturn(0); 434 } 435 436 /* solve the upper triangular system - RS is the right side and HH is 437 the upper triangular matrix - put soln in nrs */ 438 if (*HH(it,it) != 0.0) nrs[it] = *RS(it) / *HH(it,it); 439 else nrs[it] = 0.0; 440 441 for (k=it-1; k>=0; k--) { 442 tt = *RS(k); 443 for (j=k+1; j<=it; j++) tt -= *HH(k,j) * nrs[j]; 444 nrs[k] = tt / *HH(k,k); 445 } 446 447 /* Accumulate the correction to the solution of the preconditioned problem in VEC_TEMP */ 448 ierr = VecZeroEntries(VEC_TEMP);CHKERRQ(ierr); 449 ierr = VecMAXPY(VEC_TEMP,it+1,nrs,&PREVEC(0));CHKERRQ(ierr); 450 451 /* add solution to previous solution */ 452 if (vdest == vguess) { 453 ierr = VecAXPY(vdest,1.0,VEC_TEMP);CHKERRQ(ierr); 454 } else { 455 ierr = VecWAXPY(vdest,1.0,VEC_TEMP,vguess);CHKERRQ(ierr); 456 } 457 PetscFunctionReturn(0); 458 } 459 460 /* 461 462 KSPPIPEFGMRESUpdateHessenberg - Do the scalar work for the orthogonalization. 463 Return new residual. 464 465 input parameters: 466 467 . ksp - Krylov space object 468 . it - plane rotations are applied to the (it+1)th column of the 469 modified hessenberg (i.e. HH(:,it)) 470 . hapend - PETSC_FALSE not happy breakdown ending. 471 472 output parameters: 473 . res - the new residual 474 475 */ 476 /* 477 . it - column of the Hessenberg that is complete, PIPEFGMRES is actually computing two columns ahead of this 478 */ 479 static PetscErrorCode KSPPIPEFGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool *hapend,PetscReal *res) 480 { 481 PetscScalar *hh,*cc,*ss,*rs; 482 PetscInt j; 483 PetscReal hapbnd; 484 KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)(ksp->data); 485 PetscErrorCode ierr; 486 487 PetscFunctionBegin; 488 hh = HH(0,it); /* pointer to beginning of column to update */ 489 cc = CC(0); /* beginning of cosine rotations */ 490 ss = SS(0); /* beginning of sine rotations */ 491 rs = RS(0); /* right hand side of least squares system */ 492 493 /* The Hessenberg matrix is now correct through column it, save that form for possible spectral analysis */ 494 for (j=0; j<=it+1; j++) *HES(j,it) = hh[j]; 495 496 /* check for the happy breakdown */ 497 hapbnd = PetscMin(PetscAbsScalar(hh[it+1] / rs[it]),pipefgmres->haptol); 498 if (PetscAbsScalar(hh[it+1]) < hapbnd) { 499 ierr = PetscInfo4(ksp,"Detected happy breakdown, current hapbnd = %14.12e H(%D,%D) = %14.12e\n",(double)hapbnd,it+1,it,(double)PetscAbsScalar(*HH(it+1,it)));CHKERRQ(ierr); 500 *hapend = PETSC_TRUE; 501 } 502 503 /* Apply all the previously computed plane rotations to the new column 504 of the Hessenberg matrix */ 505 /* Note: this uses the rotation [conj(c) s ; -s c], c= cos(theta), s= sin(theta), 506 and some refs have [c s ; -conj(s) c] (don't be confused!) */ 507 508 for (j=0; j<it; j++) { 509 PetscScalar hhj = hh[j]; 510 hh[j] = PetscConj(cc[j])*hhj + ss[j]*hh[j+1]; 511 hh[j+1] = -ss[j] *hhj + cc[j]*hh[j+1]; 512 } 513 514 /* 515 compute the new plane rotation, and apply it to: 516 1) the right-hand-side of the Hessenberg system (RS) 517 note: it affects RS(it) and RS(it+1) 518 2) the new column of the Hessenberg matrix 519 note: it affects HH(it,it) which is currently pointed to 520 by hh and HH(it+1, it) (*(hh+1)) 521 thus obtaining the updated value of the residual... 522 */ 523 524 /* compute new plane rotation */ 525 526 if (!*hapend) { 527 PetscReal delta = PetscSqrtReal(PetscSqr(PetscAbsScalar(hh[it])) + PetscSqr(PetscAbsScalar(hh[it+1]))); 528 if (delta == 0.0) { 529 ksp->reason = KSP_DIVERGED_NULL; 530 PetscFunctionReturn(0); 531 } 532 533 cc[it] = hh[it] / delta; /* new cosine value */ 534 ss[it] = hh[it+1] / delta; /* new sine value */ 535 536 hh[it] = PetscConj(cc[it])*hh[it] + ss[it]*hh[it+1]; 537 rs[it+1] = -ss[it]*rs[it]; 538 rs[it] = PetscConj(cc[it])*rs[it]; 539 *res = PetscAbsScalar(rs[it+1]); 540 } else { /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply 541 another rotation matrix (so RH doesn't change). The new residual is 542 always the new sine term times the residual from last time (RS(it)), 543 but now the new sine rotation would be zero...so the residual should 544 be zero...so we will multiply "zero" by the last residual. This might 545 not be exactly what we want to do here -could just return "zero". */ 546 547 *res = 0.0; 548 } 549 PetscFunctionReturn(0); 550 } 551 552 /* 553 KSPBuildSolution_PIPEFGMRES 554 555 Input Parameter: 556 . ksp - the Krylov space object 557 . ptr- 558 559 Output Parameter: 560 . result - the solution 561 562 Note: this calls KSPPIPEFGMRESBuildSoln - the same function that KSPPIPEFGMRESCycle 563 calls directly. 564 565 */ 566 PetscErrorCode KSPBuildSolution_PIPEFGMRES(KSP ksp,Vec ptr,Vec *result) 567 { 568 KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data; 569 PetscErrorCode ierr; 570 571 PetscFunctionBegin; 572 if (!ptr) { 573 if (!pipefgmres->sol_temp) { 574 ierr = VecDuplicate(ksp->vec_sol,&pipefgmres->sol_temp);CHKERRQ(ierr); 575 ierr = PetscLogObjectParent((PetscObject)ksp,(PetscObject)pipefgmres->sol_temp);CHKERRQ(ierr); 576 } 577 ptr = pipefgmres->sol_temp; 578 } 579 if (!pipefgmres->nrs) { 580 /* allocate the work area */ 581 ierr = PetscMalloc1(pipefgmres->max_k,&pipefgmres->nrs);CHKERRQ(ierr); 582 ierr = PetscLogObjectMemory((PetscObject)ksp,pipefgmres->max_k*sizeof(PetscScalar));CHKERRQ(ierr); 583 } 584 585 ierr = KSPPIPEFGMRESBuildSoln(pipefgmres->nrs,ksp->vec_sol,ptr,ksp,pipefgmres->it);CHKERRQ(ierr); 586 if (result) *result = ptr; 587 PetscFunctionReturn(0); 588 } 589 590 PetscErrorCode KSPSetFromOptions_PIPEFGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp) 591 { 592 PetscErrorCode ierr; 593 KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data; 594 PetscBool flg; 595 PetscScalar shift; 596 597 PetscFunctionBegin; 598 ierr = KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);CHKERRQ(ierr); 599 ierr = PetscOptionsHead(PetscOptionsObject,"KSP pipelined FGMRES Options");CHKERRQ(ierr); 600 ierr = PetscOptionsScalar("-ksp_pipefgmres_shift","shift parameter","KSPPIPEFGMRESSetShift",pipefgmres->shift,&shift,&flg);CHKERRQ(ierr); 601 if (flg) { ierr = KSPPIPEFGMRESSetShift(ksp,shift);CHKERRQ(ierr); } 602 ierr = PetscOptionsTail();CHKERRQ(ierr); 603 PetscFunctionReturn(0); 604 } 605 606 PetscErrorCode KSPView_PIPEFGMRES(KSP ksp,PetscViewer viewer) 607 { 608 KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data; 609 PetscErrorCode ierr; 610 PetscBool iascii,isstring; 611 612 PetscFunctionBegin; 613 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 614 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 615 616 if (iascii) { 617 ierr = PetscViewerASCIIPrintf(viewer," restart=%D\n",pipefgmres->max_k);CHKERRQ(ierr); 618 ierr = PetscViewerASCIIPrintf(viewer," happy breakdown tolerance %g\n",(double)pipefgmres->haptol);CHKERRQ(ierr); 619 #if defined(PETSC_USE_COMPLEX) 620 ierr = PetscViewerASCIIPrintf(viewer," shift=%g+%gi\n",PetscRealPart(pipefgmres->shift),PetscImaginaryPart(pipefgmres->shift));CHKERRQ(ierr); 621 #else 622 ierr = PetscViewerASCIIPrintf(viewer," shift=%g\n",pipefgmres->shift);CHKERRQ(ierr); 623 #endif 624 } else if (isstring) { 625 ierr = PetscViewerStringSPrintf(viewer,"restart %D",pipefgmres->max_k);CHKERRQ(ierr); 626 #if defined(PETSC_USE_COMPLEX) 627 ierr = PetscViewerStringSPrintf(viewer," shift=%g+%gi\n",PetscRealPart(pipefgmres->shift),PetscImaginaryPart(pipefgmres->shift));CHKERRQ(ierr); 628 #else 629 ierr = PetscViewerStringSPrintf(viewer," shift=%g\n",pipefgmres->shift);CHKERRQ(ierr); 630 #endif 631 } 632 PetscFunctionReturn(0); 633 } 634 635 PetscErrorCode KSPReset_PIPEFGMRES(KSP ksp) 636 { 637 KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data; 638 PetscErrorCode ierr; 639 PetscInt i; 640 641 PetscFunctionBegin; 642 ierr = PetscFree(pipefgmres->prevecs);CHKERRQ(ierr); 643 ierr = PetscFree(pipefgmres->zvecs);CHKERRQ(ierr); 644 for (i=0; i<pipefgmres->nwork_alloc; i++) { 645 ierr = VecDestroyVecs(pipefgmres->mwork_alloc[i],&pipefgmres->prevecs_user_work[i]);CHKERRQ(ierr); 646 ierr = VecDestroyVecs(pipefgmres->mwork_alloc[i],&pipefgmres->zvecs_user_work[i]);CHKERRQ(ierr); 647 } 648 ierr = PetscFree(pipefgmres->prevecs_user_work);CHKERRQ(ierr); 649 ierr = PetscFree(pipefgmres->zvecs_user_work);CHKERRQ(ierr); 650 ierr = PetscFree(pipefgmres->redux);CHKERRQ(ierr); 651 ierr = KSPReset_GMRES(ksp);CHKERRQ(ierr); 652 PetscFunctionReturn(0); 653 } 654 655 /*MC 656 KSPPIPEFGMRES - Implements the Pipelined Generalized Minimal Residual method. 657 658 A flexible, 1-stage pipelined variant of GMRES. 659 660 Options Database Keys: 661 + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against 662 . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence) 663 . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of 664 . -ksp_pipefgmres_shift - the shift to use (defaults to 1. See KSPPIPEFGMRESSetShift() 665 vectors are allocated as needed) 666 - -ksp_gmres_krylov_monitor - plot the Krylov space generated 667 668 669 Level: intermediate 670 671 Notes: 672 673 This variant is not "explicitly normalized" like KSPPGMRES, and requires a shift parameter. 674 675 A heuristic for choosing the shift parameter is the largest eigenvalue of the preconditioned operator. 676 677 Only right preconditioning is supported (but this preconditioner may be nonlinear/variable/inexact, as with KSPFGMRES). 678 MPI configuration may be necessary for reductions to make asynchronous progress, which is important for performance of pipelined methods. 679 See the FAQ on the PETSc website for details. 680 681 Developer Notes: This class is subclassed off of KSPGMRES. 682 683 Reference: 684 P. Sanan, S.M. Schnepp, and D.A. May, 685 "Pipelined, Flexible Krylov Subspace Methods," 686 SIAM Journal on Scientific Computing 2016 38:5, C441-C470, 687 DOI: 10.1137/15M1049130 688 689 .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPLGMRES, KSPPIPECG, KSPPIPECR, KSPPGMRES, KSPFGMRES 690 KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESMonitorKrylov(), KSPPIPEFGMRESSetShift() 691 M*/ 692 693 PETSC_EXTERN PetscErrorCode KSPCreate_PIPEFGMRES(KSP ksp) 694 { 695 KSP_PIPEFGMRES *pipefgmres; 696 PetscErrorCode ierr; 697 698 PetscFunctionBegin; 699 ierr = PetscNewLog(ksp,&pipefgmres);CHKERRQ(ierr); 700 701 ksp->data = (void*)pipefgmres; 702 ksp->ops->buildsolution = KSPBuildSolution_PIPEFGMRES; 703 ksp->ops->setup = KSPSetUp_PIPEFGMRES; 704 ksp->ops->solve = KSPSolve_PIPEFGMRES; 705 ksp->ops->reset = KSPReset_PIPEFGMRES; 706 ksp->ops->destroy = KSPDestroy_PIPEFGMRES; 707 ksp->ops->view = KSPView_PIPEFGMRES; 708 ksp->ops->setfromoptions = KSPSetFromOptions_PIPEFGMRES; 709 ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES; 710 ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES; 711 712 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,3);CHKERRQ(ierr); 713 714 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);CHKERRQ(ierr); 715 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);CHKERRQ(ierr); 716 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);CHKERRQ(ierr); 717 718 pipefgmres->nextra_vecs = 1; 719 pipefgmres->haptol = 1.0e-30; 720 pipefgmres->q_preallocate = 0; 721 pipefgmres->delta_allocate = PIPEFGMRES_DELTA_DIRECTIONS; 722 pipefgmres->orthog = 0; 723 pipefgmres->nrs = 0; 724 pipefgmres->sol_temp = 0; 725 pipefgmres->max_k = PIPEFGMRES_DEFAULT_MAXK; 726 pipefgmres->Rsvd = 0; 727 pipefgmres->orthogwork = 0; 728 pipefgmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER; 729 pipefgmres->shift = 1.0; 730 PetscFunctionReturn(0); 731 } 732 733 static PetscErrorCode KSPPIPEFGMRESGetNewVectors(KSP ksp,PetscInt it) 734 { 735 KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data; 736 PetscInt nwork = pipefgmres->nwork_alloc; /* number of work vector chunks allocated */ 737 PetscInt nalloc; /* number to allocate */ 738 PetscErrorCode ierr; 739 PetscInt k; 740 741 PetscFunctionBegin; 742 nalloc = pipefgmres->delta_allocate; /* number of vectors to allocate 743 in a single chunk */ 744 745 /* Adjust the number to allocate to make sure that we don't exceed the 746 number of available slots (pipefgmres->vecs_allocated)*/ 747 if (it + VEC_OFFSET + nalloc >= pipefgmres->vecs_allocated) { 748 nalloc = pipefgmres->vecs_allocated - it - VEC_OFFSET; 749 } 750 if (!nalloc) PetscFunctionReturn(0); 751 752 pipefgmres->vv_allocated += nalloc; /* vv_allocated is the number of vectors allocated */ 753 754 /* work vectors */ 755 ierr = KSPCreateVecs(ksp,nalloc,&pipefgmres->user_work[nwork],0,NULL);CHKERRQ(ierr); 756 ierr = PetscLogObjectParents(ksp,nalloc,pipefgmres->user_work[nwork]);CHKERRQ(ierr); 757 for (k=0; k < nalloc; k++) { 758 pipefgmres->vecs[it+VEC_OFFSET+k] = pipefgmres->user_work[nwork][k]; 759 } 760 /* specify size of chunk allocated */ 761 pipefgmres->mwork_alloc[nwork] = nalloc; 762 763 /* preconditioned vectors (note we don't use VEC_OFFSET) */ 764 ierr = KSPCreateVecs(ksp,nalloc,&pipefgmres->prevecs_user_work[nwork],0,NULL);CHKERRQ(ierr); 765 ierr = PetscLogObjectParents(ksp,nalloc,pipefgmres->prevecs_user_work[nwork]);CHKERRQ(ierr); 766 for (k=0; k < nalloc; k++) { 767 pipefgmres->prevecs[it+k] = pipefgmres->prevecs_user_work[nwork][k]; 768 } 769 770 ierr = KSPCreateVecs(ksp,nalloc,&pipefgmres->zvecs_user_work[nwork],0,NULL);CHKERRQ(ierr); 771 ierr = PetscLogObjectParents(ksp,nalloc,pipefgmres->zvecs_user_work[nwork]);CHKERRQ(ierr); 772 for (k=0; k < nalloc; k++) { 773 pipefgmres->zvecs[it+k] = pipefgmres->zvecs_user_work[nwork][k]; 774 } 775 776 /* increment the number of work vector chunks */ 777 pipefgmres->nwork_alloc++; 778 PetscFunctionReturn(0); 779 } 780 /*@ 781 KSPPIPEFGMRESSetShift - Set the shift parameter for the flexible, pipelined GMRES solver. 782 783 A heuristic is to set this to be comparable to the largest eigenvalue of the preconditioned operator. This can be acheived with PETSc itself by using a few iterations of a Krylov method. See KSPComputeEigenvalues (and note the caveats there). 784 785 Logically Collective on KSP 786 787 Input Parameters: 788 + ksp - the Krylov space context 789 - shift - the shift 790 791 Level: intermediate 792 793 Options Database: 794 . -ksp_pipefgmres_shift <shift> 795 796 .seealso: KSPComputeEigenvalues() 797 @*/ 798 PetscErrorCode KSPPIPEFGMRESSetShift(KSP ksp,PetscScalar shift) 799 { 800 KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data; 801 802 PetscFunctionBegin; 803 PetscValidHeaderSpecific(ksp,KSP_CLASSID,1); 804 PetscValidLogicalCollectiveScalar(ksp,shift,2); 805 pipefgmres->shift = shift; 806 PetscFunctionReturn(0); 807 } 808