Actual source code: dgmres.c

  1: /*
  2:     Implements deflated GMRES.
  3: */

  5: #include <../src/ksp/ksp/impls/gmres/dgmres/dgmresimpl.h>

  7: PetscLogEvent KSP_DGMRESComputeDeflationData, KSP_DGMRESApplyDeflation;

  9: static PetscErrorCode KSPDGMRESGetNewVectors(KSP, PetscInt);
 10: static PetscErrorCode KSPDGMRESUpdateHessenberg(KSP, PetscInt, PetscBool, PetscReal *);
 11: static PetscErrorCode KSPDGMRESBuildSoln(PetscScalar *, Vec, Vec, KSP, PetscInt);

 13: static PetscErrorCode KSPDGMRESSetEigen(KSP ksp, PetscInt nb_eig)
 14: {
 15:   PetscFunctionBegin;
 16:   PetscTryMethod((ksp), "KSPDGMRESSetEigen_C", (KSP, PetscInt), (ksp, nb_eig));
 17:   PetscFunctionReturn(PETSC_SUCCESS);
 18: }
 19: static PetscErrorCode KSPDGMRESSetMaxEigen(KSP ksp, PetscInt max_neig)
 20: {
 21:   PetscFunctionBegin;
 22:   PetscTryMethod((ksp), "KSPDGMRESSetMaxEigen_C", (KSP, PetscInt), (ksp, max_neig));
 23:   PetscFunctionReturn(PETSC_SUCCESS);
 24: }
 25: static PetscErrorCode KSPDGMRESComputeSchurForm(KSP ksp, PetscInt *neig)
 26: {
 27:   PetscFunctionBegin;
 28:   PetscUseMethod((ksp), "KSPDGMRESComputeSchurForm_C", (KSP, PetscInt *), (ksp, neig));
 29:   PetscFunctionReturn(PETSC_SUCCESS);
 30: }
 31: PetscErrorCode KSPDGMRESComputeDeflationData(KSP ksp, PetscInt *curneigh)
 32: {
 33:   PetscFunctionBegin;
 34:   PetscUseMethod((ksp), "KSPDGMRESComputeDeflationData_C", (KSP, PetscInt *), (ksp, curneigh));
 35:   PetscFunctionReturn(PETSC_SUCCESS);
 36: }
 37: static PetscErrorCode KSPDGMRESApplyDeflation(KSP ksp, Vec x, Vec y)
 38: {
 39:   PetscFunctionBegin;
 40:   PetscUseMethod((ksp), "KSPDGMRESApplyDeflation_C", (KSP, Vec, Vec), (ksp, x, y));
 41:   PetscFunctionReturn(PETSC_SUCCESS);
 42: }

 44: static PetscErrorCode KSPDGMRESImproveEig(KSP ksp, PetscInt neig)
 45: {
 46:   PetscFunctionBegin;
 47:   PetscUseMethod((ksp), "KSPDGMRESImproveEig_C", (KSP, PetscInt), (ksp, neig));
 48:   PetscFunctionReturn(PETSC_SUCCESS);
 49: }

 51: PetscErrorCode KSPSetUp_DGMRES(KSP ksp)
 52: {
 53:   KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
 54:   PetscInt    neig   = dgmres->neig + EIG_OFFSET;
 55:   PetscInt    max_k  = dgmres->max_k + 1;

 57:   PetscFunctionBegin;
 58:   PetscCall(KSPSetUp_GMRES(ksp));
 59:   if (!dgmres->neig) PetscFunctionReturn(PETSC_SUCCESS);

 61:   /* Allocate workspace for the Schur vectors*/
 62:   PetscCall(PetscMalloc1(neig * max_k, &SR));
 63:   dgmres->wr    = NULL;
 64:   dgmres->wi    = NULL;
 65:   dgmres->perm  = NULL;
 66:   dgmres->modul = NULL;
 67:   dgmres->Q     = NULL;
 68:   dgmres->Z     = NULL;

 70:   UU   = NULL;
 71:   XX   = NULL;
 72:   MX   = NULL;
 73:   AUU  = NULL;
 74:   XMX  = NULL;
 75:   XMU  = NULL;
 76:   UMX  = NULL;
 77:   AUAU = NULL;
 78:   TT   = NULL;
 79:   TTF  = NULL;
 80:   INVP = NULL;
 81:   X1   = NULL;
 82:   X2   = NULL;
 83:   MU   = NULL;
 84:   PetscFunctionReturn(PETSC_SUCCESS);
 85: }

 87: /*
 88:  Run GMRES, possibly with restart.  Return residual history if requested.
 89:  input parameters:

 91:  .       gmres  - structure containing parameters and work areas

 93:  output parameters:
 94:  .        nres    - residuals (from preconditioned system) at each step.
 95:  If restarting, consider passing nres+it.  If null,
 96:  ignored
 97:  .        itcount - number of iterations used.  nres[0] to nres[itcount]
 98:  are defined.  If null, ignored.

100:  Notes:
101:  On entry, the value in vector VEC_VV(0) should be the initial residual
102:  (this allows shortcuts where the initial preconditioned residual is 0).
103:  */
104: static PetscErrorCode KSPDGMRESCycle(PetscInt *itcount, KSP ksp)
105: {
106:   KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
107:   PetscReal   res_norm, res, hapbnd, tt;
108:   PetscInt    it     = 0;
109:   PetscInt    max_k  = dgmres->max_k;
110:   PetscBool   hapend = PETSC_FALSE;
111:   PetscReal   res_old;
112:   PetscInt    test = 0;

114:   PetscFunctionBegin;
115:   PetscCall(VecNormalize(VEC_VV(0), &res_norm));
116:   KSPCheckNorm(ksp, res_norm);
117:   res     = res_norm;
118:   *GRS(0) = res_norm;

120:   /* check for the convergence */
121:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
122:   if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res;
123:   else ksp->rnorm = 0.0;
124:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
125:   dgmres->it = (it - 1);
126:   PetscCall(KSPLogResidualHistory(ksp, ksp->rnorm));
127:   PetscCall(KSPMonitor(ksp, ksp->its, ksp->rnorm));
128:   if (!res) {
129:     if (itcount) *itcount = 0;
130:     ksp->reason = KSP_CONVERGED_ATOL;
131:     PetscCall(PetscInfo(ksp, "Converged due to zero residual norm on entry\n"));
132:     PetscFunctionReturn(PETSC_SUCCESS);
133:   }
134:   /* record the residual norm to test if deflation is needed */
135:   res_old = res;

137:   PetscCall((*ksp->converged)(ksp, ksp->its, ksp->rnorm, &ksp->reason, ksp->cnvP));
138:   while (!ksp->reason && it < max_k && ksp->its < ksp->max_it) {
139:     if (it) {
140:       PetscCall(KSPLogResidualHistory(ksp, ksp->rnorm));
141:       PetscCall(KSPMonitor(ksp, ksp->its, ksp->rnorm));
142:     }
143:     dgmres->it = (it - 1);
144:     if (dgmres->vv_allocated <= it + VEC_OFFSET + 1) PetscCall(KSPDGMRESGetNewVectors(ksp, it + 1));
145:     if (dgmres->r > 0) {
146:       if (ksp->pc_side == PC_LEFT) {
147:         /* Apply the first preconditioner */
148:         PetscCall(KSP_PCApplyBAorAB(ksp, VEC_VV(it), VEC_TEMP, VEC_TEMP_MATOP));
149:         /* Then apply Deflation as a preconditioner */
150:         PetscCall(KSPDGMRESApplyDeflation(ksp, VEC_TEMP, VEC_VV(1 + it)));
151:       } else if (ksp->pc_side == PC_RIGHT) {
152:         PetscCall(KSPDGMRESApplyDeflation(ksp, VEC_VV(it), VEC_TEMP));
153:         PetscCall(KSP_PCApplyBAorAB(ksp, VEC_TEMP, VEC_VV(1 + it), VEC_TEMP_MATOP));
154:       }
155:     } else {
156:       PetscCall(KSP_PCApplyBAorAB(ksp, VEC_VV(it), VEC_VV(1 + it), VEC_TEMP_MATOP));
157:     }
158:     dgmres->matvecs += 1;
159:     /* update Hessenberg matrix and do Gram-Schmidt */
160:     PetscCall((*dgmres->orthog)(ksp, it));

162:     /* vv(i+1) . vv(i+1) */
163:     PetscCall(VecNormalize(VEC_VV(it + 1), &tt));
164:     /* save the magnitude */
165:     *HH(it + 1, it)  = tt;
166:     *HES(it + 1, it) = tt;

168:     /* check for the happy breakdown */
169:     hapbnd = PetscAbsScalar(tt / *GRS(it));
170:     if (hapbnd > dgmres->haptol) hapbnd = dgmres->haptol;
171:     if (tt < hapbnd) {
172:       PetscCall(PetscInfo(ksp, "Detected happy breakdown, current hapbnd = %g tt = %g\n", (double)hapbnd, (double)tt));
173:       hapend = PETSC_TRUE;
174:     }
175:     PetscCall(KSPDGMRESUpdateHessenberg(ksp, it, hapend, &res));

177:     it++;
178:     dgmres->it = (it - 1); /* For converged */
179:     ksp->its++;
180:     if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res;
181:     else ksp->rnorm = 0.0;
182:     if (ksp->reason) break;

184:     PetscCall((*ksp->converged)(ksp, ksp->its, ksp->rnorm, &ksp->reason, ksp->cnvP));

186:     /* Catch error in happy breakdown and signal convergence and break from loop */
187:     if (hapend) {
188:       if (!ksp->reason) {
189:         PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "Reached happy break down, but convergence was not indicated. Residual norm = %g", (double)res);
190:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
191:         break;
192:       }
193:     }
194:   }

