Actual source code: gmreig.c

petsc-3.3-p5 2012-12-01
  2: #include <../src/ksp/ksp/impls/gmres/gmresimpl.h>
  3: #include <petscblaslapack.h>

  7: PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP ksp,PetscReal *emax,PetscReal *emin)
  8: {
  9: #if defined(PETSC_MISSING_LAPACK_GESVD) 
 11:   /*
 12:       The Cray math libraries on T3D/T3E, and early versions of Intel Math Kernel Libraries (MKL)
 13:       for PCs do not seem to have the DGESVD() lapack routines
 14:   */
 15:   SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable\nNot able to provide singular value estimates.");
 16: #else
 17:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
 19:   PetscInt       n = gmres->it + 1,i,N = gmres->max_k + 2;
 20:   PetscBLASInt   bn, bN ,lwork, idummy,lierr;
 21:   PetscScalar    *R = gmres->Rsvd,*work = R + N*N,sdummy;
 22:   PetscReal      *realpart = gmres->Dsvd;

 25:   bn = PetscBLASIntCast(n);
 26:   bN = PetscBLASIntCast(N);
 27:   lwork = PetscBLASIntCast(5*N);
 28:   idummy = PetscBLASIntCast(N);
 29:   if (n <= 0) {
 30:     *emax = *emin = 1.0;
 31:     return(0);
 32:   }
 33:   /* copy R matrix to work space */
 34:   PetscMemcpy(R,gmres->hh_origin,(gmres->max_k+2)*(gmres->max_k+1)*sizeof(PetscScalar));

 36:   /* zero below diagonal garbage */
 37:   for (i=0; i<n; i++) {
 38:     R[i*N+i+1] = 0.0;
 39:   }
 40: 
 41:   /* compute Singular Values */
 42:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
 43: #if !defined(PETSC_USE_COMPLEX)
 44:   LAPACKgesvd_("N","N",&bn,&bn,R,&bN,realpart,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,&lierr);
 45: #else
 46:   LAPACKgesvd_("N","N",&bn,&bn,R,&bN,realpart,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,realpart+N,&lierr);
 47: #endif
 48:   if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SVD Lapack routine %d",(int)lierr);
 49:   PetscFPTrapPop();

 51:   *emin = realpart[n-1];
 52:   *emax = realpart[0];
 53: #endif
 54:   return(0);
 55: }

 57: /* ------------------------------------------------------------------------ */
 58: /* ESSL has a different calling sequence for dgeev() and zgeev() than standard LAPACK */
 61: PetscErrorCode KSPComputeEigenvalues_GMRES(KSP ksp,PetscInt nmax,PetscReal *r,PetscReal *c,PetscInt *neig)
 62: {
 63: #if defined(PETSC_HAVE_ESSL)
 64:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
 66:   PetscInt       n = gmres->it + 1,N = gmres->max_k + 1;
 67:   PetscInt       i,*perm;
 68:   PetscScalar    *R = gmres->Rsvd;
 69:   PetscScalar    *cwork = R + N*N,sdummy;
 70:   PetscReal      *work,*realpart = gmres->Dsvd ;
 71:   PetscBLASInt   zero = 0,bn,bN,idummy,lwork;

 74:   bn = PetscBLASIntCast(n);
 75:   bN = PetscBLASIntCast(N);
 76:   idummy = -1;                  /* unused */
 77:   lwork = PetscBLASIntCast(5*N);
 78:   if (nmax < n) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_SIZ,"Not enough room in work space r and c for eigenvalues");
 79:   *neig = n;

 81:   if (!n) {
 82:     return(0);
 83:   }
 84:   /* copy R matrix to work space */
 85:   PetscMemcpy(R,gmres->hes_origin,N*N*sizeof(PetscScalar));

 87:   /* compute eigenvalues */

 89:   /* for ESSL version need really cwork of length N (complex), 2N
 90:      (real); already at least 5N of space has been allocated */

 92:   PetscMalloc(lwork*sizeof(PetscReal),&work);
 93:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
 94:   LAPACKgeev_(&zero,R,&bN,cwork,&sdummy,&idummy,&idummy,&bn,work,&lwork);
 95:   PetscFPTrapPop();
 96:   PetscFree(work);

 98:   /* For now we stick with the convention of storing the real and imaginary
 99:      components of evalues separately.  But is this what we really want? */
100:   PetscMalloc(n*sizeof(PetscInt),&perm);

