static PetscErrorCode KSPSolve_AGMRES(KSP ksp) { PetscErrorCode ierr; PetscInt its; KSP_AGMRES *agmres = (KSP_AGMRES*)ksp->data; PetscBool guess_zero = ksp->guess_zero; PetscReal res_old, res; PetscInt test; PetscFunctionBegin; ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); ksp->its = 0; ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); ksp->reason = KSP_CONVERGED_ITERATING; if (!agmres->HasShifts) { /* Compute Shifts for the Newton basis */ ierr = KSPComputeShifts_DGMRES(ksp);CHKERRQ(ierr); } /* NOTE: At this step, the initial guess is not equal to zero since one cycle of the classical GMRES is performed to compute the shifts */ ierr = (*ksp->converged)(ksp,0,ksp->rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); while (!ksp->reason) { ierr = KSPInitialResidual(ksp,ksp->vec_sol,VEC_TMP,VEC_TMP_MATOP,VEC_V(0),ksp->vec_rhs);CHKERRQ(ierr); if ((ksp->pc_side == PC_LEFT) && agmres->r && agmres->DeflPrecond) { ierr = KSPDGMRESApplyDeflation_DGMRES(ksp, VEC_V(0), VEC_TMP);CHKERRQ(ierr); ierr = VecCopy(VEC_TMP, VEC_V(0));CHKERRQ(ierr); agmres->matvecs += 1; } ierr = VecNormalize(VEC_V(0),&(ksp->rnorm));CHKERRQ(ierr); KSPCheckNorm(ksp,ksp->rnorm); res_old = ksp->rnorm; /* Record the residual norm to test if deflation is needed */ ksp->ops->buildsolution = KSPBuildSolution_AGMRES; ierr = KSPAGMRESCycle(&its,ksp);CHKERRQ(ierr); if (ksp->its >= ksp->max_it) { if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS; break; } /* compute the eigenvectors to augment the subspace : use an adaptive strategy */ res = ksp->rnorm; if (!ksp->reason && agmres->neig > 0) { test = agmres->max_k * PetscLogReal(ksp->rtol/res) / PetscLogReal(res/res_old); /* estimate the remaining number of steps */ if ((test > agmres->smv*(ksp->max_it-ksp->its)) || agmres->force) { if (!agmres->force && ((test > agmres->bgv*(ksp->max_it-ksp->its)) && ((agmres->r + 1) < agmres->max_neig))) { agmres->neig += 1; /* Augment the number of eigenvalues to deflate if the convergence is too slow */ } ierr = KSPDGMRESComputeDeflationData_DGMRES(ksp,&agmres->neig);CHKERRQ(ierr); } } ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */ } ksp->guess_zero = guess_zero; /* restore if user has provided nonzero initial guess */ PetscFunctionReturn(0); }
PetscErrorCode KSPComputeShifts_DGMRES(KSP ksp) { PetscErrorCode ierr; KSP_AGMRES *agmres = (KSP_AGMRES*)(ksp->data); PetscInt max_k = agmres->max_k; /* size of the (non augmented) Krylov subspace */ PetscInt Neig = 0; PetscInt max_it = ksp->max_it; /* Perform one cycle of dgmres to find the eigenvalues and compute the first approximations of the eigenvectors */ PetscFunctionBegin; ierr = PetscLogEventBegin(KSP_AGMRESComputeShifts, ksp, 0,0,0);CHKERRQ(ierr); /* Send the size of the augmented basis to DGMRES */ ksp->max_it = max_k; /* set this to have DGMRES performing only one cycle */ ksp->ops->buildsolution = KSPBuildSolution_DGMRES; ierr = KSPSolve_DGMRES(ksp); ksp->guess_zero = PETSC_FALSE; if (ksp->reason == KSP_CONVERGED_RTOL) { ierr = PetscLogEventEnd(KSP_AGMRESComputeShifts, ksp, 0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } else ksp->reason = KSP_CONVERGED_ITERATING; if ((agmres->r == 0) && (agmres->neig > 0)) { /* Compute the eigenvalues for the shifts and the eigenvectors (to augment the Newton basis) */ agmres->HasSchur = PETSC_FALSE; ierr = KSPDGMRESComputeDeflationData_DGMRES (ksp, &Neig);CHKERRQ (ierr); Neig = max_k; } else { /* From DGMRES, compute only the eigenvalues needed as Shifts for the Newton Basis */ ierr = KSPDGMRESComputeSchurForm_DGMRES(ksp, &Neig);CHKERRQ(ierr); } /* It may happen that the Ritz values from one cycle of GMRES are not accurate enough to provide a good stability. In this case, another cycle of GMRES is performed. The two sets of values thus generated are sorted and the most accurate are kept as shifts */ PetscBool flg; ierr = PetscOptionsHasName(NULL, "-ksp_agmres_ImproveShifts", &flg);CHKERRQ(ierr); if (!flg) { ierr = KSPAGMRESLejaOrdering(agmres->wr, agmres->wi, agmres->Rshift, agmres->Ishift, max_k);CHKERRQ(ierr); } else { /* Perform another cycle of DGMRES to find another set of eigenvalues */ PetscInt i; PetscScalar *wr, *wi,*Rshift, *Ishift; ierr = PetscMalloc4(2*max_k, &wr, 2*max_k, &wi, 2*max_k, &Rshift, 2*max_k, &Ishift);CHKERRQ(ierr); for (i = 0; i < max_k; i++) { wr[i] = agmres->wr[i]; wi[i] = agmres->wi[i]; } ierr = KSPSolve_DGMRES(ksp); ksp->guess_zero = PETSC_FALSE; if (ksp->reason == KSP_CONVERGED_RTOL) PetscFunctionReturn(0); else ksp->reason = KSP_CONVERGED_ITERATING; if (agmres->neig > 0) { /* Compute the eigenvalues for the shifts) and the eigenvectors (to augment the Newton basis */ agmres->HasSchur = PETSC_FALSE; ierr = KSPDGMRESComputeDeflationData_DGMRES(ksp, &Neig);CHKERRQ(ierr); Neig = max_k; } else { /* From DGMRES, compute only the eigenvalues needed as Shifts for the Newton Basis */ ierr = KSPDGMRESComputeSchurForm_DGMRES(ksp, &Neig);CHKERRQ(ierr); } for (i = 0; i < max_k; i++) { wr[max_k+i] = agmres->wr[i]; wi[max_k+i] = agmres->wi[i]; } ierr = KSPAGMRESLejaOrdering(wr, wi, Rshift, Ishift, 2*max_k);CHKERRQ(ierr); for (i = 0; i< max_k; i++) { agmres->Rshift[i] = Rshift[i]; agmres->Ishift[i] = Ishift[i]; } ierr = PetscFree(Rshift);CHKERRQ(ierr); ierr = PetscFree(wr);CHKERRQ(ierr); ierr = PetscFree(Ishift);CHKERRQ(ierr); ierr = PetscFree(wi);CHKERRQ(ierr); } agmres->HasShifts = PETSC_TRUE; ksp->max_it = max_it; ierr = PetscLogEventEnd(KSP_AGMRESComputeShifts, ksp, 0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); }