196:   if (itcount) *itcount = it;

198:   /*
199:    Down here we have to solve for the "best" coefficients of the Krylov
200:    columns, add the solution values together, and possibly unwind the
201:    preconditioning from the solution
202:    */
203:   /* Form the solution (or the solution so far) */
204:   PetscCall(KSPDGMRESBuildSoln(GRS(0), ksp->vec_sol, ksp->vec_sol, ksp, it - 1));

206:   /* Monitor if we know that we will not return for a restart */
207:   if (it && (ksp->reason || ksp->its >= ksp->max_it)) {
208:     PetscCall(KSPLogResidualHistory(ksp, ksp->rnorm));
209:     PetscCall(KSPMonitor(ksp, ksp->its, ksp->rnorm));
210:   }

212:   /* Compute data for the deflation to be used during the next restart */
213:   if (!ksp->reason && ksp->its < ksp->max_it) {
214:     test = max_k * PetscLogReal(ksp->rtol / res) / PetscLogReal(res / res_old);
215:     /* Compute data for the deflation if the residual rtol will not be reached in the remaining number of steps allowed  */
216:     if ((test > dgmres->smv * (ksp->max_it - ksp->its)) || dgmres->force) PetscCall(KSPDGMRESComputeDeflationData(ksp, NULL));
217:   }
218:   PetscFunctionReturn(PETSC_SUCCESS);
219: }

221: PetscErrorCode KSPSolve_DGMRES(KSP ksp)
222: {
223:   PetscInt    i, its = 0, itcount;
224:   KSP_DGMRES *dgmres     = (KSP_DGMRES *)ksp->data;
225:   PetscBool   guess_zero = ksp->guess_zero;

227:   PetscFunctionBegin;
228:   PetscCheck(!ksp->calc_sings || dgmres->Rsvd, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ORDER, "Must call KSPSetComputeSingularValues() before KSPSetUp() is called");

230:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
231:   ksp->its        = 0;
232:   dgmres->matvecs = 0;
233:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));

235:   itcount = 0;
236:   while (!ksp->reason) {
237:     PetscCall(KSPInitialResidual(ksp, ksp->vec_sol, VEC_TEMP, VEC_TEMP_MATOP, VEC_VV(0), ksp->vec_rhs));
238:     if (ksp->pc_side == PC_LEFT) {
239:       dgmres->matvecs += 1;
240:       if (dgmres->r > 0) {
241:         PetscCall(KSPDGMRESApplyDeflation(ksp, VEC_VV(0), VEC_TEMP));
242:         PetscCall(VecCopy(VEC_TEMP, VEC_VV(0)));
243:       }
244:     }

246:     PetscCall(KSPDGMRESCycle(&its, ksp));
247:     itcount += its;
248:     if (itcount >= ksp->max_it) {
249:       if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
250:       break;
251:     }
252:     ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
253:   }
254:   ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */

256:   for (i = 0; i < dgmres->r; i++) PetscCall(VecViewFromOptions(UU[i], (PetscObject)ksp, "-ksp_dgmres_view_deflation_vecs"));
257:   PetscFunctionReturn(PETSC_SUCCESS);
258: }

260: PetscErrorCode KSPDestroy_DGMRES(KSP ksp)
261: {
262:   KSP_DGMRES *dgmres   = (KSP_DGMRES *)ksp->data;
263:   PetscInt    neig1    = dgmres->neig + EIG_OFFSET;
264:   PetscInt    max_neig = dgmres->max_neig;

266:   PetscFunctionBegin;
267:   if (dgmres->r) {
268:     PetscCall(VecDestroyVecs(max_neig, &UU));
269:     PetscCall(VecDestroyVecs(max_neig, &MU));
270:     if (XX) {
271:       PetscCall(VecDestroyVecs(neig1, &XX));
272:       PetscCall(VecDestroyVecs(neig1, &MX));
273:     }
274:     PetscCall(PetscFree(TT));
275:     PetscCall(PetscFree(TTF));
276:     PetscCall(PetscFree(INVP));
277:     PetscCall(PetscFree(XMX));
278:     PetscCall(PetscFree(UMX));
279:     PetscCall(PetscFree(XMU));
280:     PetscCall(PetscFree(X1));
281:     PetscCall(PetscFree(X2));
282:     PetscCall(PetscFree(dgmres->work));
283:     PetscCall(PetscFree(dgmres->iwork));
284:     PetscCall(PetscFree(dgmres->wr));
285:     PetscCall(PetscFree(dgmres->wi));
286:     PetscCall(PetscFree(dgmres->modul));
287:     PetscCall(PetscFree(dgmres->Q));
288:     PetscCall(PetscFree(ORTH));
289:     PetscCall(PetscFree(AUAU));
290:     PetscCall(PetscFree(AUU));
291:     PetscCall(PetscFree(SR2));
292:   }
293:   PetscCall(PetscFree(SR));
294:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetEigen_C", NULL));
295:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetMaxEigen_C", NULL));
296:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetRatio_C", NULL));
297:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESForce_C", NULL));
298:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESComputeSchurForm_C", NULL));
299:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESComputeDeflationData_C", NULL));
300:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESApplyDeflation_C", NULL));
301:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESImproveEig_C", NULL));
302:   PetscCall(KSPDestroy_GMRES(ksp));
303:   PetscFunctionReturn(PETSC_SUCCESS);
304: }

