Actual source code: gmreig.c

  1: /*$Id: gmreig.c,v 1.28 2001/08/07 03:03:51 balay Exp $*/

 3:  #include src/sles/ksp/impls/gmres/gmresp.h
 4:  #include petscblaslapack.h

  6: int KSPComputeExtremeSingularValues_GMRES(KSP ksp,PetscReal *emax,PetscReal *emin)
  7: {
  8:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
  9:   int       n = gmres->it + 1,N = gmres->max_k + 2,ierr,lwork = 5*N,idummy = N,i;
 10:   PetscScalar    *R = gmres->Rsvd,*work = R + N*N,sdummy;
 11:   PetscReal *realpart = gmres->Dsvd;

 14:   if (!n) {
 15:     *emax = *emin = 1.0;
 16:     return(0);
 17:   }
 18:   /* copy R matrix to work space */
 19:   PetscMemcpy(R,gmres->hh_origin,N*N*sizeof(PetscScalar));

 21:   /* zero below diagonal garbage */
 22:   for (i=0; i<n; i++) {
 23:     R[i*N+i+1] = 0.0;
 24:   }
 25: 
 26:   /* compute Singular Values */
 27:   /*
 28:       The Cray math libraries on T3D/T3E, and Intel Math Kernel Libraries (MKL) for PCs do not 
 29:       seem to have the DGESVD() lapack routines
 30:   */
 31: #if defined(PETSC_MISSING_LAPACK_GESVD) 
 32:   SETERRQ(PETSC_ERR_SUP,"GESVD - Lapack routine is unavilablenNot able to provide singular value estimates.");
 33: #else
 34: #if !defined(PETSC_USE_COMPLEX)
 35:   LAgesvd_("N","N",&n,&n,R,&N,realpart,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,&ierr);
 36: #else
 37:   LAgesvd_("N","N",&n,&n,R,&N,realpart,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,realpart+N,&ierr);
 38: #endif
 39:   if (ierr) SETERRQ(PETSC_ERR_LIB,"Error in SVD Lapack routine");

 41:   *emin = realpart[n-1];
 42:   *emax = realpart[0];

 44:   return(0);
 45: #endif
 46: }
 47: /* ------------------------------------------------------------------------ */
 48: /* ESSL has a different calling sequence for dgeev() and zgeev() than standard LAPACK */
 49: #if defined(PETSC_HAVE_ESSL)
 50: int KSPComputeEigenvalues_GMRES(KSP ksp,int nmax,PetscReal *r,PetscReal *c,int *neig)
 51: {
 52:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
 53:   int       n = gmres->it + 1,N = gmres->max_k + 1,ierr,lwork = 5*N;
 54:   int       idummy = N,i,*perm,clen,zero;
 55:   PetscScalar    *R = gmres->Rsvd;
 56:   PetscScalar    *cwork = R + N*N,sdummy;
 57:   PetscReal *work,*realpart = gmres->Dsvd,*imagpart = realpart + N ;

 60:   if (nmax < n) SETERRQ(PETSC_ERR_ARG_SIZ,"Not enough room in work space r and c for eigenvalues");
 61:   *neig = n;

 63:   if (!n) {
 64:     return(0);
 65:   }
 66:   /* copy R matrix to work space */
 67:   PetscMemcpy(R,gmres->hes_origin,N*N*sizeof(PetscScalar));

 69:   /* compute eigenvalues */

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

 74:   PetscMalloc(lwork*sizeof(PetscReal),&work);
 75:   zero = 0;
 76:   LAgeev_(&zero,R,&N,cwork,&sdummy,&idummy,&idummy,&n,work,&lwork);
 77:   PetscFree(work);

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

 83: #if !defined(PETSC_USE_COMPLEX)
 84:   for (i=0; i<n; i++) {
 85:     realpart[i] = cwork[2*i];
 86:     perm[i]     = i;
 87:   }
 88:   PetscSortRealWithPermutation(n,realpart,perm);
 89:   for (i=0; i<n; i++) {
 90:     r[i] = cwork[2*perm[i]];
 91:     c[i] = cwork[2*perm[i]+1];
 92:   }
 93: #else
 94:   for (i=0; i<n; i++) {
 95:     realpart[i] = PetscRealPart(cwork[i]);
 96:     perm[i]     = i;
 97:   }
 98:   PetscSortRealWithPermutation(n,realpart,perm);
 99:   for (i=0; i<n; i++) {
100:     r[i] = PetscRealPart(cwork[perm[i]]);
101:     c[i] = PetscImaginaryPart(cwork[perm[i]]);
102:   }
103: #endif
104:   PetscFree(perm);
105:   return(0);
106: }
107: #elif !defined(PETSC_USE_COMPLEX)
108: int KSPComputeEigenvalues_GMRES(KSP ksp,int nmax,PetscReal *r,PetscReal *c,int *neig)
109: {
110:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
111:   int       n = gmres->it + 1,N = gmres->max_k + 1,ierr,lwork = 5*N,idummy = N,i,*perm;
112:   PetscScalar    *R = gmres->Rsvd,*work = R + N*N;
113:   PetscScalar    *realpart = gmres->Dsvd,*imagpart = realpart + N,sdummy;

116:   if (nmax < n) SETERRQ(PETSC_ERR_ARG_SIZ,"Not enough room in work space r and c for eigenvalues");
117:   *neig = n;

119:   if (!n) {
120:     return(0);
121:   }

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

126:   /* compute eigenvalues */
127: #if defined(PETSC_MISSING_LAPACK_GEEV) 
128:   SETERRQ(PETSC_ERR_SUP,"GEEV - Lapack routine is unavilablenNot able to provide eigen values.");
129: #else
130:   LAgeev_("N","N",&n,R,&N,realpart,imagpart,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,&ierr);
131: #endif
132:   if (ierr) SETERRQ(PETSC_ERR_LIB,"Error in LAPACK routine");
133:   PetscMalloc(n*sizeof(int),&perm);
134:   for (i=0; i<n; i++) { perm[i] = i;}
135:   PetscSortRealWithPermutation(n,realpart,perm);
136:   for (i=0; i<n; i++) {
137:     r[i] = realpart[perm[i]];
138:     c[i] = imagpart[perm[i]];
139:   }
140:   PetscFree(perm);
141:   return(0);
142: }
143: #else
144: int KSPComputeEigenvalues_GMRES(KSP ksp,int nmax,PetscReal *r,PetscReal *c,int *neig)
145: {
146:   KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
147:   int       n = gmres->it + 1,N = gmres->max_k + 1,ierr,lwork = 5*N,idummy = N,i,*perm;
148:   PetscScalar    *R = gmres->Rsvd,*work = R + N*N,*eigs = work + 5*N,sdummy;

151:   if (nmax < n) SETERRQ(PETSC_ERR_ARG_SIZ,"Not enough room in work space r and c for eigenvalues");
152:   *neig = n;

154:   if (!n) {
155:     return(0);
156:   }
157:   /* copy R matrix to work space */
158:   PetscMemcpy(R,gmres->hes_origin,N*N*sizeof(PetscScalar));

160:   /* compute eigenvalues */
161: #if defined(PETSC_MISSING_LAPACK_GEEV) 
162:   SETERRQ(PETSC_ERR_SUP,"GEEV - Lapack routine is unavilablenNot able to provide eigen values.");
163: #else
164:   LAgeev_("N","N",&n,R,&N,eigs,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,gmres->Dsvd,&ierr);
165: #endif
166:   if (ierr) SETERRQ(PETSC_ERR_LIB,"Error in LAPACK routine");
167:   PetscMalloc(n*sizeof(int),&perm);
168:   for (i=0; i<n; i++) { perm[i] = i;}
169:   for (i=0; i<n; i++) { r[i]    = PetscRealPart(eigs[i]);}
170:   PetscSortRealWithPermutation(n,r,perm);
171:   for (i=0; i<n; i++) {
172:     r[i] = PetscRealPart(eigs[perm[i]]);
173:     c[i] = PetscImaginaryPart(eigs[perm[i]]);
174:   }
175:   PetscFree(perm);
176:   return(0);
177: }
178: #endif