Actual source code: agmresleja.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: #define PETSCKSP_DLL
  2: /*
  3:  * Functions in this file reorder the Ritz values in the (modified) Leja order.
  4:  *
  5:  * References : [1] Bai, Zhaojun and  Hu, D. and Reichel, L. A Newton basis GMRES implementation. IMA J. Numer. Anal. 14 (1994), no. 4, 563-581.
  6:  *
  7:  */
  8:  #include <../src/ksp/ksp/impls/gmres/agmres/agmresimpl.h>

 10: static PetscErrorCode KSPAGMRESLejafmaxarray(PetscScalar *re, PetscInt pt, PetscInt n,PetscInt *pos)
 11: {
 12:   PetscInt    i;
 13:   PetscScalar mx;

 16:   mx   = re[0];
 17:   *pos = 0;
 18:   for (i = pt; i < n; i++) {
 19:     if (mx < re[i]) {
 20:       mx   = re[i];
 21:       *pos = i;
 22:     }
 23:   }
 24:   return(0);
 25: }

 27: static PetscErrorCode KSPAGMRESLejaCfpdMax(PetscScalar *rm, PetscScalar *im, PetscInt *spos, PetscInt nbre, PetscInt n, PetscInt *rpos)
 28: {
 29:   PetscScalar rd, id, pd, max;
 30:   PetscInt    i, j;

 33:   pd    = 1.0;
 34:   max   = 0.0;
 35:   *rpos = 0;
 36:   for (i = 0; i < n; i++) {
 37:     for (j = 0; j < nbre; j++) {
 38:       rd = rm[i] - rm[spos[j]];
 39:       id = im[i] - im[spos[j]];
 40:       pd = pd * PetscSqrtReal(rd*rd + id*id);
 41:     }
 42:     if (max < pd) {
 43:       *rpos = i;
 44:       max   = pd;
 45:     }
 46:     pd = 1.0;
 47:   }
 48:   return(0);
 49: }

 51: PetscErrorCode KSPAGMRESLejaOrdering(PetscScalar *re, PetscScalar *im, PetscScalar *rre, PetscScalar *rim, PetscInt m)
 52: {
 53:   PetscInt       *spos;
 54:   PetscScalar    *n_cmpl,temp;
 56:   PetscInt       i, pos, j;

 59:   PetscMalloc1(m, &n_cmpl);
 60:   PetscMalloc1(m, &spos);
 61:   /* Check the proper order of complex conjugate pairs */
 62:   j = 0;
 63:   while (j  < m) {
 64:     if (im[j] != 0.0) { /* complex eigenvalue */
 65:       if (im[j] < 0.0) { /* change the order */
 66:         temp    = im[j+1];
 67:         im[j+1] = im[j];
 68:         im[j]   = temp;
 69:       }
 70:       j += 2;
 71:     } else j++;
 72:   }

 74:   for (i = 0; i < m; i++) n_cmpl[i] = PetscSqrtReal(re[i]*re[i]+im[i]*im[i]);
 75:   KSPAGMRESLejafmaxarray(n_cmpl, 0, m, &pos);
 76:   j = 0;
 77:   if (im[pos] >= 0.0) {
 78:     rre[0] = re[pos];
 79:     rim[0] = im[pos];
 80:     j++;
 81:     spos[0] = pos;
 82:   }
 83:   while (j < (m)) {
 84:     if (im[pos] > 0) {
 85:       rre[j]  = re[pos+1];
 86:       rim[j]  = im[pos+1];
 87:       spos[j] = pos + 1;
 88:       j++;
 89:     }
 90:     KSPAGMRESLejaCfpdMax(re, im, spos, j, m, &pos);
 91:     if (im[pos] < 0) pos--;

 93:     if ((im[pos] >= 0) && (j < m)) {
 94:       rre[j]  = re[pos];
 95:       rim[j]  = im[pos];
 96:       spos[j] = pos;
 97:       j++;
 98:     }
 99:   }
100:   PetscFree(spos);
101:   PetscFree(n_cmpl);
102:   return(0);
103: }