306: /*
307:  KSPDGMRESBuildSoln - create the solution from the starting vector and the
308:  current iterates.

310:  Input parameters:
311:  nrs - work area of size it + 1.
312:  vs  - index of initial guess
313:  vdest - index of result.  Note that vs may == vdest (replace
314:  guess with the solution).

316:  This is an internal routine that knows about the GMRES internals.
317:  */
318: static PetscErrorCode KSPDGMRESBuildSoln(PetscScalar *nrs, Vec vs, Vec vdest, KSP ksp, PetscInt it)
319: {
320:   PetscScalar tt;
321:   PetscInt    ii, k, j;
322:   KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;

324:   /* Solve for solution vector that minimizes the residual */

326:   PetscFunctionBegin;
327:   /* If it is < 0, no gmres steps have been performed */
328:   if (it < 0) {
329:     PetscCall(VecCopy(vs, vdest)); /* VecCopy() is smart, exists immediately if vguess == vdest */
330:     PetscFunctionReturn(PETSC_SUCCESS);
331:   }
332:   PetscCheck(*HH(it, it) != 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_CONV_FAILED, "Likely your matrix is the zero operator. HH(it,it) is identically zero; it = %" PetscInt_FMT " GRS(it) = %g", it, (double)PetscAbsScalar(*GRS(it)));
333:   if (*HH(it, it) != 0.0) nrs[it] = *GRS(it) / *HH(it, it);
334:   else nrs[it] = 0.0;

336:   for (ii = 1; ii <= it; ii++) {
337:     k  = it - ii;
338:     tt = *GRS(k);
339:     for (j = k + 1; j <= it; j++) tt = tt - *HH(k, j) * nrs[j];
340:     PetscCheck(*HH(k, k) != 0.0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_CONV_FAILED, "Likely your matrix is singular. HH(k,k) is identically zero; it = %" PetscInt_FMT " k = %" PetscInt_FMT, it, k);
341:     nrs[k] = tt / *HH(k, k);
342:   }

344:   /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
345:   PetscCall(VecMAXPBY(VEC_TEMP, it + 1, nrs, 0, &VEC_VV(0)));

347:   /* Apply deflation */
348:   if (ksp->pc_side == PC_RIGHT && dgmres->r > 0) {
349:     PetscCall(KSPDGMRESApplyDeflation(ksp, VEC_TEMP, VEC_TEMP_MATOP));
350:     PetscCall(VecCopy(VEC_TEMP_MATOP, VEC_TEMP));
351:   }
352:   PetscCall(KSPUnwindPreconditioner(ksp, VEC_TEMP, VEC_TEMP_MATOP));

354:   /* add solution to previous solution */
355:   if (vdest != vs) PetscCall(VecCopy(vs, vdest));
356:   PetscCall(VecAXPY(vdest, 1.0, VEC_TEMP));
357:   PetscFunctionReturn(PETSC_SUCCESS);
358: }

360: /*
361:  Do the scalar work for the orthogonalization.  Return new residual norm.
362:  */
363: static PetscErrorCode KSPDGMRESUpdateHessenberg(KSP ksp, PetscInt it, PetscBool hapend, PetscReal *res)
364: {
365:   PetscScalar *hh, *cc, *ss, tt;
366:   PetscInt     j;
367:   KSP_DGMRES  *dgmres = (KSP_DGMRES *)ksp->data;

369:   PetscFunctionBegin;
370:   hh = HH(0, it);
371:   cc = CC(0);
372:   ss = SS(0);

374:   /* Apply all the previously computed plane rotations to the new column
375:    of the Hessenberg matrix */
376:   for (j = 1; j <= it; j++) {
377:     tt  = *hh;
378:     *hh = PetscConj(*cc) * tt + *ss * *(hh + 1);
379:     hh++;
380:     *hh = *cc++ * *hh - (*ss++ * tt);
381:   }

383:   /*
384:    compute the new plane rotation, and apply it to:
385:    1) the right-hand side of the Hessenberg system
386:    2) the new column of the Hessenberg matrix
387:    thus obtaining the updated value of the residual
388:    */
389:   if (!hapend) {
390:     tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh + 1)) * *(hh + 1));
391:     if (tt == 0.0) {
392:       ksp->reason = KSP_DIVERGED_NULL;
393:       PetscFunctionReturn(PETSC_SUCCESS);
394:     }
395:     *cc          = *hh / tt;
396:     *ss          = *(hh + 1) / tt;
397:     *GRS(it + 1) = -(*ss * *GRS(it));
398:     *GRS(it)     = PetscConj(*cc) * *GRS(it);
399:     *hh          = PetscConj(*cc) * *hh + *ss * *(hh + 1);
400:     *res         = PetscAbsScalar(*GRS(it + 1));
401:   } else {
402:     /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
403:      another rotation matrix (so RH doesn't change).  The new residual is
404:      always the new sine term times the residual from last time (GRS(it)),
405:      but now the new sine rotation would be zero...so the residual should
406:      be zero...so we will multiply "zero" by the last residual.  This might
407:      not be exactly what we want to do here -could just return "zero". */
408:     *res = 0.0;
409:   }
410:   PetscFunctionReturn(PETSC_SUCCESS);
411: }

413: /*
414:   Allocates more work vectors, starting from VEC_VV(it).
415: */
416: static PetscErrorCode KSPDGMRESGetNewVectors(KSP ksp, PetscInt it)
417: {
418:   KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
419:   PetscInt    nwork  = dgmres->nwork_alloc, k, nalloc;

421:   PetscFunctionBegin;
422:   nalloc = PetscMin(ksp->max_it, dgmres->delta_allocate);
423:   /* Adjust the number to allocate to make sure that we don't exceed the
424:    number of available slots */
425:   if (it + VEC_OFFSET + nalloc >= dgmres->vecs_allocated) nalloc = dgmres->vecs_allocated - it - VEC_OFFSET;
426:   if (!nalloc) PetscFunctionReturn(PETSC_SUCCESS);

428:   dgmres->vv_allocated += nalloc;

430:   PetscCall(KSPCreateVecs(ksp, nalloc, &dgmres->user_work[nwork], 0, NULL));

432:   dgmres->mwork_alloc[nwork] = nalloc;
433:   for (k = 0; k < nalloc; k++) dgmres->vecs[it + VEC_OFFSET + k] = dgmres->user_work[nwork][k];
434:   dgmres->nwork_alloc++;
435:   PetscFunctionReturn(PETSC_SUCCESS);
436: }

438: PetscErrorCode KSPBuildSolution_DGMRES(KSP ksp, Vec ptr, Vec *result)
439: {
440:   KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;

442:   PetscFunctionBegin;
443:   if (!ptr) {
444:     if (!dgmres->sol_temp) PetscCall(VecDuplicate(ksp->vec_sol, &dgmres->sol_temp));
445:     ptr = dgmres->sol_temp;
446:   }
447:   if (!dgmres->nrs) {
448:     /* allocate the work area */
449:     PetscCall(PetscMalloc1(dgmres->max_k, &dgmres->nrs));
450:   }
451:   PetscCall(KSPDGMRESBuildSoln(dgmres->nrs, ksp->vec_sol, ptr, ksp, dgmres->it));
452:   if (result) *result = ptr;
453:   PetscFunctionReturn(PETSC_SUCCESS);
454: }

