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: /* Monitor if we know that we will not return for a restart */
197: if (it && (ksp->reason || ksp->its >= ksp->max_it)) {
198: PetscCall(KSPLogResidualHistory(ksp, ksp->rnorm));
199: PetscCall(KSPMonitor(ksp, ksp->its, ksp->rnorm));
200: }
201: if (itcount) *itcount = it;
203: /*
204: Down here we have to solve for the "best" coefficients of the Krylov
205: columns, add the solution values together, and possibly unwind the
206: preconditioning from the solution
207: */
208: /* Form the solution (or the solution so far) */
209: PetscCall(KSPDGMRESBuildSoln(GRS(0), ksp->vec_sol, ksp->vec_sol, ksp, it - 1));
211: /* Compute data for the deflation to be used during the next restart */
212: if (!ksp->reason && ksp->its < ksp->max_it) {
213: test = max_k * PetscLogReal(ksp->rtol / res) / PetscLogReal(res / res_old);
214: /* Compute data for the deflation if the residual rtol will not be reached in the remaining number of steps allowed */
215: if ((test > dgmres->smv * (ksp->max_it - ksp->its)) || dgmres->force) PetscCall(KSPDGMRESComputeDeflationData(ksp, NULL));
216: }
217: PetscFunctionReturn(PETSC_SUCCESS);
218: }
220: PetscErrorCode KSPSolve_DGMRES(KSP ksp)
221: {
222: PetscInt i, its = 0, itcount;
223: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
224: PetscBool guess_zero = ksp->guess_zero;
226: PetscFunctionBegin;
227: PetscCheck(!ksp->calc_sings || dgmres->Rsvd, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ORDER, "Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
229: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)ksp));
230: ksp->its = 0;
231: dgmres->matvecs = 0;
232: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)ksp));
234: itcount = 0;
235: while (!ksp->reason) {
236: PetscCall(KSPInitialResidual(ksp, ksp->vec_sol, VEC_TEMP, VEC_TEMP_MATOP, VEC_VV(0), ksp->vec_rhs));
237: if (ksp->pc_side == PC_LEFT) {
238: dgmres->matvecs += 1;
239: if (dgmres->r > 0) {
240: PetscCall(KSPDGMRESApplyDeflation(ksp, VEC_VV(0), VEC_TEMP));
241: PetscCall(VecCopy(VEC_TEMP, VEC_VV(0)));
242: }
243: }
245: PetscCall(KSPDGMRESCycle(&its, ksp));
246: itcount += its;
247: if (itcount >= ksp->max_it) {
248: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
249: break;
250: }
251: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
252: }
253: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
255: for (i = 0; i < dgmres->r; i++) PetscCall(VecViewFromOptions(UU[i], (PetscObject)ksp, "-ksp_dgmres_view_deflation_vecs"));
256: PetscFunctionReturn(PETSC_SUCCESS);
257: }
259: PetscErrorCode KSPDestroy_DGMRES(KSP ksp)
260: {
261: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
262: PetscInt neig1 = dgmres->neig + EIG_OFFSET;
263: PetscInt max_neig = dgmres->max_neig;
265: PetscFunctionBegin;
266: if (dgmres->r) {
267: PetscCall(VecDestroyVecs(max_neig, &UU));
268: PetscCall(VecDestroyVecs(max_neig, &MU));
269: if (XX) {
270: PetscCall(VecDestroyVecs(neig1, &XX));
271: PetscCall(VecDestroyVecs(neig1, &MX));
272: }
273: PetscCall(PetscFree(TT));
274: PetscCall(PetscFree(TTF));
275: PetscCall(PetscFree(INVP));
276: PetscCall(PetscFree(XMX));
277: PetscCall(PetscFree(UMX));
278: PetscCall(PetscFree(XMU));
279: PetscCall(PetscFree(X1));
280: PetscCall(PetscFree(X2));
281: PetscCall(PetscFree(dgmres->work));
282: PetscCall(PetscFree(dgmres->iwork));
283: PetscCall(PetscFree(dgmres->wr));
284: PetscCall(PetscFree(dgmres->wi));
285: PetscCall(PetscFree(dgmres->modul));
286: PetscCall(PetscFree(dgmres->Q));
287: PetscCall(PetscFree(ORTH));
288: PetscCall(PetscFree(AUAU));
289: PetscCall(PetscFree(AUU));
290: PetscCall(PetscFree(SR2));
291: }
292: PetscCall(PetscFree(SR));
293: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetEigen_C", NULL));
294: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetMaxEigen_C", NULL));
295: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetRatio_C", NULL));
296: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESForce_C", NULL));
297: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESComputeSchurForm_C", NULL));
298: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESComputeDeflationData_C", NULL));
299: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESApplyDeflation_C", NULL));
300: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESImproveEig_C", NULL));
301: PetscCall(KSPDestroy_GMRES(ksp));
302: PetscFunctionReturn(PETSC_SUCCESS);
303: }
305: /*
306: KSPDGMRESBuildSoln - create the solution from the starting vector and the
307: current iterates.
309: Input parameters:
310: nrs - work area of size it + 1.
311: vs - index of initial guess
312: vdest - index of result. Note that vs may == vdest (replace
313: guess with the solution).
315: This is an internal routine that knows about the GMRES internals.
316: */
317: static PetscErrorCode KSPDGMRESBuildSoln(PetscScalar *nrs, Vec vs, Vec vdest, KSP ksp, PetscInt it)
318: {
319: PetscScalar tt;
320: PetscInt ii, k, j;
321: KSP_DGMRES *dgmres = (KSP_DGMRES *)(ksp->data);
323: /* Solve for solution vector that minimizes the residual */
325: PetscFunctionBegin;
326: /* If it is < 0, no gmres steps have been performed */
327: if (it < 0) {
328: PetscCall(VecCopy(vs, vdest)); /* VecCopy() is smart, exists immediately if vguess == vdest */
329: PetscFunctionReturn(PETSC_SUCCESS);
330: }
331: 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)));
332: if (*HH(it, it) != 0.0) nrs[it] = *GRS(it) / *HH(it, it);
333: else nrs[it] = 0.0;
335: for (ii = 1; ii <= it; ii++) {
336: k = it - ii;
337: tt = *GRS(k);
338: for (j = k + 1; j <= it; j++) tt = tt - *HH(k, j) * nrs[j];
339: 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);
340: nrs[k] = tt / *HH(k, k);
341: }
343: /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
344: PetscCall(VecMAXPBY(VEC_TEMP, it + 1, nrs, 0, &VEC_VV(0)));
346: /* Apply deflation */
347: if (ksp->pc_side == PC_RIGHT && dgmres->r > 0) {
348: PetscCall(KSPDGMRESApplyDeflation(ksp, VEC_TEMP, VEC_TEMP_MATOP));
349: PetscCall(VecCopy(VEC_TEMP_MATOP, VEC_TEMP));
350: }
351: PetscCall(KSPUnwindPreconditioner(ksp, VEC_TEMP, VEC_TEMP_MATOP));
353: /* add solution to previous solution */
354: if (vdest != vs) PetscCall(VecCopy(vs, vdest));
355: PetscCall(VecAXPY(vdest, 1.0, VEC_TEMP));
356: PetscFunctionReturn(PETSC_SUCCESS);
357: }
359: /*
360: Do the scalar work for the orthogonalization. Return new residual norm.
361: */
362: static PetscErrorCode KSPDGMRESUpdateHessenberg(KSP ksp, PetscInt it, PetscBool hapend, PetscReal *res)
363: {
364: PetscScalar *hh, *cc, *ss, tt;
365: PetscInt j;
366: KSP_DGMRES *dgmres = (KSP_DGMRES *)(ksp->data);
368: PetscFunctionBegin;
369: hh = HH(0, it);
370: cc = CC(0);
371: ss = SS(0);
373: /* Apply all the previously computed plane rotations to the new column
374: of the Hessenberg matrix */
375: for (j = 1; j <= it; j++) {
376: tt = *hh;
377: *hh = PetscConj(*cc) * tt + *ss * *(hh + 1);
378: hh++;
379: *hh = *cc++ * *hh - (*ss++ * tt);
380: }
382: /*
383: compute the new plane rotation, and apply it to:
384: 1) the right-hand-side of the Hessenberg system
385: 2) the new column of the Hessenberg matrix
386: thus obtaining the updated value of the residual
387: */
388: if (!hapend) {
389: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh + 1)) * *(hh + 1));
390: if (tt == 0.0) {
391: ksp->reason = KSP_DIVERGED_NULL;
392: PetscFunctionReturn(PETSC_SUCCESS);
393: }
394: *cc = *hh / tt;
395: *ss = *(hh + 1) / tt;
396: *GRS(it + 1) = -(*ss * *GRS(it));
397: *GRS(it) = PetscConj(*cc) * *GRS(it);
398: *hh = PetscConj(*cc) * *hh + *ss * *(hh + 1);
399: *res = PetscAbsScalar(*GRS(it + 1));
400: } else {
401: /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
402: another rotation matrix (so RH doesn't change). The new residual is
403: always the new sine term times the residual from last time (GRS(it)),
404: but now the new sine rotation would be zero...so the residual should
405: be zero...so we will multiply "zero" by the last residual. This might
406: not be exactly what we want to do here -could just return "zero". */
407: *res = 0.0;
408: }
409: PetscFunctionReturn(PETSC_SUCCESS);
410: }
412: /*
413: Allocates more work vectors, starting from VEC_VV(it).
414: */
415: static PetscErrorCode KSPDGMRESGetNewVectors(KSP ksp, PetscInt it)
416: {
417: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
418: PetscInt nwork = dgmres->nwork_alloc, k, nalloc;
420: PetscFunctionBegin;
421: nalloc = PetscMin(ksp->max_it, dgmres->delta_allocate);
422: /* Adjust the number to allocate to make sure that we don't exceed the
423: number of available slots */
424: if (it + VEC_OFFSET + nalloc >= dgmres->vecs_allocated) nalloc = dgmres->vecs_allocated - it - VEC_OFFSET;
425: if (!nalloc) PetscFunctionReturn(PETSC_SUCCESS);
427: dgmres->vv_allocated += nalloc;
429: PetscCall(KSPCreateVecs(ksp, nalloc, &dgmres->user_work[nwork], 0, NULL));
431: dgmres->mwork_alloc[nwork] = nalloc;
432: for (k = 0; k < nalloc; k++) dgmres->vecs[it + VEC_OFFSET + k] = dgmres->user_work[nwork][k];
433: dgmres->nwork_alloc++;
434: PetscFunctionReturn(PETSC_SUCCESS);
435: }
437: PetscErrorCode KSPBuildSolution_DGMRES(KSP ksp, Vec ptr, Vec *result)
438: {
439: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
441: PetscFunctionBegin;
442: if (!ptr) {
443: if (!dgmres->sol_temp) PetscCall(VecDuplicate(ksp->vec_sol, &dgmres->sol_temp));
444: ptr = dgmres->sol_temp;
445: }
446: if (!dgmres->nrs) {
447: /* allocate the work area */
448: PetscCall(PetscMalloc1(dgmres->max_k, &dgmres->nrs));
449: }
450: PetscCall(KSPDGMRESBuildSoln(dgmres->nrs, ksp->vec_sol, ptr, ksp, dgmres->it));
451: if (result) *result = ptr;
452: PetscFunctionReturn(PETSC_SUCCESS);
453: }
455: static PetscErrorCode KSPView_DGMRES(KSP ksp, PetscViewer viewer)
456: {
457: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
458: PetscBool iascii, isharmonic;
460: PetscFunctionBegin;
461: PetscCall(KSPView_GMRES(ksp, viewer));
462: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
463: if (iascii) {
464: PetscCall(PetscViewerASCIIPrintf(viewer, " Adaptive strategy is used: %s\n", PetscBools[dgmres->force]));
465: PetscCall(PetscOptionsHasName(((PetscObject)ksp)->options, ((PetscObject)ksp)->prefix, "-ksp_dgmres_harmonic_ritz", &isharmonic));
466: if (isharmonic) {
467: PetscCall(PetscViewerASCIIPrintf(viewer, " Frequency of extracted eigenvalues = %" PetscInt_FMT " using Harmonic Ritz values \n", dgmres->neig));
468: } else {
469: PetscCall(PetscViewerASCIIPrintf(viewer, " Frequency of extracted eigenvalues = %" PetscInt_FMT " using Ritz values \n", dgmres->neig));
470: }
471: PetscCall(PetscViewerASCIIPrintf(viewer, " Total number of extracted eigenvalues = %" PetscInt_FMT "\n", dgmres->r));
472: PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum number of eigenvalues set to be extracted = %" PetscInt_FMT "\n", dgmres->max_neig));
473: PetscCall(PetscViewerASCIIPrintf(viewer, " relaxation parameter for the adaptive strategy(smv) = %g\n", (double)dgmres->smv));
474: PetscCall(PetscViewerASCIIPrintf(viewer, " Number of matvecs : %" PetscInt_FMT "\n", dgmres->matvecs));
475: }
476: PetscFunctionReturn(PETSC_SUCCESS);
477: }
479: PetscErrorCode KSPDGMRESSetEigen_DGMRES(KSP ksp, PetscInt neig)
480: {
481: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
483: PetscFunctionBegin;
484: 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 ");
485: dgmres->neig = neig;
486: PetscFunctionReturn(PETSC_SUCCESS);
487: }
489: static PetscErrorCode KSPDGMRESSetMaxEigen_DGMRES(KSP ksp, PetscInt max_neig)
490: {
491: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
493: PetscFunctionBegin;
494: 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 ");
495: dgmres->max_neig = max_neig;
496: PetscFunctionReturn(PETSC_SUCCESS);
497: }
499: static PetscErrorCode KSPDGMRESSetRatio_DGMRES(KSP ksp, PetscReal ratio)
500: {
501: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
503: PetscFunctionBegin;
504: PetscCheck(ratio > 0, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "The relaxation parameter value must be positive");
505: dgmres->smv = ratio;
506: PetscFunctionReturn(PETSC_SUCCESS);
507: }
509: static PetscErrorCode KSPDGMRESForce_DGMRES(KSP ksp, PetscBool force)
510: {
511: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
513: PetscFunctionBegin;
514: dgmres->force = force;
515: PetscFunctionReturn(PETSC_SUCCESS);
516: }
518: PetscErrorCode KSPSetFromOptions_DGMRES(KSP ksp, PetscOptionItems *PetscOptionsObject)
519: {
520: PetscInt neig;
521: PetscInt max_neig;
522: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
523: PetscBool flg;
525: PetscFunctionBegin;
526: PetscCall(KSPSetFromOptions_GMRES(ksp, PetscOptionsObject));
527: PetscOptionsHeadBegin(PetscOptionsObject, "KSP DGMRES Options");
528: PetscCall(PetscOptionsInt("-ksp_dgmres_eigen", "Number of smallest eigenvalues to extract at each restart", "KSPDGMRESSetEigen", dgmres->neig, &neig, &flg));
529: if (flg) PetscCall(KSPDGMRESSetEigen(ksp, neig));
530: PetscCall(PetscOptionsInt("-ksp_dgmres_max_eigen", "Maximum Number of smallest eigenvalues to extract ", "KSPDGMRESSetMaxEigen", dgmres->max_neig, &max_neig, &flg));
531: if (flg) PetscCall(KSPDGMRESSetMaxEigen(ksp, max_neig));
532: PetscCall(PetscOptionsReal("-ksp_dgmres_ratio", "Relaxation parameter for the smaller number of matrix-vectors product allowed", "KSPDGMRESSetRatio", dgmres->smv, &dgmres->smv, NULL));
533: 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));
534: 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));
535: PetscOptionsHeadEnd();
536: PetscFunctionReturn(PETSC_SUCCESS);
537: }
539: PetscErrorCode KSPDGMRESComputeDeflationData_DGMRES(KSP ksp, PetscInt *ExtrNeig)
540: {
541: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
542: PetscInt i, j, k;
543: PetscBLASInt nr, bmax;
544: PetscInt r = dgmres->r;
545: PetscInt neig; /* number of eigenvalues to extract at each restart */
546: PetscInt neig1 = dgmres->neig + EIG_OFFSET; /* max number of eig that can be extracted at each restart */
547: PetscInt max_neig = dgmres->max_neig; /* Max number of eigenvalues to extract during the iterative process */
548: PetscInt N = dgmres->max_k + 1;
549: PetscInt n = dgmres->it + 1;
550: PetscReal alpha;
552: PetscFunctionBegin;
553: PetscCall(PetscLogEventBegin(KSP_DGMRESComputeDeflationData, ksp, 0, 0, 0));
554: if (dgmres->neig == 0 || (max_neig < (r + neig1) && !dgmres->improve)) {
555: PetscCall(PetscLogEventEnd(KSP_DGMRESComputeDeflationData, ksp, 0, 0, 0));
556: PetscFunctionReturn(PETSC_SUCCESS);
557: }
559: PetscCall(KSPDGMRESComputeSchurForm(ksp, &neig));
560: /* Form the extended Schur vectors X=VV*Sr */
561: if (!XX) PetscCall(VecDuplicateVecs(VEC_VV(0), neig1, &XX));
562: for (j = 0; j < neig; j++) PetscCall(VecMAXPBY(XX[j], n, &SR[j * N], 0, &VEC_VV(0)));
564: /* Orthogonalize X against U */
565: if (!ORTH) PetscCall(PetscMalloc1(max_neig, &ORTH));
566: if (r > 0) {
567: /* modified Gram-Schmidt */
568: for (j = 0; j < neig; j++) {
569: for (i = 0; i < r; i++) {
570: /* First, compute U'*X[j] */
571: PetscCall(VecDot(XX[j], UU[i], &alpha));
572: /* Then, compute X(j)=X(j)-U*U'*X(j) */
573: PetscCall(VecAXPY(XX[j], -alpha, UU[i]));
574: }
575: }
576: }
577: /* Compute MX = M^{-1}*A*X */
578: if (!MX) PetscCall(VecDuplicateVecs(VEC_VV(0), neig1, &MX));
579: for (j = 0; j < neig; j++) PetscCall(KSP_PCApplyBAorAB(ksp, XX[j], MX[j], VEC_TEMP_MATOP));
580: dgmres->matvecs += neig;
582: if ((r + neig1) > max_neig && dgmres->improve) { /* Improve the approximate eigenvectors in X by solving a new generalized eigenvalue -- expensive to do this */
583: PetscCall(KSPDGMRESImproveEig(ksp, neig));
584: PetscCall(PetscLogEventEnd(KSP_DGMRESComputeDeflationData, ksp, 0, 0, 0));
585: PetscFunctionReturn(PETSC_SUCCESS); /* We return here since data for M have been improved in KSPDGMRESImproveEig()*/
586: }
588: /* Compute XMX = X'*M^{-1}*A*X -- size (neig, neig) */
589: if (!XMX) PetscCall(PetscMalloc1(neig1 * neig1, &XMX));
590: for (j = 0; j < neig; j++) PetscCall(VecMDot(MX[j], neig, XX, &(XMX[j * neig1])));
592: if (r > 0) {
593: /* Compute UMX = U'*M^{-1}*A*X -- size (r, neig) */
594: if (!UMX) PetscCall(PetscMalloc1(max_neig * neig1, &UMX));
595: for (j = 0; j < neig; j++) PetscCall(VecMDot(MX[j], r, UU, &(UMX[j * max_neig])));
596: /* Compute XMU = X'*M^{-1}*A*U -- size(neig, r) */
597: if (!XMU) PetscCall(PetscMalloc1(max_neig * neig1, &XMU));
598: for (j = 0; j < r; j++) PetscCall(VecMDot(MU[j], neig, XX, &(XMU[j * neig1])));
599: }
601: /* Form the new matrix T = [T UMX; XMU XMX]; */
602: if (!TT) PetscCall(PetscMalloc1(max_neig * max_neig, &TT));
603: if (r > 0) {
604: /* Add XMU to T */
605: for (j = 0; j < r; j++) PetscCall(PetscArraycpy(&(TT[max_neig * j + r]), &(XMU[neig1 * j]), neig));
606: /* Add [UMX; XMX] to T */
607: for (j = 0; j < neig; j++) {
608: k = r + j;
609: PetscCall(PetscArraycpy(&(TT[max_neig * k]), &(UMX[max_neig * j]), r));
610: PetscCall(PetscArraycpy(&(TT[max_neig * k + r]), &(XMX[neig1 * j]), neig));
611: }
612: } else { /* Add XMX to T */
613: for (j = 0; j < neig; j++) PetscCall(PetscArraycpy(&(TT[max_neig * j]), &(XMX[neig1 * j]), neig));
614: }
616: dgmres->r += neig;
617: r = dgmres->r;
618: PetscCall(PetscBLASIntCast(r, &nr));
619: /*LU Factorize T with Lapack xgetrf routine */
621: PetscCall(PetscBLASIntCast(max_neig, &bmax));
622: if (!TTF) PetscCall(PetscMalloc1(bmax * bmax, &TTF));
623: PetscCall(PetscArraycpy(TTF, TT, bmax * r));
624: if (!INVP) PetscCall(PetscMalloc1(bmax, &INVP));
625: {
626: PetscBLASInt info;
627: PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&nr, &nr, TTF, &bmax, INVP, &info));
628: PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XGETRF INFO=%d", (int)info);
629: }
631: /* Save X in U and MX in MU for the next cycles and increase the size of the invariant subspace */
632: if (!UU) {
633: PetscCall(VecDuplicateVecs(VEC_VV(0), max_neig, &UU));
634: PetscCall(VecDuplicateVecs(VEC_VV(0), max_neig, &MU));
635: }
636: for (j = 0; j < neig; j++) {
637: PetscCall(VecCopy(XX[j], UU[r - neig + j]));
638: PetscCall(VecCopy(MX[j], MU[r - neig + j]));
639: }
640: PetscCall(PetscLogEventEnd(KSP_DGMRESComputeDeflationData, ksp, 0, 0, 0));
641: PetscFunctionReturn(PETSC_SUCCESS);
642: }
644: PetscErrorCode KSPDGMRESComputeSchurForm_DGMRES(KSP ksp, PetscInt *neig)
645: {
646: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
647: PetscInt N = dgmres->max_k + 1, n = dgmres->it + 1;
648: PetscBLASInt bn;
649: PetscReal *A;
650: PetscBLASInt ihi;
651: PetscBLASInt ldA = 0; /* leading dimension of A */
652: PetscBLASInt ldQ; /* leading dimension of Q */
653: PetscReal *Q; /* orthogonal matrix of (left) Schur vectors */
654: PetscReal *work; /* working vector */
655: PetscBLASInt lwork; /* size of the working vector */
656: PetscInt *perm; /* Permutation vector to sort eigenvalues */
657: PetscInt i, j;
658: PetscBLASInt NbrEig; /* Number of eigenvalues really extracted */
659: PetscReal *wr, *wi, *modul; /* Real and imaginary part and modulus of the eigenvalues of A */
660: PetscBLASInt *select;
661: PetscBLASInt *iwork;
662: PetscBLASInt liwork;
663: PetscScalar *Ht; /* Transpose of the Hessenberg matrix */
664: PetscScalar *t; /* Store the result of the solution of H^T*t=h_{m+1,m}e_m */
665: PetscBLASInt *ipiv; /* Permutation vector to be used in LAPACK */
666: PetscBool flag; /* determine whether to use Ritz vectors or harmonic Ritz vectors */
668: PetscFunctionBegin;
669: PetscCall(PetscBLASIntCast(n, &bn));
670: PetscCall(PetscBLASIntCast(N, &ldA));
671: ihi = ldQ = bn;
672: PetscCall(PetscBLASIntCast(5 * N, &lwork));
674: #if defined(PETSC_USE_COMPLEX)
675: SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "No support for complex numbers.");
676: #endif
678: PetscCall(PetscMalloc1(ldA * ldA, &A));
679: PetscCall(PetscMalloc1(ldQ * n, &Q));
680: PetscCall(PetscMalloc1(lwork, &work));
681: if (!dgmres->wr) {
682: PetscCall(PetscMalloc1(n, &dgmres->wr));
683: PetscCall(PetscMalloc1(n, &dgmres->wi));
684: }
685: wr = dgmres->wr;
686: wi = dgmres->wi;
687: PetscCall(PetscMalloc1(n, &modul));
688: PetscCall(PetscMalloc1(n, &perm));
689: /* copy the Hessenberg matrix to work space */
690: PetscCall(PetscArraycpy(A, dgmres->hes_origin, ldA * ldA));
691: PetscCall(PetscOptionsHasName(((PetscObject)ksp)->options, ((PetscObject)ksp)->prefix, "-ksp_dgmres_harmonic_ritz", &flag));
692: if (flag) {
693: /* Compute the matrix H + H^{-T}*h^2_{m+1,m}e_m*e_m^T */
694: /* Transpose the Hessenberg matrix */
695: PetscCall(PetscMalloc1(bn * bn, &Ht));
696: for (i = 0; i < bn; i++) {
697: for (j = 0; j < bn; j++) Ht[i * bn + j] = dgmres->hes_origin[j * ldA + i];
698: }
700: /* Solve the system H^T*t = h_{m+1,m}e_m */
701: PetscCall(PetscCalloc1(bn, &t));
702: t[bn - 1] = dgmres->hes_origin[(bn - 1) * ldA + bn]; /* Pick the last element H(m+1,m) */
703: PetscCall(PetscMalloc1(bn, &ipiv));
704: /* Call the LAPACK routine dgesv to solve the system Ht^-1 * t */
705: {
706: PetscBLASInt info;
707: PetscBLASInt nrhs = 1;
708: PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&bn, &nrhs, Ht, &bn, ipiv, t, &bn, &info));
709: PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error while calling the Lapack routine DGESV");
710: }
711: /* Now form H + H^{-T}*h^2_{m+1,m}e_m*e_m^T */
712: for (i = 0; i < bn; i++) A[(bn - 1) * bn + i] += t[i];
713: PetscCall(PetscFree(t));
714: PetscCall(PetscFree(Ht));
715: }
716: /* Compute eigenvalues with the Schur form */
717: {
718: PetscBLASInt info = 0;
719: PetscBLASInt ilo = 1;
720: PetscCallBLAS("LAPACKhseqr", LAPACKhseqr_("S", "I", &bn, &ilo, &ihi, A, &ldA, wr, wi, Q, &ldQ, work, &lwork, &info));
721: PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XHSEQR %d", (int)info);
722: }
723: PetscCall(PetscFree(work));
725: /* sort the eigenvalues */
726: for (i = 0; i < n; i++) modul[i] = PetscSqrtReal(wr[i] * wr[i] + wi[i] * wi[i]);
727: for (i = 0; i < n; i++) perm[i] = i;
729: PetscCall(PetscSortRealWithPermutation(n, modul, perm));
730: /* save the complex modulus of the largest eigenvalue in magnitude */
731: if (dgmres->lambdaN < modul[perm[n - 1]]) dgmres->lambdaN = modul[perm[n - 1]];
732: /* count the number of extracted eigenvalues (with complex conjugates) */
733: NbrEig = 0;
734: while (NbrEig < dgmres->neig) {
735: if (wi[perm[NbrEig]] != 0) NbrEig += 2;
736: else NbrEig += 1;
737: }
738: /* Reorder the Schur decomposition so that the cluster of smallest eigenvalues appears in the leading diagonal blocks of A */
740: PetscCall(PetscCalloc1(n, &select));
742: if (!dgmres->GreatestEig) {
743: for (j = 0; j < NbrEig; j++) select[perm[j]] = 1;
744: } else {
745: for (j = 0; j < NbrEig; j++) select[perm[n - j - 1]] = 1;
746: }
747: /* call Lapack dtrsen */
748: lwork = PetscMax(1, 4 * NbrEig * (bn - NbrEig));
749: liwork = PetscMax(1, 2 * NbrEig * (bn - NbrEig));
750: PetscCall(PetscMalloc1(lwork, &work));
751: PetscCall(PetscMalloc1(liwork, &iwork));
752: {
753: PetscBLASInt info = 0;
754: PetscReal CondEig; /* lower bound on the reciprocal condition number for the selected cluster of eigenvalues */
755: PetscReal CondSub; /* estimated reciprocal condition number of the specified invariant subspace. */
756: PetscCallBLAS("LAPACKtrsen", LAPACKtrsen_("B", "V", select, &bn, A, &ldA, Q, &ldQ, wr, wi, &NbrEig, &CondEig, &CondSub, work, &lwork, iwork, &liwork, &info));
757: PetscCheck(info != 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Unable to reorder the eigenvalues with the LAPACK routine: ILL-CONDITIONED PROBLEM");
758: }
759: PetscCall(PetscFree(select));
761: /* Extract the Schur vectors */
762: for (j = 0; j < NbrEig; j++) PetscCall(PetscArraycpy(&SR[j * N], &(Q[j * ldQ]), n));
763: *neig = NbrEig;
764: PetscCall(PetscFree(A));
765: PetscCall(PetscFree(work));
766: PetscCall(PetscFree(perm));
767: PetscCall(PetscFree(work));
768: PetscCall(PetscFree(iwork));
769: PetscCall(PetscFree(modul));
770: PetscCall(PetscFree(Q));
771: PetscFunctionReturn(PETSC_SUCCESS);
772: }
774: PetscErrorCode KSPDGMRESApplyDeflation_DGMRES(KSP ksp, Vec x, Vec y)
775: {
776: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
777: PetscInt i, r = dgmres->r;
778: PetscReal alpha = 1.0;
779: PetscInt max_neig = dgmres->max_neig;
780: PetscBLASInt br, bmax;
781: PetscReal lambda = dgmres->lambdaN;
783: PetscFunctionBegin;
784: PetscCall(PetscBLASIntCast(r, &br));
785: PetscCall(PetscBLASIntCast(max_neig, &bmax));
786: PetscCall(PetscLogEventBegin(KSP_DGMRESApplyDeflation, ksp, 0, 0, 0));
787: if (!r) {
788: PetscCall(VecCopy(x, y));
789: PetscFunctionReturn(PETSC_SUCCESS);
790: }
791: /* Compute U'*x */
792: if (!X1) {
793: PetscCall(PetscMalloc1(bmax, &X1));
794: PetscCall(PetscMalloc1(bmax, &X2));
795: }
796: PetscCall(VecMDot(x, r, UU, X1));
798: /* Solve T*X1=X2 for X1*/
799: PetscCall(PetscArraycpy(X2, X1, br));
800: {
801: PetscBLASInt info;
802: PetscBLASInt nrhs = 1;
803: PetscCallBLAS("LAPACKgetrs", LAPACKgetrs_("N", &br, &nrhs, TTF, &bmax, INVP, X1, &bmax, &info));
804: PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XGETRS %d", (int)info);
805: }
806: /* Iterative refinement -- is it really necessary ?? */
807: if (!WORK) {
808: PetscCall(PetscMalloc1(3 * bmax, &WORK));
809: PetscCall(PetscMalloc1(bmax, &IWORK));
810: }
811: {
812: PetscBLASInt info;
813: PetscReal berr, ferr;
814: PetscBLASInt nrhs = 1;
815: PetscCallBLAS("LAPACKgerfs", LAPACKgerfs_("N", &br, &nrhs, TT, &bmax, TTF, &bmax, INVP, X2, &bmax, X1, &bmax, &ferr, &berr, WORK, IWORK, &info));
816: PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XGERFS %d", (int)info);
817: }
819: for (i = 0; i < r; i++) X2[i] = X1[i] / lambda - X2[i];
821: /* Compute X2=U*X2 */
822: PetscCall(VecMAXPBY(y, r, X2, 0, UU));
823: PetscCall(VecAXPY(y, alpha, x));
825: PetscCall(PetscLogEventEnd(KSP_DGMRESApplyDeflation, ksp, 0, 0, 0));
826: PetscFunctionReturn(PETSC_SUCCESS);
827: }
829: static PetscErrorCode KSPDGMRESImproveEig_DGMRES(KSP ksp, PetscInt neig)
830: {
831: KSP_DGMRES *dgmres = (KSP_DGMRES *)ksp->data;
832: PetscInt j, r_old, r = dgmres->r;
833: PetscBLASInt i = 0;
834: PetscInt neig1 = dgmres->neig + EIG_OFFSET;
835: PetscInt bmax = dgmres->max_neig;
836: PetscInt aug = r + neig; /* actual size of the augmented invariant basis */
837: PetscInt aug1 = bmax + neig1; /* maximum size of the augmented invariant basis */
838: PetscBLASInt ldA; /* leading dimension of AUAU and AUU*/
839: PetscBLASInt N; /* size of AUAU */
840: PetscReal *Q; /* orthogonal matrix of (left) schur vectors */
841: PetscReal *Z; /* orthogonal matrix of (right) schur vectors */
842: PetscReal *work; /* working vector */
843: PetscBLASInt lwork; /* size of the working vector */
844: PetscInt *perm; /* Permutation vector to sort eigenvalues */
845: PetscReal *wr, *wi, *beta, *modul; /* Real and imaginary part and modulus of the eigenvalues of A*/
846: PetscBLASInt NbrEig = 0, nr, bm;
847: PetscBLASInt *select;
848: PetscBLASInt liwork, *iwork;
850: PetscFunctionBegin;
851: /* Block construction of the matrices AUU=(AU)'*U and (AU)'*AU*/
852: if (!AUU) {
853: PetscCall(PetscMalloc1(aug1 * aug1, &AUU));
854: PetscCall(PetscMalloc1(aug1 * aug1, &AUAU));
855: }
856: /* AUU = (AU)'*U = [(MU)'*U (MU)'*X; (MX)'*U (MX)'*X]
857: * Note that MU and MX have been computed previously either in ComputeDataDeflation() or down here in a previous call to this function */
858: /* (MU)'*U size (r x r) -- store in the <r> first columns of AUU*/
859: for (j = 0; j < r; j++) PetscCall(VecMDot(UU[j], r, MU, &AUU[j * aug1]));
860: /* (MU)'*X size (r x neig) -- store in AUU from the column <r>*/
861: for (j = 0; j < neig; j++) PetscCall(VecMDot(XX[j], r, MU, &AUU[(r + j) * aug1]));
862: /* (MX)'*U size (neig x r) -- store in the <r> first columns of AUU from the row <r>*/
863: for (j = 0; j < r; j++) PetscCall(VecMDot(UU[j], neig, MX, &AUU[j * aug1 + r]));
864: /* (MX)'*X size (neig neig) -- store in AUU from the column <r> and the row <r>*/
865: for (j = 0; j < neig; j++) PetscCall(VecMDot(XX[j], neig, MX, &AUU[(r + j) * aug1 + r]));
867: /* AUAU = (AU)'*AU = [(MU)'*MU (MU)'*MX; (MX)'*MU (MX)'*MX] */
868: /* (MU)'*MU size (r x r) -- store in the <r> first columns of AUAU*/
869: for (j = 0; j < r; j++) PetscCall(VecMDot(MU[j], r, MU, &AUAU[j * aug1]));
870: /* (MU)'*MX size (r x neig) -- store in AUAU from the column <r>*/
871: for (j = 0; j < neig; j++) PetscCall(VecMDot(MX[j], r, MU, &AUAU[(r + j) * aug1]));
872: /* (MX)'*MU size (neig x r) -- store in the <r> first columns of AUAU from the row <r>*/
873: for (j = 0; j < r; j++) PetscCall(VecMDot(MU[j], neig, MX, &AUAU[j * aug1 + r]));
874: /* (MX)'*MX size (neig neig) -- store in AUAU from the column <r> and the row <r>*/
875: for (j = 0; j < neig; j++) PetscCall(VecMDot(MX[j], neig, MX, &AUAU[(r + j) * aug1 + r]));
877: /* Computation of the eigenvectors */
878: PetscCall(PetscBLASIntCast(aug1, &ldA));
879: PetscCall(PetscBLASIntCast(aug, &N));
880: lwork = 8 * N + 20; /* sizeof the working space */
881: PetscCall(PetscMalloc1(N, &wr));
882: PetscCall(PetscMalloc1(N, &wi));
883: PetscCall(PetscMalloc1(N, &beta));
884: PetscCall(PetscMalloc1(N, &modul));
885: PetscCall(PetscMalloc1(N, &perm));
886: PetscCall(PetscMalloc1(N * N, &Q));
887: PetscCall(PetscMalloc1(N * N, &Z));
888: PetscCall(PetscMalloc1(lwork, &work));
889: {
890: PetscBLASInt info = 0;
891: PetscCallBLAS("LAPACKgges", LAPACKgges_("V", "V", "N", NULL, &N, AUAU, &ldA, AUU, &ldA, &i, wr, wi, beta, Q, &N, Z, &N, work, &lwork, NULL, &info));
892: PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XGGES %d", (int)info);
893: }
894: for (i = 0; i < N; i++) {
895: if (beta[i] != 0.0) {
896: wr[i] /= beta[i];
897: wi[i] /= beta[i];
898: }
899: }
900: /* sort the eigenvalues */
901: for (i = 0; i < N; i++) modul[i] = PetscSqrtReal(wr[i] * wr[i] + wi[i] * wi[i]);
902: for (i = 0; i < N; i++) perm[i] = i;
903: PetscCall(PetscSortRealWithPermutation(N, modul, perm));
904: /* Save the norm of the largest eigenvalue */
905: if (dgmres->lambdaN < modul[perm[N - 1]]) dgmres->lambdaN = modul[perm[N - 1]];
906: /* Allocate space to extract the first r schur vectors */
907: if (!SR2) PetscCall(PetscMalloc1(aug1 * bmax, &SR2));
908: /* count the number of extracted eigenvalues (complex conjugates count as 2) */
909: while (NbrEig < bmax) {
910: if (wi[perm[NbrEig]] == 0) NbrEig += 1;
911: else NbrEig += 2;
912: }
913: if (NbrEig > bmax) NbrEig = bmax - 1;
914: r_old = r; /* previous size of r */
915: dgmres->r = r = NbrEig;
917: /* Select the eigenvalues to reorder */
918: PetscCall(PetscCalloc1(N, &select));
919: if (!dgmres->GreatestEig) {
920: for (j = 0; j < NbrEig; j++) select[perm[j]] = 1;
921: } else {
922: for (j = 0; j < NbrEig; j++) select[perm[N - j - 1]] = 1;
923: }
924: /* Reorder and extract the new <r> schur vectors */
925: lwork = PetscMax(4 * N + 16, 2 * NbrEig * (N - NbrEig));
926: liwork = PetscMax(N + 6, 2 * NbrEig * (N - NbrEig));
927: PetscCall(PetscFree(work));
928: PetscCall(PetscMalloc1(lwork, &work));
929: PetscCall(PetscMalloc1(liwork, &iwork));
930: {
931: PetscBLASInt info = 0;
932: PetscReal Dif[2];
933: PetscBLASInt ijob = 2;
934: PetscBLASInt wantQ = 1, wantZ = 1;
935: 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));
936: PetscCheck(info != 1, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Unable to reorder the eigenvalues with the LAPACK routine: ill-conditioned problem.");
937: }
938: PetscCall(PetscFree(select));
940: for (j = 0; j < r; j++) PetscCall(PetscArraycpy(&SR2[j * aug1], &(Z[j * N]), N));
942: /* Multiply the Schur vectors SR2 by U (and X) to get a new U
943: -- save it temporarily in MU */
944: for (j = 0; j < r; j++) {
945: PetscCall(VecMAXPBY(MU[j], r_old, &SR2[j * aug1], 0, UU));
946: PetscCall(VecMAXPY(MU[j], neig, &SR2[j * aug1 + r_old], XX));
947: }
948: /* Form T = U'*MU*U */
949: for (j = 0; j < r; j++) {
950: PetscCall(VecCopy(MU[j], UU[j]));
951: PetscCall(KSP_PCApplyBAorAB(ksp, UU[j], MU[j], VEC_TEMP_MATOP));
952: }
953: dgmres->matvecs += r;
954: for (j = 0; j < r; j++) PetscCall(VecMDot(MU[j], r, UU, &TT[j * bmax]));
955: /* Factorize T */
956: PetscCall(PetscArraycpy(TTF, TT, bmax * r));
957: PetscCall(PetscBLASIntCast(r, &nr));
958: PetscCall(PetscBLASIntCast(bmax, &bm));
959: {
960: PetscBLASInt info;
961: PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&nr, &nr, TTF, &bm, INVP, &info));
962: PetscCheck(!info, PetscObjectComm((PetscObject)ksp), PETSC_ERR_LIB, "Error in LAPACK routine XGETRF INFO=%d", (int)info);
963: }
964: /* Free Memory */
965: PetscCall(PetscFree(wr));
966: PetscCall(PetscFree(wi));
967: PetscCall(PetscFree(beta));
968: PetscCall(PetscFree(modul));
969: PetscCall(PetscFree(perm));
970: PetscCall(PetscFree(Q));
971: PetscCall(PetscFree(Z));
972: PetscCall(PetscFree(work));
973: PetscCall(PetscFree(iwork));
974: PetscFunctionReturn(PETSC_SUCCESS);
975: }
977: /*MC
978: KSPDGMRES - Implements the deflated GMRES as defined in {cite}`erhel1996restarted` and {cite}`wakam2013memory`
980: Options Database Keys:
981: GMRES Options (inherited):
982: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
983: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
984: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially
985: (otherwise groups of vectors are allocated as needed)
986: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against
987: the Krylov space (fast) (the default)
988: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
989: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
990: stability of the classical Gram-Schmidt orthogonalization.
991: - -ksp_gmres_krylov_monitor - plot the Krylov space generated
993: DGMRES Options Database Keys:
994: + -ksp_dgmres_eigen <neig> - number of smallest eigenvalues to extract at each restart
995: . -ksp_dgmres_max_eigen <max_neig> - maximum number of eigenvalues that can be extracted during the iterative process
996: . -ksp_dgmres_force - use the deflation at each restart; switch off the adaptive strategy.
997: - -ksp_dgmres_view_deflation_vecs <viewerspec> - View the deflation vectors, where viewerspec is a key that can be
998: parsed by `PetscOptionsGetViewer()`. If neig > 1, viewerspec should
999: end with ":append". No vectors will be viewed if the adaptive
1000: strategy chooses not to deflate, so -ksp_dgmres_force should also
1001: be given.
1002: The deflation vectors span a subspace that may be a good
1003: approximation of the subspace of smallest eigenvectors of the
1004: preconditioned operator, so this option can aid in understanding
1005: the performance of a preconditioner.
1007: Level: beginner
1009: Notes:
1010: Left and right preconditioning are supported, but not symmetric preconditioning. Complex arithmetic is not supported
1012: In this implementation, the adaptive strategy allows switching to deflated GMRES when the stagnation occurs.
1014: Contributed by:
1015: Desire NUENTSA WAKAM, INRIA
1017: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPFGMRES`, `KSPLGMRES`,
1018: `KSPGMRESSetRestart()`, `KSPGMRESSetHapTol()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
1019: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`,
1020: `KSPGMRESCGSRefinementType`, `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESMonitorKrylov()`, `KSPSetPCSide()`
1021: M*/
1023: PETSC_EXTERN PetscErrorCode KSPCreate_DGMRES(KSP ksp)
1024: {
1025: KSP_DGMRES *dgmres;
1027: PetscFunctionBegin;
1028: PetscCall(PetscNew(&dgmres));
1029: ksp->data = (void *)dgmres;
1031: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3));
1032: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2));
1033: PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1));
1035: ksp->ops->buildsolution = KSPBuildSolution_DGMRES;
1036: ksp->ops->setup = KSPSetUp_DGMRES;
1037: ksp->ops->solve = KSPSolve_DGMRES;
1038: ksp->ops->destroy = KSPDestroy_DGMRES;
1039: ksp->ops->view = KSPView_DGMRES;
1040: ksp->ops->setfromoptions = KSPSetFromOptions_DGMRES;
1041: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
1042: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
1044: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetPreAllocateVectors_C", KSPGMRESSetPreAllocateVectors_GMRES));
1045: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetOrthogonalization_C", KSPGMRESSetOrthogonalization_GMRES));
1046: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetRestart_C", KSPGMRESSetRestart_GMRES));
1047: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetHapTol_C", KSPGMRESSetHapTol_GMRES));
1048: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPGMRESSetCGSRefinementType_C", KSPGMRESSetCGSRefinementType_GMRES));
1049: /* -- New functions defined in DGMRES -- */
1050: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetEigen_C", KSPDGMRESSetEigen_DGMRES));
1051: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetMaxEigen_C", KSPDGMRESSetMaxEigen_DGMRES));
1052: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESSetRatio_C", KSPDGMRESSetRatio_DGMRES));
1053: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESForce_C", KSPDGMRESForce_DGMRES));
1054: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESComputeSchurForm_C", KSPDGMRESComputeSchurForm_DGMRES));
1055: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESComputeDeflationData_C", KSPDGMRESComputeDeflationData_DGMRES));
1056: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESApplyDeflation_C", KSPDGMRESApplyDeflation_DGMRES));
1057: PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPDGMRESImproveEig_C", KSPDGMRESImproveEig_DGMRES));
1059: PetscCall(PetscLogEventRegister("DGMRESCompDefl", KSP_CLASSID, &KSP_DGMRESComputeDeflationData));
1060: PetscCall(PetscLogEventRegister("DGMRESApplyDefl", KSP_CLASSID, &KSP_DGMRESApplyDeflation));
1062: dgmres->haptol = 1.0e-30;
1063: dgmres->q_preallocate = 0;
1064: dgmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
1065: dgmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
1066: dgmres->nrs = NULL;
1067: dgmres->sol_temp = NULL;
1068: dgmres->max_k = GMRES_DEFAULT_MAXK;
1069: dgmres->Rsvd = NULL;
1070: dgmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
1071: dgmres->orthogwork = NULL;
1073: /* Default values for the deflation */
1074: dgmres->r = 0;
1075: dgmres->neig = DGMRES_DEFAULT_EIG;
1076: dgmres->max_neig = DGMRES_DEFAULT_MAXEIG - 1;
1077: dgmres->lambdaN = 0.0;
1078: dgmres->smv = SMV;
1079: dgmres->matvecs = 0;
1080: dgmres->GreatestEig = PETSC_FALSE; /* experimental */
1081: dgmres->HasSchur = PETSC_FALSE;
1082: PetscFunctionReturn(PETSC_SUCCESS);
1083: }