102: #if !defined(PETSC_USE_COMPLEX)
103:   for (i=0; i<n; i++) {
104:     realpart[i] = cwork[2*i];
105:     perm[i]     = i;
106:   }
107:   PetscSortRealWithPermutation(n,realpart,perm);
108:   for (i=0; i<n; i++) {
109:     r[i] = cwork[2*perm[i]];
110:     c[i] = cwork[2*perm[i]+1];
111:   }
112: #else
113:   for (i=0; i<n; i++) {
114:     realpart[i] = PetscRealPart(cwork[i]);
115:     perm[i]     = i;
116:   }
117:   PetscSortRealWithPermutation(n,realpart,perm);
118:   for (i=0; i<n; i++) {
119:     r[i] = PetscRealPart(cwork[perm[i]]);
120:     c[i] = PetscImaginaryPart(cwork[perm[i]]);
121:   }
122: #endif
123:   PetscFree(perm);
124: #elif defined(PETSC_MISSING_LAPACK_GEEV) 
126:   SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
127: #elif !defined(PETSC_USE_COMPLEX)
128:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
130:   PetscInt       n = gmres->it + 1,N = gmres->max_k + 1,i,*perm;
131:   PetscBLASInt   bn, bN, lwork, idummy, lierr;
132:   PetscScalar    *R = gmres->Rsvd,*work = R + N*N;
133:   PetscScalar    *realpart = gmres->Dsvd,*imagpart = realpart + N,sdummy;

136:   bn = PetscBLASIntCast(n);
137:   bN = PetscBLASIntCast(N);
138:   lwork = PetscBLASIntCast(5*N);
139:   idummy = PetscBLASIntCast(N);
140:   if (nmax < n) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_SIZ,"Not enough room in work space r and c for eigenvalues");
141:   *neig = n;

143:   if (!n) {
144:     return(0);
145:   }

147:   /* copy R matrix to work space */
148:   PetscMemcpy(R,gmres->hes_origin,N*N*sizeof(PetscScalar));

150:   /* compute eigenvalues */
151:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
152:   LAPACKgeev_("N","N",&bn,R,&bN,realpart,imagpart,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,&lierr);
153:   if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in LAPACK routine %d",(int)lierr);
154:   PetscFPTrapPop();
155:   PetscMalloc(n*sizeof(PetscInt),&perm);
156:   for (i=0; i<n; i++) { perm[i] = i;}
157:   PetscSortRealWithPermutation(n,realpart,perm);
158:   for (i=0; i<n; i++) {
159:     r[i] = realpart[perm[i]];
160:     c[i] = imagpart[perm[i]];
161:   }
162:   PetscFree(perm);
163: #else
164:   KSP_GMRES      *gmres = (KSP_GMRES*)ksp->data;
166:   PetscInt       n = gmres->it + 1,N = gmres->max_k + 1,i,*perm;
167:   PetscScalar    *R = gmres->Rsvd,*work = R + N*N,*eigs = work + 5*N,sdummy;
168:   PetscBLASInt   bn,bN,lwork,idummy,lierr;

171:   bn = PetscBLASIntCast(n);
172:   bN = PetscBLASIntCast(N);
173:   lwork = PetscBLASIntCast(5*N);
174:   idummy = PetscBLASIntCast(N);
175:   if (nmax < n) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_SIZ,"Not enough room in work space r and c for eigenvalues");
176:   *neig = n;

178:   if (!n) {
179:     return(0);
180:   }
181:   /* copy R matrix to work space */
182:   PetscMemcpy(R,gmres->hes_origin,N*N*sizeof(PetscScalar));

184:   /* compute eigenvalues */
185:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
186:   LAPACKgeev_("N","N",&bn,R,&bN,eigs,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,gmres->Dsvd,&lierr);
187:   if (lierr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in LAPACK routine");
188:   PetscFPTrapPop();
189:   PetscMalloc(n*sizeof(PetscInt),&perm);
190:   for (i=0; i<n; i++) { perm[i] = i;}
191:   for (i=0; i<n; i++) { r[i]    = PetscRealPart(eigs[i]);}
192:   PetscSortRealWithPermutation(n,r,perm);
193:   for (i=0; i<n; i++) {
194:     r[i] = PetscRealPart(eigs[perm[i]]);
195:     c[i] = PetscImaginaryPart(eigs[perm[i]]);
196:   }
197:   PetscFree(perm);
198: #endif
199:   return(0);
200: }