456: static PetscErrorCode KSPView_DGMRES(KSP ksp, PetscViewer viewer)
457: {
458:   KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
459:   PetscBool   iascii, isharmonic;

461:   PetscFunctionBegin;
462:   PetscCall(KSPView_GMRES(ksp, viewer));
463:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
464:   if (iascii) {
465:     PetscCall(PetscViewerASCIIPrintf(viewer, "    Adaptive strategy is used: %s\n", PetscBools[dgmres->force]));
466:     PetscCall(PetscOptionsHasName(((PetscObject)ksp)->options, ((PetscObject)ksp)->prefix, "-ksp_dgmres_harmonic_ritz", &isharmonic));
467:     if (isharmonic) {
468:       PetscCall(PetscViewerASCIIPrintf(viewer, "   Frequency of extracted eigenvalues = %" PetscInt_FMT " using Harmonic Ritz values \n", dgmres->neig));
469:     } else {
470:       PetscCall(PetscViewerASCIIPrintf(viewer, "   Frequency of extracted eigenvalues = %" PetscInt_FMT " using Ritz values \n", dgmres->neig));
471:     }
472:     PetscCall(PetscViewerASCIIPrintf(viewer, "   Total number of extracted eigenvalues = %" PetscInt_FMT "\n", dgmres->r));
473:     PetscCall(PetscViewerASCIIPrintf(viewer, "   Maximum number of eigenvalues set to be extracted = %" PetscInt_FMT "\n", dgmres->max_neig));
474:     PetscCall(PetscViewerASCIIPrintf(viewer, "   relaxation parameter for the adaptive strategy(smv)  = %g\n", (double)dgmres->smv));
475:     PetscCall(PetscViewerASCIIPrintf(viewer, "   Number of matvecs : %" PetscInt_FMT "\n", dgmres->matvecs));
476:   }
477:   PetscFunctionReturn(PETSC_SUCCESS);
478: }

480: PetscErrorCode KSPDGMRESSetEigen_DGMRES(KSP ksp, PetscInt neig)
481: {
482:   KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;

484:   PetscFunctionBegin;
485:   PetscCheck(neig >= 0 && neig <= dgmres->max_k, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "The value of neig must be positive and less than the restart value ");
486:   dgmres->neig = neig;
487:   PetscFunctionReturn(PETSC_SUCCESS);
488: }

490: static PetscErrorCode KSPDGMRESSetMaxEigen_DGMRES(KSP ksp, PetscInt max_neig)
491: {
492:   KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;

494:   PetscFunctionBegin;
495:   PetscCheck(max_neig >= 0 && max_neig <= dgmres->max_k, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "The value of max_neig must be positive and less than the restart value ");
496:   dgmres->max_neig = max_neig;
497:   PetscFunctionReturn(PETSC_SUCCESS);
498: }

500: static PetscErrorCode KSPDGMRESSetRatio_DGMRES(KSP ksp, PetscReal ratio)
501: {
502:   KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;

504:   PetscFunctionBegin;
505:   PetscCheck(ratio > 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "The relaxation parameter value must be positive");
506:   dgmres->smv = ratio;
507:   PetscFunctionReturn(PETSC_SUCCESS);
508: }

510: static PetscErrorCode KSPDGMRESForce_DGMRES(KSP ksp, PetscBool force)
511: {
512:   KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;

514:   PetscFunctionBegin;
515:   dgmres->force = force;
516:   PetscFunctionReturn(PETSC_SUCCESS);
517: }

519: PetscErrorCode KSPSetFromOptions_DGMRES(KSP ksp, PetscOptionItems *PetscOptionsObject)
520: {
521:   PetscInt    neig;
522:   PetscInt    max_neig;
523:   KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
524:   PetscBool   flg;

526:   PetscFunctionBegin;
527:   PetscCall(KSPSetFromOptions_GMRES(ksp, PetscOptionsObject));
528:   PetscOptionsHeadBegin(PetscOptionsObject, "KSP DGMRES Options");
529:   PetscCall(PetscOptionsInt("-ksp_dgmres_eigen", "Number of smallest eigenvalues to extract at each restart", "KSPDGMRESSetEigen", dgmres->neig, &neig, &flg));
530:   if (flg) PetscCall(KSPDGMRESSetEigen(ksp, neig));
531:   PetscCall(PetscOptionsInt("-ksp_dgmres_max_eigen", "Maximum Number of smallest eigenvalues to extract ", "KSPDGMRESSetMaxEigen", dgmres->max_neig, &max_neig, &flg));
532:   if (flg) PetscCall(KSPDGMRESSetMaxEigen(ksp, max_neig));
533:   PetscCall(PetscOptionsReal("-ksp_dgmres_ratio", "Relaxation parameter for the smaller number of matrix-vectors product allowed", "KSPDGMRESSetRatio", dgmres->smv, &dgmres->smv, NULL));
534:   PetscCall(PetscOptionsBool("-ksp_dgmres_improve", "Improve the computation of eigenvalues by solving a new generalized eigenvalue problem (experimental - not stable at this time)", NULL, dgmres->improve, &dgmres->improve, NULL));
535:   PetscCall(PetscOptionsBool("-ksp_dgmres_force", "Sets DGMRES always at restart active, i.e do not use the adaptive strategy", "KSPDGMRESForce", dgmres->force, &dgmres->force, NULL));
536:   PetscOptionsHeadEnd();
537:   PetscFunctionReturn(PETSC_SUCCESS);
538: }

