1 2 /* 3 This file implements PGMRES (a Pipelined Generalized Minimal Residual method) 4 */ 5 6 #include <../src/ksp/ksp/impls/gmres/pgmres/pgmresimpl.h> /*I "petscksp.h" I*/ 7 #define PGMRES_DELTA_DIRECTIONS 10 8 #define PGMRES_DEFAULT_MAXK 30 9 10 static PetscErrorCode KSPPGMRESUpdateHessenberg(KSP,PetscInt,PetscBool*,PetscReal*); 11 static PetscErrorCode KSPPGMRESBuildSoln(PetscScalar*,Vec,Vec,KSP,PetscInt); 12 13 /* 14 15 KSPSetUp_PGMRES - Sets up the workspace needed by pgmres. 16 17 This is called once, usually automatically by KSPSolve() or KSPSetUp(), 18 but can be called directly by KSPSetUp(). 19 20 */ 21 #undef __FUNCT__ 22 #define __FUNCT__ "KSPSetUp_PGMRES" 23 static PetscErrorCode KSPSetUp_PGMRES(KSP ksp) 24 { 25 PetscErrorCode ierr; 26 27 PetscFunctionBegin; 28 ierr = KSPSetUp_GMRES(ksp);CHKERRQ(ierr); 29 PetscFunctionReturn(0); 30 } 31 32 /* 33 34 KSPPGMRESCycle - Run pgmres, possibly with restart. Return residual 35 history if requested. 36 37 input parameters: 38 . pgmres - structure containing parameters and work areas 39 40 output parameters: 41 . itcount - number of iterations used. If null, ignored. 42 . converged - 0 if not converged 43 44 Notes: 45 On entry, the value in vector VEC_VV(0) should be 46 the initial residual. 47 48 49 */ 50 #undef __FUNCT__ 51 #define __FUNCT__ "KSPPGMRESCycle" 52 static PetscErrorCode KSPPGMRESCycle(PetscInt *itcount,KSP ksp) 53 { 54 KSP_PGMRES *pgmres = (KSP_PGMRES*)(ksp->data); 55 PetscReal res_norm,res,newnorm; 56 PetscErrorCode ierr; 57 PetscInt it = 0,j,k; 58 PetscBool hapend = PETSC_FALSE; 59 60 PetscFunctionBegin; 61 if (itcount) *itcount = 0; 62 ierr = VecNormalize(VEC_VV(0),&res_norm);CHKERRQ(ierr); 63 KSPCheckNorm(ksp,res_norm); 64 res = res_norm; 65 *RS(0) = res_norm; 66 67 /* check for the convergence */ 68 ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); 69 ksp->rnorm = res; 70 ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); 71 pgmres->it = it-2; 72 ierr = KSPLogResidualHistory(ksp,res);CHKERRQ(ierr); 73 ierr = KSPMonitor(ksp,ksp->its,res);CHKERRQ(ierr); 74 if (!res) { 75 ksp->reason = KSP_CONVERGED_ATOL; 76 ierr = PetscInfo(ksp,"Converged due to zero residual norm on entry\n");CHKERRQ(ierr); 77 PetscFunctionReturn(0); 78 } 79 80 ierr = (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); 81 for (; !ksp->reason; it++) { 82 Vec Zcur,Znext; 83 if (pgmres->vv_allocated <= it + VEC_OFFSET + 1) { 84 ierr = KSPGMRESGetNewVectors(ksp,it+1);CHKERRQ(ierr); 85 } 86 /* VEC_VV(it-1) is orthogonal, it will be normalized once the VecNorm arrives. */ 87 Zcur = VEC_VV(it); /* Zcur is not yet orthogonal, but the VecMDot to orthogonalize it has been started. */ 88 Znext = VEC_VV(it+1); /* This iteration will compute Znext, update with a deferred correction once we know how 89 * Zcur relates to the previous vectors, and start the reduction to orthogonalize it. */ 90 91 if (it < pgmres->max_k+1 && ksp->its+1 < PetscMax(2,ksp->max_it)) { /* We don't know whether what we have computed is enough, so apply the matrix. */ 92 ierr = KSP_PCApplyBAorAB(ksp,Zcur,Znext,VEC_TEMP_MATOP);CHKERRQ(ierr); 93 } 94 95 if (it > 1) { /* Complete the pending reduction */ 96 ierr = VecNormEnd(VEC_VV(it-1),NORM_2,&newnorm);CHKERRQ(ierr); 97 *HH(it-1,it-2) = newnorm; 98 } 99 if (it > 0) { /* Finish the reduction computing the latest column of H */ 100 ierr = VecMDotEnd(Zcur,it,&(VEC_VV(0)),HH(0,it-1));CHKERRQ(ierr); 101 } 102 103 if (it > 1) { 104 /* normalize the base vector from two iterations ago, basis is complete up to here */ 105 ierr = VecScale(VEC_VV(it-1),1./ *HH(it-1,it-2));CHKERRQ(ierr); 106 107 ierr = KSPPGMRESUpdateHessenberg(ksp,it-2,&hapend,&res);CHKERRQ(ierr); 108 pgmres->it = it-2; 109 ksp->its++; 110 ksp->rnorm = res; 111 112 ierr = (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); 113 if (it < pgmres->max_k+1 || ksp->reason || ksp->its == ksp->max_it) { /* Monitor if we are done or still iterating, but not before a restart. */ 114 ierr = KSPLogResidualHistory(ksp,res);CHKERRQ(ierr); 115 ierr = KSPMonitor(ksp,ksp->its,res);CHKERRQ(ierr); 116 } 117 if (ksp->reason) break; 118 /* Catch error in happy breakdown and signal convergence and break from loop */ 119 if (hapend) { 120 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); 121 else { 122 ksp->reason = KSP_DIVERGED_BREAKDOWN; 123 break; 124 } 125 } 126 127 if (!(it < pgmres->max_k+1 && ksp->its < ksp->max_it)) break; 128 129 /* The it-2 column of H was not scaled when we computed Zcur, apply correction */ 130 ierr = VecScale(Zcur,1./ *HH(it-1,it-2));CHKERRQ(ierr); 131 /* And Znext computed in this iteration was computed using the under-scaled Zcur */ 132 ierr = VecScale(Znext,1./ *HH(it-1,it-2));CHKERRQ(ierr); 133 134 /* In the previous iteration, we projected an unnormalized Zcur against the Krylov basis, so we need to fix the column of H resulting from that projection. */ 135 for (k=0; k<it; k++) *HH(k,it-1) /= *HH(it-1,it-2); 136 /* When Zcur was projected against the Krylov basis, VV(it-1) was still not normalized, so fix that too. This 137 * column is complete except for HH(it,it-1) which we won't know until the next iteration. */ 138 *HH(it-1,it-1) /= *HH(it-1,it-2); 139 } 140 141 if (it > 0) { 142 PetscScalar *work; 143 if (!pgmres->orthogwork) {ierr = PetscMalloc1(pgmres->max_k + 2,&pgmres->orthogwork);CHKERRQ(ierr);} 144 work = pgmres->orthogwork; 145 /* Apply correction computed by the VecMDot in the last iteration to Znext. The original form is 146 * 147 * Znext -= sum_{j=0}^{i-1} Z[j+1] * H[j,i-1] 148 * 149 * where 150 * 151 * Z[j] = sum_{k=0}^j V[k] * H[k,j-1] 152 * 153 * substituting 154 * 155 * Znext -= sum_{j=0}^{i-1} sum_{k=0}^{j+1} V[k] * H[k,j] * H[j,i-1] 156 * 157 * rearranging the iteration space from row-column to column-row 158 * 159 * Znext -= sum_{k=0}^i sum_{j=k-1}^{i-1} V[k] * H[k,j] * H[j,i-1] 160 * 161 * Note that column it-1 of HH is correct. For all previous columns, we must look at HES because HH has already 162 * been transformed to upper triangular form. 163 */ 164 for (k=0; k<it+1; k++) { 165 work[k] = 0; 166 for (j=PetscMax(0,k-1); j<it-1; j++) work[k] -= *HES(k,j) * *HH(j,it-1); 167 } 168 ierr = VecMAXPY(Znext,it+1,work,&VEC_VV(0));CHKERRQ(ierr); 169 ierr = VecAXPY(Znext,-*HH(it-1,it-1),Zcur);CHKERRQ(ierr); 170 171 /* Orthogonalize Zcur against existing basis vectors. */ 172 for (k=0; k<it; k++) work[k] = -*HH(k,it-1); 173 ierr = VecMAXPY(Zcur,it,work,&VEC_VV(0));CHKERRQ(ierr); 174 /* Zcur is now orthogonal, and will be referred to as VEC_VV(it) again, though it is still not normalized. */ 175 /* Begin computing the norm of the new vector, will be normalized after the MatMult in the next iteration. */ 176 ierr = VecNormBegin(VEC_VV(it),NORM_2,&newnorm);CHKERRQ(ierr); 177 } 178 179 /* Compute column of H (to the diagonal, but not the subdiagonal) to be able to orthogonalize the newest vector. */ 180 ierr = VecMDotBegin(Znext,it+1,&VEC_VV(0),HH(0,it));CHKERRQ(ierr); 181 182 /* Start an asynchronous split-mode reduction, the result of the MDot and Norm will be collected on the next iteration. */ 183 ierr = PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)Znext));CHKERRQ(ierr); 184 } 185 186 if (itcount) *itcount = it-1; /* Number of iterations actually completed. */ 187 188 /* 189 Down here we have to solve for the "best" coefficients of the Krylov 190 columns, add the solution values together, and possibly unwind the 191 preconditioning from the solution 192 */ 193 /* Form the solution (or the solution so far) */ 194 ierr = KSPPGMRESBuildSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,it-2);CHKERRQ(ierr); 195 PetscFunctionReturn(0); 196 } 197 198 /* 199 KSPSolve_PGMRES - This routine applies the PGMRES method. 200 201 202 Input Parameter: 203 . ksp - the Krylov space object that was set to use pgmres 204 205 Output Parameter: 206 . outits - number of iterations used 207 208 */ 209 #undef __FUNCT__ 210 #define __FUNCT__ "KSPSolve_PGMRES" 211 static PetscErrorCode KSPSolve_PGMRES(KSP ksp) 212 { 213 PetscErrorCode ierr; 214 PetscInt its,itcount; 215 KSP_PGMRES *pgmres = (KSP_PGMRES*)ksp->data; 216 PetscBool guess_zero = ksp->guess_zero; 217 218 PetscFunctionBegin; 219 if (ksp->calc_sings && !pgmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called"); 220 ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); 221 ksp->its = 0; 222 ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); 223 224 itcount = 0; 225 ksp->reason = KSP_CONVERGED_ITERATING; 226 while (!ksp->reason) { 227 ierr = KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);CHKERRQ(ierr); 228 ierr = KSPPGMRESCycle(&its,ksp);CHKERRQ(ierr); 229 itcount += its; 230 if (itcount >= ksp->max_it) { 231 if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS; 232 break; 233 } 234 ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */ 235 } 236 ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */ 237 PetscFunctionReturn(0); 238 } 239 240 #undef __FUNCT__ 241 #define __FUNCT__ "KSPDestroy_PGMRES" 242 static PetscErrorCode KSPDestroy_PGMRES(KSP ksp) 243 { 244 PetscErrorCode ierr; 245 246 PetscFunctionBegin; 247 ierr = KSPDestroy_GMRES(ksp);CHKERRQ(ierr); 248 PetscFunctionReturn(0); 249 } 250 251 /* 252 KSPPGMRESBuildSoln - create the solution from the starting vector and the 253 current iterates. 254 255 Input parameters: 256 nrs - work area of size it + 1. 257 vguess - index of initial guess 258 vdest - index of result. Note that vguess may == vdest (replace 259 guess with the solution). 260 it - HH upper triangular part is a block of size (it+1) x (it+1) 261 262 This is an internal routine that knows about the PGMRES internals. 263 */ 264 #undef __FUNCT__ 265 #define __FUNCT__ "KSPPGMRESBuildSoln" 266 static PetscErrorCode KSPPGMRESBuildSoln(PetscScalar *nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it) 267 { 268 PetscScalar tt; 269 PetscErrorCode ierr; 270 PetscInt k,j; 271 KSP_PGMRES *pgmres = (KSP_PGMRES*)(ksp->data); 272 273 PetscFunctionBegin; 274 /* Solve for solution vector that minimizes the residual */ 275 276 if (it < 0) { /* no pgmres steps have been performed */ 277 ierr = VecCopy(vguess,vdest);CHKERRQ(ierr); /* VecCopy() is smart, exits immediately if vguess == vdest */ 278 PetscFunctionReturn(0); 279 } 280 281 /* solve the upper triangular system - RS is the right side and HH is 282 the upper triangular matrix - put soln in nrs */ 283 if (*HH(it,it) != 0.0) nrs[it] = *RS(it) / *HH(it,it); 284 else nrs[it] = 0.0; 285 286 for (k=it-1; k>=0; k--) { 287 tt = *RS(k); 288 for (j=k+1; j<=it; j++) tt -= *HH(k,j) * nrs[j]; 289 nrs[k] = tt / *HH(k,k); 290 } 291 292 /* Accumulate the correction to the solution of the preconditioned problem in TEMP */ 293 ierr = VecZeroEntries(VEC_TEMP);CHKERRQ(ierr); 294 ierr = VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));CHKERRQ(ierr); 295 ierr = KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);CHKERRQ(ierr); 296 /* add solution to previous solution */ 297 if (vdest == vguess) { 298 ierr = VecAXPY(vdest,1.0,VEC_TEMP);CHKERRQ(ierr); 299 } else { 300 ierr = VecWAXPY(vdest,1.0,VEC_TEMP,vguess);CHKERRQ(ierr); 301 } 302 PetscFunctionReturn(0); 303 } 304 305 /* 306 307 KSPPGMRESUpdateHessenberg - Do the scalar work for the orthogonalization. 308 Return new residual. 309 310 input parameters: 311 312 . ksp - Krylov space object 313 . it - plane rotations are applied to the (it+1)th column of the 314 modified hessenberg (i.e. HH(:,it)) 315 . hapend - PETSC_FALSE not happy breakdown ending. 316 317 output parameters: 318 . res - the new residual 319 320 */ 321 #undef __FUNCT__ 322 #define __FUNCT__ "KSPPGMRESUpdateHessenberg" 323 /* 324 . it - column of the Hessenberg that is complete, PGMRES is actually computing two columns ahead of this 325 */ 326 static PetscErrorCode KSPPGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool *hapend,PetscReal *res) 327 { 328 PetscScalar *hh,*cc,*ss,*rs; 329 PetscInt j; 330 PetscReal hapbnd; 331 KSP_PGMRES *pgmres = (KSP_PGMRES*)(ksp->data); 332 PetscErrorCode ierr; 333 334 PetscFunctionBegin; 335 hh = HH(0,it); /* pointer to beginning of column to update */ 336 cc = CC(0); /* beginning of cosine rotations */ 337 ss = SS(0); /* beginning of sine rotations */ 338 rs = RS(0); /* right hand side of least squares system */ 339 340 /* The Hessenberg matrix is now correct through column it, save that form for possible spectral analysis */ 341 for (j=0; j<=it+1; j++) *HES(j,it) = hh[j]; 342 343 /* check for the happy breakdown */ 344 hapbnd = PetscMin(PetscAbsScalar(hh[it+1] / rs[it]),pgmres->haptol); 345 if (PetscAbsScalar(hh[it+1]) < hapbnd) { 346 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); 347 *hapend = PETSC_TRUE; 348 } 349 350 /* Apply all the previously computed plane rotations to the new column 351 of the Hessenberg matrix */ 352 /* Note: this uses the rotation [conj(c) s ; -s c], c= cos(theta), s= sin(theta), 353 and some refs have [c s ; -conj(s) c] (don't be confused!) */ 354 355 for (j=0; j<it; j++) { 356 PetscScalar hhj = hh[j]; 357 hh[j] = PetscConj(cc[j])*hhj + ss[j]*hh[j+1]; 358 hh[j+1] = -ss[j] *hhj + cc[j]*hh[j+1]; 359 } 360 361 /* 362 compute the new plane rotation, and apply it to: 363 1) the right-hand-side of the Hessenberg system (RS) 364 note: it affects RS(it) and RS(it+1) 365 2) the new column of the Hessenberg matrix 366 note: it affects HH(it,it) which is currently pointed to 367 by hh and HH(it+1, it) (*(hh+1)) 368 thus obtaining the updated value of the residual... 369 */ 370 371 /* compute new plane rotation */ 372 373 if (!*hapend) { 374 PetscReal delta = PetscSqrtReal(PetscSqr(PetscAbsScalar(hh[it])) + PetscSqr(PetscAbsScalar(hh[it+1]))); 375 if (delta == 0.0) { 376 ksp->reason = KSP_DIVERGED_NULL; 377 PetscFunctionReturn(0); 378 } 379 380 cc[it] = hh[it] / delta; /* new cosine value */ 381 ss[it] = hh[it+1] / delta; /* new sine value */ 382 383 hh[it] = PetscConj(cc[it])*hh[it] + ss[it]*hh[it+1]; 384 rs[it+1] = -ss[it]*rs[it]; 385 rs[it] = PetscConj(cc[it])*rs[it]; 386 *res = PetscAbsScalar(rs[it+1]); 387 } else { /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply 388 another rotation matrix (so RH doesn't change). The new residual is 389 always the new sine term times the residual from last time (RS(it)), 390 but now the new sine rotation would be zero...so the residual should 391 be zero...so we will multiply "zero" by the last residual. This might 392 not be exactly what we want to do here -could just return "zero". */ 393 394 *res = 0.0; 395 } 396 PetscFunctionReturn(0); 397 } 398 399 /* 400 KSPBuildSolution_PGMRES 401 402 Input Parameter: 403 . ksp - the Krylov space object 404 . ptr- 405 406 Output Parameter: 407 . result - the solution 408 409 Note: this calls KSPPGMRESBuildSoln - the same function that KSPPGMRESCycle 410 calls directly. 411 412 */ 413 #undef __FUNCT__ 414 #define __FUNCT__ "KSPBuildSolution_PGMRES" 415 PetscErrorCode KSPBuildSolution_PGMRES(KSP ksp,Vec ptr,Vec *result) 416 { 417 KSP_PGMRES *pgmres = (KSP_PGMRES*)ksp->data; 418 PetscErrorCode ierr; 419 420 PetscFunctionBegin; 421 if (!ptr) { 422 if (!pgmres->sol_temp) { 423 ierr = VecDuplicate(ksp->vec_sol,&pgmres->sol_temp);CHKERRQ(ierr); 424 ierr = PetscLogObjectParent((PetscObject)ksp,(PetscObject)pgmres->sol_temp);CHKERRQ(ierr); 425 } 426 ptr = pgmres->sol_temp; 427 } 428 if (!pgmres->nrs) { 429 /* allocate the work area */ 430 ierr = PetscMalloc1(pgmres->max_k,&pgmres->nrs);CHKERRQ(ierr); 431 ierr = PetscLogObjectMemory((PetscObject)ksp,pgmres->max_k*sizeof(PetscScalar));CHKERRQ(ierr); 432 } 433 434 ierr = KSPPGMRESBuildSoln(pgmres->nrs,ksp->vec_sol,ptr,ksp,pgmres->it);CHKERRQ(ierr); 435 if (result) *result = ptr; 436 PetscFunctionReturn(0); 437 } 438 439 #undef __FUNCT__ 440 #define __FUNCT__ "KSPSetFromOptions_PGMRES" 441 PetscErrorCode KSPSetFromOptions_PGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp) 442 { 443 PetscErrorCode ierr; 444 445 PetscFunctionBegin; 446 ierr = KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);CHKERRQ(ierr); 447 ierr = PetscOptionsHead(PetscOptionsObject,"KSP pipelined GMRES Options");CHKERRQ(ierr); 448 ierr = PetscOptionsTail();CHKERRQ(ierr); 449 PetscFunctionReturn(0); 450 } 451 452 #undef __FUNCT__ 453 #define __FUNCT__ "KSPReset_PGMRES" 454 PetscErrorCode KSPReset_PGMRES(KSP ksp) 455 { 456 PetscErrorCode ierr; 457 458 PetscFunctionBegin; 459 ierr = KSPReset_GMRES(ksp);CHKERRQ(ierr); 460 PetscFunctionReturn(0); 461 } 462 463 /*MC 464 KSPPGMRES - Implements the Pipelined Generalized Minimal Residual method. 465 466 Options Database Keys: 467 + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against 468 . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence) 469 . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of 470 vectors are allocated as needed) 471 . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default) 472 . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower) 473 . -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the 474 stability of the classical Gram-Schmidt orthogonalization. 475 - -ksp_gmres_krylov_monitor - plot the Krylov space generated 476 477 Level: beginner 478 479 Notes: 480 MPI configuration may be necessary for reductions to make asynchronous progress, which is important for performance of pipelined methods. 481 See the FAQ on the PETSc website for details. 482 483 Reference: 484 Ghysels, Ashby, Meerbergen, Vanroose, Hiding global communication latencies in the GMRES algorithm on massively parallel machines, 2012. 485 486 Developer Notes: This object is subclassed off of KSPGMRES 487 488 .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES, KSPPIPECG, KSPPIPECR, 489 KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 490 KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(), 491 KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov() 492 M*/ 493 494 #undef __FUNCT__ 495 #define __FUNCT__ "KSPCreate_PGMRES" 496 PETSC_EXTERN PetscErrorCode KSPCreate_PGMRES(KSP ksp) 497 { 498 KSP_PGMRES *pgmres; 499 PetscErrorCode ierr; 500 501 PetscFunctionBegin; 502 ierr = PetscNewLog(ksp,&pgmres);CHKERRQ(ierr); 503 504 ksp->data = (void*)pgmres; 505 ksp->ops->buildsolution = KSPBuildSolution_PGMRES; 506 ksp->ops->setup = KSPSetUp_PGMRES; 507 ksp->ops->solve = KSPSolve_PGMRES; 508 ksp->ops->reset = KSPReset_PGMRES; 509 ksp->ops->destroy = KSPDestroy_PGMRES; 510 ksp->ops->view = KSPView_GMRES; 511 ksp->ops->setfromoptions = KSPSetFromOptions_PGMRES; 512 ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES; 513 ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES; 514 515 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);CHKERRQ(ierr); 516 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);CHKERRQ(ierr); 517 518 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);CHKERRQ(ierr); 519 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);CHKERRQ(ierr); 520 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);CHKERRQ(ierr); 521 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);CHKERRQ(ierr); 522 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);CHKERRQ(ierr); 523 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);CHKERRQ(ierr); 524 ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);CHKERRQ(ierr); 525 526 pgmres->nextra_vecs = 1; 527 pgmres->haptol = 1.0e-30; 528 pgmres->q_preallocate = 0; 529 pgmres->delta_allocate = PGMRES_DELTA_DIRECTIONS; 530 pgmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization; 531 pgmres->nrs = 0; 532 pgmres->sol_temp = 0; 533 pgmres->max_k = PGMRES_DEFAULT_MAXK; 534 pgmres->Rsvd = 0; 535 pgmres->orthogwork = 0; 536 pgmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER; 537 PetscFunctionReturn(0); 538 } 539