540: PetscErrorCode KSPDGMRESComputeDeflationData_DGMRES(KSP ksp, PetscInt *ExtrNeig)
541: {
542:   KSP_DGMRES  *dgmres = (KSP_DGMRES *)ksp->data;
543:   PetscInt     i, j, k;
544:   PetscBLASInt nr, bmax;
545:   PetscInt     r = dgmres->r;
546:   PetscInt     neig;                                 /* number of eigenvalues to extract at each restart */
547:   PetscInt     neig1    = dgmres->neig + EIG_OFFSET; /* max number of eig that can be extracted at each restart */
548:   PetscInt     max_neig = dgmres->max_neig;          /* Max number of eigenvalues to extract during the iterative process */
549:   PetscInt     N        = dgmres->max_k + 1;
550:   PetscInt     n        = dgmres->it + 1;
551:   PetscReal    alpha;

553:   PetscFunctionBegin;
554:   PetscCall(PetscLogEventBegin(KSP_DGMRESComputeDeflationData, ksp, 0, 0, 0));
555:   if (dgmres->neig == 0 || (max_neig < (r + neig1) && !dgmres->improve)) {
556:     PetscCall(PetscLogEventEnd(KSP_DGMRESComputeDeflationData, ksp, 0, 0, 0));
557:     PetscFunctionReturn(PETSC_SUCCESS);
558:   }

560:   PetscCall(KSPDGMRESComputeSchurForm(ksp, &neig));
561:   /* Form the extended Schur vectors X=VV*Sr */
562:   if (!XX) PetscCall(VecDuplicateVecs(VEC_VV(0), neig1, &XX));
563:   for (j = 0; j < neig; j++) PetscCall(VecMAXPBY(XX[j], n, &SR[j * N], 0, &VEC_VV(0)));

565:   /* Orthogonalize X against U */
566:   if (!ORTH) PetscCall(PetscMalloc1(max_neig, &ORTH));
567:   if (r > 0) {
568:     /* modified Gram-Schmidt */
569:     for (j = 0; j < neig; j++) {
570:       for (i = 0; i < r; i++) {
571:         /* First, compute U'*X[j] */
572:         PetscCall(VecDot(XX[j], UU[i], &alpha));
573:         /* Then, compute X(j)=X(j)-U*U'*X(j) */
574:         PetscCall(VecAXPY(XX[j], -alpha, UU[i]));
575:       }
576:     }
577:   }
578:   /* Compute MX = M^{-1}*A*X */
579:   if (!MX) PetscCall(VecDuplicateVecs(VEC_VV(0), neig1, &MX));
580:   for (j = 0; j < neig; j++) PetscCall(KSP_PCApplyBAorAB(ksp, XX[j], MX[j], VEC_TEMP_MATOP));
581:   dgmres->matvecs += neig;

583:   if ((r + neig1) > max_neig && dgmres->improve) { /* Improve the approximate eigenvectors in X by solving a new generalized eigenvalue -- expensive to do this */
584:     PetscCall(KSPDGMRESImproveEig(ksp, neig));
585:     PetscCall(PetscLogEventEnd(KSP_DGMRESComputeDeflationData, ksp, 0, 0, 0));
586:     PetscFunctionReturn(PETSC_SUCCESS); /* We return here since data for M have been improved in KSPDGMRESImproveEig()*/
587:   }

589:   /* Compute XMX = X'*M^{-1}*A*X -- size (neig, neig) */
590:   if (!XMX) PetscCall(PetscMalloc1(neig1 * neig1, &XMX));
591:   for (j = 0; j < neig; j++) PetscCall(VecMDot(MX[j], neig, XX, &XMX[j * neig1]));

593:   if (r > 0) {
594:     /* Compute UMX = U'*M^{-1}*A*X -- size (r, neig) */
595:     if (!UMX) PetscCall(PetscMalloc1(max_neig * neig1, &UMX));
596:     for (j = 0; j < neig; j++) PetscCall(VecMDot(MX[j], r, UU, &UMX[j * max_neig]));
597:     /* Compute XMU = X'*M^{-1}*A*U -- size(neig, r) */
598:     if (!XMU) PetscCall(PetscMalloc1(max_neig * neig1, &XMU));
599:     for (j = 0; j < r; j++) PetscCall(VecMDot(MU[j], neig, XX, &XMU[j * neig1]));
600:   }

602:   /* Form the new matrix T = [T UMX; XMU XMX]; */
603:   if (!TT) PetscCall(PetscMalloc1(max_neig * max_neig, &TT));
604:   if (r > 0) {
605:     /* Add XMU to T */
606:     for (j = 0; j < r; j++) PetscCall(PetscArraycpy(&TT[max_neig * j + r], &XMU[neig1 * j], neig));
607:     /* Add [UMX; XMX] to T */
608:     for (j = 0; j < neig; j++) {
609:       k = r + j;
610:       PetscCall(PetscArraycpy(&TT[max_neig * k], &UMX[max_neig * j], r));
611:       PetscCall(PetscArraycpy(&TT[max_neig * k + r], &XMX[neig1 * j], neig));
612:     }
613:   } else { /* Add XMX to T */
614:     for (j = 0; j < neig; j++) PetscCall(PetscArraycpy(&TT[max_neig * j], &XMX[neig1 * j], neig));
615:   }

617:   dgmres->r += neig;
618:   r = dgmres->r;
619:   PetscCall(PetscBLASIntCast(r, &nr));
620:   /*LU Factorize T with Lapack xgetrf routine */

622:   PetscCall(PetscBLASIntCast(max_neig, &bmax));
623:   if (!TTF) PetscCall(PetscMalloc1(bmax * bmax, &TTF));
624:   PetscCall(PetscArraycpy(TTF, TT, bmax * r));
625:   if (!INVP) PetscCall(PetscMalloc1(bmax, &INVP));
626:   {
627:     PetscBLASInt info;
628:     PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&nr, &nr, TTF, &bmax, INVP, &info));
629:     PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XGETRF INFO=%d", (int)info);
630:   }

632:   /* Save X in U and MX in MU for the next cycles and increase the size of the invariant subspace */
633:   if (!UU) {
634:     PetscCall(VecDuplicateVecs(VEC_VV(0), max_neig, &UU));
635:     PetscCall(VecDuplicateVecs(VEC_VV(0), max_neig, &MU));
636:   }
637:   for (j = 0; j < neig; j++) {
638:     PetscCall(VecCopy(XX[j], UU[r - neig + j]));
639:     PetscCall(VecCopy(MX[j], MU[r - neig + j]));
640:   }
641:   PetscCall(PetscLogEventEnd(KSP_DGMRESComputeDeflationData, ksp, 0, 0, 0));
642:   PetscFunctionReturn(PETSC_SUCCESS);
643: }

645: PetscErrorCode KSPDGMRESComputeSchurForm_DGMRES(KSP ksp, PetscInt *neig)
646: {
647:   KSP_DGMRES   *dgmres = (KSP_DGMRES *)ksp->data;
648:   PetscInt      N = dgmres->max_k + 1, n = dgmres->it + 1;
649:   PetscBLASInt  bn;
650:   PetscReal    *A;
651:   PetscBLASInt  ihi;
652:   PetscBLASInt  ldA = 0; /* leading dimension of A */
653:   PetscBLASInt  ldQ;     /* leading dimension of Q */
654:   PetscReal    *Q;       /*  orthogonal matrix of  (left) Schur vectors */
655:   PetscReal    *work;    /* working vector */
656:   PetscBLASInt  lwork;   /* size of the working vector */
657:   PetscInt     *perm;    /* Permutation vector to sort eigenvalues */
658:   PetscInt      i, j;
659:   PetscBLASInt  NbrEig;          /* Number of eigenvalues really extracted */
660:   PetscReal    *wr, *wi, *modul; /* Real and imaginary part and modulus of the eigenvalues of A */
661:   PetscBLASInt *select;
662:   PetscBLASInt *iwork;
663:   PetscBLASInt  liwork;
664:   PetscScalar  *Ht;   /* Transpose of the Hessenberg matrix */
665:   PetscScalar  *t;    /* Store the result of the solution of H^T*t=h_{m+1,m}e_m */
666:   PetscBLASInt *ipiv; /* Permutation vector to be used in LAPACK */
667:   PetscBool     flag; /* determine whether to use Ritz vectors or harmonic Ritz vectors */

669:   PetscFunctionBegin;
670:   PetscCall(PetscBLASIntCast(n, &bn));
671:   PetscCall(PetscBLASIntCast(N, &ldA));
672:   ihi = ldQ = bn;
673:   PetscCall(PetscBLASIntCast(5 * N, &lwork));

675: #if defined(PETSC_USE_COMPLEX)
676:   SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "No support for complex numbers.");
677: #endif

679:   PetscCall(PetscMalloc1(ldA * ldA, &A));
680:   PetscCall(PetscMalloc1(ldQ * n, &Q));
681:   PetscCall(PetscMalloc1(lwork, &work));
682:   if (!dgmres->wr) {
683:     PetscCall(PetscMalloc1(n, &dgmres->wr));
684:     PetscCall(PetscMalloc1(n, &dgmres->wi));
685:   }
686:   wr = dgmres->wr;
687:   wi = dgmres->wi;
688:   PetscCall(PetscMalloc1(n, &modul));
689:   PetscCall(PetscMalloc1(n, &perm));
690:   /* copy the Hessenberg matrix to work space */
691:   PetscCall(PetscArraycpy(A, dgmres->hes_origin, ldA * ldA));
692:   PetscCall(PetscOptionsHasName(((PetscObject)ksp)->options, ((PetscObject)ksp)->prefix, "-ksp_dgmres_harmonic_ritz", &flag));
693:   if (flag) {
694:     /* Compute the matrix H + H^{-T}*h^2_{m+1,m}e_m*e_m^T */
695:     /* Transpose the Hessenberg matrix */
696:     PetscCall(PetscMalloc1(bn * bn, &Ht));
697:     for (i = 0; i < bn; i++) {
698:       for (j = 0; j < bn; j++) Ht[i * bn + j] = dgmres->hes_origin[j * ldA + i];
699:     }

701:     /* Solve the system H^T*t = h_{m+1,m}e_m */
702:     PetscCall(PetscCalloc1(bn, &t));
703:     t[bn - 1] = dgmres->hes_origin[(bn - 1) * ldA + bn]; /* Pick the last element H(m+1,m) */
704:     PetscCall(PetscMalloc1(bn, &ipiv));
705:     /* Call the LAPACK routine dgesv to solve the system Ht^-1 * t */
706:     {
707:       PetscBLASInt info;
708:       PetscBLASInt nrhs = 1;
709:       PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&bn, &nrhs, Ht, &bn, ipiv, t, &bn, &info));
710:       PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error while calling the Lapack routine DGESV");
711:     }
712:     /* Now form H + H^{-T}*h^2_{m+1,m}e_m*e_m^T */
713:     for (i = 0; i < bn; i++) A[(bn - 1) * bn + i] += t[i];
714:     PetscCall(PetscFree(t));
715:     PetscCall(PetscFree(Ht));
716:   }
717:   /* Compute eigenvalues with the Schur form */
718:   {
719:     PetscBLASInt info = 0;
720:     PetscBLASInt ilo  = 1;
721:     PetscCallBLAS("LAPACKhseqr", LAPACKhseqr_("S", "I", &bn, &ilo, &ihi, A, &ldA, wr, wi, Q, &ldQ, work, &lwork, &info));
722:     PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XHSEQR %d", (int)info);
723:   }
724:   PetscCall(PetscFree(work));

726:   /* sort the eigenvalues */
727:   for (i = 0; i < n; i++) modul[i] = PetscSqrtReal(wr[i] * wr[i] + wi[i] * wi[i]);
728:   for (i = 0; i < n; i++) perm[i] = i;

730:   PetscCall(PetscSortRealWithPermutation(n, modul, perm));
731:   /* save the complex modulus of the largest eigenvalue in magnitude */
732:   if (dgmres->lambdaN < modul[perm[n - 1]]) dgmres->lambdaN = modul[perm[n - 1]];
733:   /* count the number of extracted eigenvalues (with complex conjugates) */
734:   NbrEig = 0;
735:   while (NbrEig < dgmres->neig) {
736:     if (wi[perm[NbrEig]] != 0) NbrEig += 2;
737:     else NbrEig += 1;
738:   }
739:   /* Reorder the Schur decomposition so that the cluster of smallest eigenvalues appears in the leading diagonal blocks of A */

741:   PetscCall(PetscCalloc1(n, &select));

743:   if (!dgmres->GreatestEig) {
744:     for (j = 0; j < NbrEig; j++) select[perm[j]] = 1;
745:   } else {
746:     for (j = 0; j < NbrEig; j++) select[perm[n - j - 1]] = 1;
747:   }
748:   /* call Lapack dtrsen */
749:   lwork  = PetscMax(1, 4 * NbrEig * (bn - NbrEig));
750:   liwork = PetscMax(1, 2 * NbrEig * (bn - NbrEig));
751:   PetscCall(PetscMalloc1(lwork, &work));
752:   PetscCall(PetscMalloc1(liwork, &iwork));
753:   {
754:     PetscBLASInt info = 0;
755:     PetscReal    CondEig; /* lower bound on the reciprocal condition number for the selected cluster of eigenvalues */
756:     PetscReal    CondSub; /* estimated reciprocal condition number of the specified invariant subspace. */
757:     PetscCallBLAS("LAPACKtrsen", LAPACKtrsen_("B", "V", select, &bn, A, &ldA, Q, &ldQ, wr, wi, &NbrEig, &CondEig, &CondSub, work, &lwork, iwork, &liwork, &info));
758:     PetscCheck(info != 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Unable to reorder the eigenvalues with the LAPACK routine: ILL-CONDITIONED PROBLEM");
759:   }
760:   PetscCall(PetscFree(select));

762:   /* Extract the Schur vectors */
763:   for (j = 0; j < NbrEig; j++) PetscCall(PetscArraycpy(&SR[j * N], &Q[j * ldQ], n));
764:   *neig = NbrEig;
765:   PetscCall(PetscFree(A));
766:   PetscCall(PetscFree(work));
767:   PetscCall(PetscFree(perm));
768:   PetscCall(PetscFree(work));
769:   PetscCall(PetscFree(iwork));
770:   PetscCall(PetscFree(modul));
771:   PetscCall(PetscFree(Q));
772:   PetscFunctionReturn(PETSC_SUCCESS);
773: }

775: PetscErrorCode KSPDGMRESApplyDeflation_DGMRES(KSP ksp, Vec x, Vec y)
776: {
777:   KSP_DGMRES  *dgmres = (KSP_DGMRES *)ksp->data;
778:   PetscInt     i, r = dgmres->r;
779:   PetscReal    alpha    = 1.0;
780:   PetscInt     max_neig = dgmres->max_neig;
781:   PetscBLASInt br, bmax;
782:   PetscReal    lambda = dgmres->lambdaN;

784:   PetscFunctionBegin;
785:   PetscCall(PetscBLASIntCast(r, &br));
786:   PetscCall(PetscBLASIntCast(max_neig, &bmax));
787:   PetscCall(PetscLogEventBegin(KSP_DGMRESApplyDeflation, ksp, 0, 0, 0));
788:   if (!r) {
789:     PetscCall(VecCopy(x, y));
790:     PetscFunctionReturn(PETSC_SUCCESS);
791:   }
792:   /* Compute U'*x */
793:   if (!X1) {
794:     PetscCall(PetscMalloc1(bmax, &X1));
795:     PetscCall(PetscMalloc1(bmax, &X2));
796:   }
797:   PetscCall(VecMDot(x, r, UU, X1));

799:   /* Solve T*X1=X2 for X1*/
800:   PetscCall(PetscArraycpy(X2, X1, br));
801:   {
802:     PetscBLASInt info;
803:     PetscBLASInt nrhs = 1;
804:     PetscCallBLAS("LAPACKgetrs", LAPACKgetrs_("N", &br, &nrhs, TTF, &bmax, INVP, X1, &bmax, &info));
805:     PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XGETRS %d", (int)info);
806:   }
807:   /* Iterative refinement -- is it really necessary ?? */
808:   if (!WORK) {
809:     PetscCall(PetscMalloc1(3 * bmax, &WORK));
810:     PetscCall(PetscMalloc1(bmax, &IWORK));
811:   }
812:   {
813:     PetscBLASInt info;
814:     PetscReal    berr, ferr;
815:     PetscBLASInt nrhs = 1;
816:     PetscCallBLAS("LAPACKgerfs", LAPACKgerfs_("N", &br, &nrhs, TT, &bmax, TTF, &bmax, INVP, X2, &bmax, X1, &bmax, &ferr, &berr, WORK, IWORK, &info));
817:     PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XGERFS %d", (int)info);
818:   }

820:   for (i = 0; i < r; i++) X2[i] = X1[i] / lambda - X2[i];

822:   /* Compute X2=U*X2 */
823:   PetscCall(VecMAXPBY(y, r, X2, 0, UU));
824:   PetscCall(VecAXPY(y, alpha, x));

826:   PetscCall(PetscLogEventEnd(KSP_DGMRESApplyDeflation, ksp, 0, 0, 0));
827:   PetscFunctionReturn(PETSC_SUCCESS);
828: }

830: static PetscErrorCode KSPDGMRESImproveEig_DGMRES(KSP ksp, PetscInt neig)
831: {
832:   KSP_DGMRES   *dgmres = (KSP_DGMRES *)ksp->data;
833:   PetscInt      j, r_old, r = dgmres->r;
834:   PetscBLASInt  i     = 0;
835:   PetscInt      neig1 = dgmres->neig + EIG_OFFSET;
836:   PetscInt      bmax  = dgmres->max_neig;
837:   PetscInt      aug   = r + neig;       /* actual size of the augmented invariant basis */
838:   PetscInt      aug1  = bmax + neig1;   /* maximum size of the augmented invariant basis */
839:   PetscBLASInt  ldA;                    /* leading dimension of AUAU and AUU*/
840:   PetscBLASInt  N;                      /* size of AUAU */
841:   PetscReal    *Q;                      /*  orthogonal matrix of  (left) schur vectors */
842:   PetscReal    *Z;                      /*  orthogonal matrix of  (right) schur vectors */
843:   PetscReal    *work;                   /* working vector */
844:   PetscBLASInt  lwork;                  /* size of the working vector */
845:   PetscInt     *perm;                   /* Permutation vector to sort eigenvalues */
846:   PetscReal    *wr, *wi, *beta, *modul; /* Real and imaginary part and modulus of the eigenvalues of A*/
847:   PetscBLASInt  NbrEig = 0, nr, bm;
848:   PetscBLASInt *select;
849:   PetscBLASInt  liwork, *iwork;

851:   PetscFunctionBegin;
852:   /* Block construction of the matrices AUU=(AU)'*U and (AU)'*AU*/
853:   if (!AUU) {
854:     PetscCall(PetscMalloc1(aug1 * aug1, &AUU));
855:     PetscCall(PetscMalloc1(aug1 * aug1, &AUAU));
856:   }
857:   /* AUU = (AU)'*U = [(MU)'*U (MU)'*X; (MX)'*U (MX)'*X]
858:    * Note that MU and MX have been computed previously either in ComputeDataDeflation() or down here in a previous call to this function */
859:   /* (MU)'*U size (r x r) -- store in the <r> first columns of AUU*/
860:   for (j = 0; j < r; j++) PetscCall(VecMDot(UU[j], r, MU, &AUU[j * aug1]));
861:   /* (MU)'*X size (r x neig) -- store in AUU from the column <r>*/
862:   for (j = 0; j < neig; j++) PetscCall(VecMDot(XX[j], r, MU, &AUU[(r + j) * aug1]));
863:   /* (MX)'*U size (neig x r) -- store in the <r> first columns of AUU from the row <r>*/
864:   for (j = 0; j < r; j++) PetscCall(VecMDot(UU[j], neig, MX, &AUU[j * aug1 + r]));
865:   /* (MX)'*X size (neig neig) --  store in AUU from the column <r> and the row <r>*/
866:   for (j = 0; j < neig; j++) PetscCall(VecMDot(XX[j], neig, MX, &AUU[(r + j) * aug1 + r]));

868:   /* AUAU = (AU)'*AU = [(MU)'*MU (MU)'*MX; (MX)'*MU (MX)'*MX] */
869:   /* (MU)'*MU size (r x r) -- store in the <r> first columns of AUAU*/
870:   for (j = 0; j < r; j++) PetscCall(VecMDot(MU[j], r, MU, &AUAU[j * aug1]));
871:   /* (MU)'*MX size (r x neig) -- store in AUAU from the column <r>*/
872:   for (j = 0; j < neig; j++) PetscCall(VecMDot(MX[j], r, MU, &AUAU[(r + j) * aug1]));
873:   /* (MX)'*MU size (neig x r) -- store in the <r> first columns of AUAU from the row <r>*/
874:   for (j = 0; j < r; j++) PetscCall(VecMDot(MU[j], neig, MX, &AUAU[j * aug1 + r]));
875:   /* (MX)'*MX size (neig neig) --  store in AUAU from the column <r> and the row <r>*/
876:   for (j = 0; j < neig; j++) PetscCall(VecMDot(MX[j], neig, MX, &AUAU[(r + j) * aug1 + r]));

878:   /* Computation of the eigenvectors */
879:   PetscCall(PetscBLASIntCast(aug1, &ldA));
880:   PetscCall(PetscBLASIntCast(aug, &N));
881:   lwork = 8 * N + 20; /* sizeof the working space */
882:   PetscCall(PetscMalloc1(N, &wr));
883:   PetscCall(PetscMalloc1(N, &wi));
884:   PetscCall(PetscMalloc1(N, &beta));
885:   PetscCall(PetscMalloc1(N, &modul));
886:   PetscCall(PetscMalloc1(N, &perm));
887:   PetscCall(PetscMalloc1(N * N, &Q));
888:   PetscCall(PetscMalloc1(N * N, &Z));
889:   PetscCall(PetscMalloc1(lwork, &work));
890:   {
891:     PetscBLASInt info = 0;
892:     PetscCallBLAS("LAPACKgges", LAPACKgges_("V", "V", "N", NULL, &N, AUAU, &ldA, AUU, &ldA, &i, wr, wi, beta, Q, &N, Z, &N, work, &lwork, NULL, &info));
893:     PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XGGES %d", (int)info);
894:   }
895:   for (i = 0; i < N; i++) {
896:     if (beta[i] != 0.0) {
897:       wr[i] /= beta[i];
898:       wi[i] /= beta[i];
899:     }
900:   }
901:   /* sort the eigenvalues */
902:   for (i = 0; i < N; i++) modul[i] = PetscSqrtReal(wr[i] * wr[i] + wi[i] * wi[i]);
903:   for (i = 0; i < N; i++) perm[i] = i;
904:   PetscCall(PetscSortRealWithPermutation(N, modul, perm));
905:   /* Save the norm of the largest eigenvalue */
906:   if (dgmres->lambdaN < modul[perm[N - 1]]) dgmres->lambdaN = modul[perm[N - 1]];
907:   /* Allocate space to extract the first r schur vectors   */
908:   if (!SR2) PetscCall(PetscMalloc1(aug1 * bmax, &SR2));
909:   /* count the number of extracted eigenvalues (complex conjugates count as 2) */
910:   while (NbrEig < bmax) {
911:     if (wi[perm[NbrEig]] == 0) NbrEig += 1;
912:     else NbrEig += 2;
913:   }
914:   if (NbrEig > bmax) NbrEig = bmax - 1;
915:   r_old     = r; /* previous size of r */
916:   dgmres->r = r = NbrEig;

918:   /* Select the eigenvalues to reorder */
919:   PetscCall(PetscCalloc1(N, &select));
920:   if (!dgmres->GreatestEig) {
921:     for (j = 0; j < NbrEig; j++) select[perm[j]] = 1;
922:   } else {
923:     for (j = 0; j < NbrEig; j++) select[perm[N - j - 1]] = 1;
924:   }
925:   /* Reorder and extract the new <r> schur vectors */
926:   lwork  = PetscMax(4 * N + 16, 2 * NbrEig * (N - NbrEig));
927:   liwork = PetscMax(N + 6, 2 * NbrEig * (N - NbrEig));
928:   PetscCall(PetscFree(work));
929:   PetscCall(PetscMalloc1(lwork, &work));
930:   PetscCall(PetscMalloc1(liwork, &iwork));
931:   {
932:     PetscBLASInt info = 0;
933:     PetscReal    Dif[2];
934:     PetscBLASInt ijob  = 2;
935:     PetscBLASInt wantQ = 1, wantZ = 1;
936:     PetscCallBLAS("LAPACKtgsen", LAPACKtgsen_(&ijob, &wantQ, &wantZ, select, &N, AUAU, &ldA, AUU, &ldA, wr, wi, beta, Q, &N, Z, &N, &NbrEig, NULL, NULL, &Dif[0], work, &lwork, iwork, &liwork, &info));
937:     PetscCheck(info != 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Unable to reorder the eigenvalues with the LAPACK routine: ill-conditioned problem.");
938:   }
939:   PetscCall(PetscFree(select));

941:   for (j = 0; j < r; j++) PetscCall(PetscArraycpy(&SR2[j * aug1], &Z[j * N], N));

943:   /* Multiply the Schur vectors SR2 by U (and X)  to get a new U
944:    -- save it temporarily in MU */
945:   for (j = 0; j < r; j++) {
946:     PetscCall(VecMAXPBY(MU[j], r_old, &SR2[j * aug1], 0, UU));
947:     PetscCall(VecMAXPY(MU[j], neig, &SR2[j * aug1 + r_old], XX));
948:   }
949:   /* Form T = U'*MU*U */
950:   for (j = 0; j < r; j++) {
951:     PetscCall(VecCopy(MU[j], UU[j]));
952:     PetscCall(KSP_PCApplyBAorAB(ksp, UU[j], MU[j], VEC_TEMP_MATOP));
953:   }
954:   dgmres->matvecs += r;
955:   for (j = 0; j < r; j++) PetscCall(VecMDot(MU[j], r, UU, &TT[j * bmax]));
956:   /* Factorize T */
957:   PetscCall(PetscArraycpy(TTF, TT, bmax * r));
958:   PetscCall(PetscBLASIntCast(r, &nr));
959:   PetscCall(PetscBLASIntCast(bmax, &bm));
960:   {
961:     PetscBLASInt info;
962:     PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&nr, &nr, TTF, &bm, INVP, &info));
963:     PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XGETRF INFO=%d", (int)info);
964:   }
965:   /* Free Memory */
966:   PetscCall(PetscFree(wr));
967:   PetscCall(PetscFree(wi));
968:   PetscCall(PetscFree(beta));
969:   PetscCall(PetscFree(modul));
970:   PetscCall(PetscFree(perm));
971:   PetscCall(PetscFree(Q));
972:   PetscCall(PetscFree(Z));
973:   PetscCall(PetscFree(work));
974:   PetscCall(PetscFree(iwork));
975:   PetscFunctionReturn(PETSC_SUCCESS);
976: }

978: /*MC
979:    KSPDGMRES - Implements the deflated GMRES as defined in {cite}`erhel1996restarted` and {cite}`wakam2013memory`

981:    Options Database Keys:
982:    GMRES Options (inherited):
983: +   -ksp_gmres_restart <restart>                                                - the number of Krylov directions to orthogonalize against
984: .   -ksp_gmres_haptol <tol>                                                     - sets the tolerance for "happy ending" (exact convergence)
985: .   -ksp_gmres_preallocate                                                      - preallocate all the Krylov search directions initially
986:                                                                                 (otherwise groups of vectors are allocated as needed)
987: .   -ksp_gmres_classicalgramschmidt                                             - use classical (unmodified) Gram-Schmidt to orthogonalize against
988:                                                                                 the Krylov space (fast) (the default)
989: .   -ksp_gmres_modifiedgramschmidt                                              - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
990: .   -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
991:                                                                                 stability of the classical Gram-Schmidt orthogonalization.
992: -   -ksp_gmres_krylov_monitor                                                   - plot the Krylov space generated

994:    DGMRES Options Database Keys:
995: +   -ksp_dgmres_eigen <neig>                     - number of smallest eigenvalues to extract at each restart
996: .   -ksp_dgmres_max_eigen <max_neig>             - maximum number of eigenvalues that can be extracted during the iterative process
997: .   -ksp_dgmres_force                            - use the deflation at each restart; switch off the adaptive strategy.
998: -   -ksp_dgmres_view_deflation_vecs <viewerspec> - View the deflation vectors, where viewerspec is a key that can be
999:                                                    parsed by `PetscOptionsGetViewer()`.  If neig > 1, viewerspec should
1000:                                                    end with ":append".  No vectors will be viewed if the adaptive
1001:                                                    strategy chooses not to deflate, so -ksp_dgmres_force should also
1002:                                                    be given.
1003:                                                    The deflation vectors span a subspace that may be a good
1004:                                                    approximation of the subspace of smallest eigenvectors of the
1005:                                                    preconditioned operator, so this option can aid in understanding
1006:                                                    the performance of a preconditioner.

1008:    Level: beginner

1010:    Notes:
1011:    Left and right preconditioning are supported, but not symmetric preconditioning. Complex arithmetic is not supported

1013:    In this implementation, the adaptive strategy allows switching to deflated GMRES when the stagnation occurs.

1015:    Contributed by:
1016:    Desire NUENTSA WAKAM, INRIA

1018: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPFGMRES`, `KSPLGMRES`,
1019:            `KSPGMRESSetRestart()`, `KSPGMRESSetHapTol()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
1020:            `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`,
1021:            `KSPGMRESCGSRefinementType`, `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESMonitorKrylov()`, `KSPSetPCSide()`
1022:  M*/

1024: PETSC_EXTERN PetscErrorCode KSPCreate_DGMRES(KSP ksp)
1025: {
1026:   KSP_DGMRES *dgmres;

1028:   PetscFunctionBegin;
1029:   PetscCall(PetscNew(&dgmres));
1030:   ksp->data = (void *)dgmres;

1032:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3));
1033:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2));
1034:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1));

1036:   ksp->ops->buildsolution                = KSPBuildSolution_DGMRES;
1037:   ksp->ops->setup                        = KSPSetUp_DGMRES;
1038:   ksp->ops->solve                        = KSPSolve_DGMRES;
1039:   ksp->ops->destroy                      = KSPDestroy_DGMRES;
1040:   ksp->ops->view                         = KSPView_DGMRES;
1041:   ksp->ops->setfromoptions               = KSPSetFromOptions_DGMRES;
1042:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
1043:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

1045:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetPreAllocateVectors_C", KSPGMRESSetPreAllocateVectors_GMRES));
1046:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetOrthogonalization_C", KSPGMRESSetOrthogonalization_GMRES));
1047:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetRestart_C", KSPGMRESSetRestart_GMRES));
1048:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetHapTol_C", KSPGMRESSetHapTol_GMRES));
1049:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetCGSRefinementType_C", KSPGMRESSetCGSRefinementType_GMRES));
1050:   /* -- New functions defined in DGMRES -- */
1051:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetEigen_C", KSPDGMRESSetEigen_DGMRES));
1052:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetMaxEigen_C", KSPDGMRESSetMaxEigen_DGMRES));
1053:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetRatio_C", KSPDGMRESSetRatio_DGMRES));
1054:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESForce_C", KSPDGMRESForce_DGMRES));
1055:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESComputeSchurForm_C", KSPDGMRESComputeSchurForm_DGMRES));
1056:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESComputeDeflationData_C", KSPDGMRESComputeDeflationData_DGMRES));
1057:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESApplyDeflation_C", KSPDGMRESApplyDeflation_DGMRES));
1058:   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESImproveEig_C", KSPDGMRESImproveEig_DGMRES));

1060:   PetscCall(PetscLogEventRegister("DGMRESCompDefl", KSP_CLASSID, &KSP_DGMRESComputeDeflationData));
1061:   PetscCall(PetscLogEventRegister("DGMRESApplyDefl", KSP_CLASSID, &KSP_DGMRESApplyDeflation));

1063:   dgmres->haptol         = 1.0e-30;
1064:   dgmres->q_preallocate  = 0;
1065:   dgmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
1066:   dgmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
1067:   dgmres->nrs            = NULL;
1068:   dgmres->sol_temp       = NULL;
1069:   dgmres->max_k          = GMRES_DEFAULT_MAXK;
1070:   dgmres->Rsvd           = NULL;
1071:   dgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
1072:   dgmres->orthogwork     = NULL;

1074:   /* Default values for the deflation */
1075:   dgmres->r           = 0;
1076:   dgmres->neig        = DGMRES_DEFAULT_EIG;
1077:   dgmres->max_neig    = DGMRES_DEFAULT_MAXEIG - 1;
1078:   dgmres->lambdaN     = 0.0;
1079:   dgmres->smv         = SMV;
1080:   dgmres->matvecs     = 0;
1081:   dgmres->GreatestEig = PETSC_FALSE; /* experimental */
1082:   dgmres->HasSchur    = PETSC_FALSE;
1083:   PetscFunctionReturn(PETSC_SUCCESS);
1084: }