Actual source code: eige.c

  1: #define PETSCKSP_DLL

 3:  #include private/kspimpl.h
 4:  #include petscblaslapack.h

  8: /*@
  9:     KSPComputeExplicitOperator - Computes the explicit preconditioned operator.  

 11:     Collective on KSP

 13:     Input Parameter:
 14: .   ksp - the Krylov subspace context

 16:     Output Parameter:
 17: .   mat - the explict preconditioned operator

 19:     Notes:
 20:     This computation is done by applying the operators to columns of the 
 21:     identity matrix.

 23:     Currently, this routine uses a dense matrix format when 1 processor
 24:     is used and a sparse format otherwise.  This routine is costly in general,
 25:     and is recommended for use only with relatively small systems.

 27:     Level: advanced
 28:    
 29: .keywords: KSP, compute, explicit, operator

 31: .seealso: KSPComputeEigenvaluesExplicitly(), PCComputeExplicitOperator()
 32: @*/
 33: PetscErrorCode  KSPComputeExplicitOperator(KSP ksp,Mat *mat)
 34: {
 35:   Vec            in,out;
 37:   PetscMPIInt    size;
 38:   PetscInt       i,M,m,*rows,start,end;
 39:   Mat            A;
 40:   MPI_Comm       comm;
 41:   PetscScalar    *array,one = 1.0;

 46:   comm = ((PetscObject)ksp)->comm;

 48:   MPI_Comm_size(comm,&size);

 50:   VecDuplicate(ksp->vec_sol,&in);
 51:   VecDuplicate(ksp->vec_sol,&out);
 52:   VecGetSize(in,&M);
 53:   VecGetLocalSize(in,&m);
 54:   VecGetOwnershipRange(in,&start,&end);
 55:   PetscMalloc(m*sizeof(PetscInt),&rows);
 56:   for (i=0; i<m; i++) {rows[i] = start + i;}

 58:   MatCreate(comm,mat);
 59:   MatSetSizes(*mat,m,m,M,M);
 60:   if (size == 1) {
 61:     MatSetType(*mat,MATSEQDENSE);
 62:     MatSeqDenseSetPreallocation(*mat,PETSC_NULL);
 63:   } else {
 64:     MatSetType(*mat,MATMPIAIJ);
 65:     MatMPIAIJSetPreallocation(*mat,0,PETSC_NULL,0,PETSC_NULL);
 66:   }
 67: 
 68:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
 69:   PCGetOperators(ksp->pc,&A,PETSC_NULL,PETSC_NULL);

 71:   for (i=0; i<M; i++) {

 73:     VecSet(in,0.0);
 74:     VecSetValues(in,1,&i,&one,INSERT_VALUES);
 75:     VecAssemblyBegin(in);
 76:     VecAssemblyEnd(in);

 78:     KSP_MatMult(ksp,A,in,out);
 79:     KSP_PCApply(ksp,out,in);
 80: 
 81:     VecGetArray(in,&array);
 82:     MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
 83:     VecRestoreArray(in,&array);

 85:   }
 86:   PetscFree(rows);
 87:   VecDestroy(in);
 88:   VecDestroy(out);
 89:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
 90:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
 91:   return(0);
 92: }

 96: /*@
 97:    KSPComputeEigenvaluesExplicitly - Computes all of the eigenvalues of the 
 98:    preconditioned operator using LAPACK.  

100:    Collective on KSP

102:    Input Parameter:
103: +  ksp - iterative context obtained from KSPCreate()
104: -  n - size of arrays r and c

106:    Output Parameters:
107: +  r - real part of computed eigenvalues
108: -  c - complex part of computed eigenvalues

110:    Notes:
111:    This approach is very slow but will generally provide accurate eigenvalue
112:    estimates.  This routine explicitly forms a dense matrix representing 
113:    the preconditioned operator, and thus will run only for relatively small
114:    problems, say n < 500.

116:    Many users may just want to use the monitoring routine
117:    KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
118:    to print the singular values at each iteration of the linear solve.

120:    The preconditoner operator, rhs vector, solution vectors should be
121:    set before this routine is called. i.e use KSPSetOperators(),KSPSolve() or
122:    KSPSetOperators()

124:    Level: advanced

126: .keywords: KSP, compute, eigenvalues, explicitly

128: .seealso: KSPComputeEigenvalues(), KSPMonitorSingularValue(), KSPComputeExtremeSingularValues(), KSPSetOperators(), KSPSolve()
129: @*/
130: PetscErrorCode  KSPComputeEigenvaluesExplicitly(KSP ksp,PetscInt nmax,PetscReal *r,PetscReal *c)
131: {
132:   Mat                BA;
133:   PetscErrorCode     ierr;
134:   PetscMPIInt        size,rank;
135:   MPI_Comm           comm = ((PetscObject)ksp)->comm;
136:   PetscScalar        *array;
137:   Mat                A;
138:   PetscInt           m,row,nz,i,n,dummy;
139:   const PetscInt     *cols;
140:   const PetscScalar  *vals;

143:   KSPComputeExplicitOperator(ksp,&BA);
144:   MPI_Comm_size(comm,&size);
145:   MPI_Comm_rank(comm,&rank);

147:   MatGetSize(BA,&n,&n);
148:   if (size > 1) { /* assemble matrix on first processor */
149:     MatCreate(((PetscObject)ksp)->comm,&A);
150:     if (!rank) {
151:       MatSetSizes(A,n,n,n,n);
152:     } else {
153:       MatSetSizes(A,0,0,n,n);
154:     }
155:     MatSetType(A,MATMPIDENSE);
156:     MatMPIDenseSetPreallocation(A,PETSC_NULL);
157:     PetscLogObjectParent(BA,A);

159:     MatGetOwnershipRange(BA,&row,&dummy);
160:     MatGetLocalSize(BA,&m,&dummy);
161:     for (i=0; i<m; i++) {
162:       MatGetRow(BA,row,&nz,&cols,&vals);
163:       MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);
164:       MatRestoreRow(BA,row,&nz,&cols,&vals);
165:       row++;
166:     }

168:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
169:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
170:     MatGetArray(A,&array);
171:   } else {
172:     MatGetArray(BA,&array);
173:   }

175: #if defined(PETSC_HAVE_ESSL)
176:   /* ESSL has a different calling sequence for dgeev() and zgeev() than standard LAPACK */
177:   if (!rank) {
178:     PetscScalar  sdummy,*cwork;
179:     PetscReal    *work,*realpart;
180:     PetscBLASInt clen,idummy,lwork,*perm,zero = 0;

182: #if !defined(PETSC_USE_COMPLEX)
183:     clen = n;
184: #else
185:     clen = 2*n;
186: #endif
187:     PetscMalloc(clen*sizeof(PetscScalar),&cwork);
188:     idummy = n;
189:     lwork  = 5*n;
190:     PetscMalloc(lwork*sizeof(PetscReal),&work);
191:     PetscMalloc(n*sizeof(PetscReal),&realpart);
192:     LAPACKgeev_(&zero,array,&n,cwork,&sdummy,&idummy,&idummy,&n,work,&lwork);
193:     PetscFree(work);

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

199: #if !defined(PETSC_USE_COMPLEX)
200:     for (i=0; i<n; i++) {
201:       realpart[i] = cwork[2*i];
202:       perm[i]     = i;
203:     }
204:     PetscSortRealWithPermutation(n,realpart,perm);
205:     for (i=0; i<n; i++) {
206:       r[i] = cwork[2*perm[i]];
207:       c[i] = cwork[2*perm[i]+1];
208:     }
209: #else
210:     for (i=0; i<n; i++) {
211:       realpart[i] = PetscRealPart(cwork[i]);
212:       perm[i]     = i;
213:     }
214:     PetscSortRealWithPermutation(n,realpart,perm);
215:     for (i=0; i<n; i++) {
216:       r[i] = PetscRealPart(cwork[perm[i]]);
217:       c[i] = PetscImaginaryPart(cwork[perm[i]]);
218:     }
219: #endif
220:     PetscFree(perm);
221:     PetscFree(realpart);
222:     PetscFree(cwork);
223:   }
224: #elif !defined(PETSC_USE_COMPLEX)
225:   if (!rank) {
226:     PetscScalar  *work;
227:     PetscReal    *realpart,*imagpart;
228:     PetscBLASInt idummy,lwork;
229:     PetscInt     *perm;

231:     idummy   = n;
232:     lwork    = 5*n;
233:     PetscMalloc(2*n*sizeof(PetscReal),&realpart);
234:     imagpart = realpart + n;
235:     PetscMalloc(5*n*sizeof(PetscReal),&work);
236: #if defined(PETSC_MISSING_LAPACK_GEEV) 
237:     SETERRQ(PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
238: #else
239:     {
240:       PetscBLASInt lierr;
241:       PetscScalar sdummy;
242:       PetscBLASInt bn = PetscBLASIntCast(n);
243:       LAPACKgeev_("N","N",&bn,array,&bn,realpart,imagpart,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,&lierr);
244:       if (lierr) SETERRQ1(PETSC_ERR_LIB,"Error in LAPACK routine %d",(int)lierr);
245:     }
246: #endif
247:     PetscFree(work);
248:     PetscMalloc(n*sizeof(PetscInt),&perm);
249:     for (i=0; i<n; i++) { perm[i] = i;}
250:     PetscSortRealWithPermutation(n,realpart,perm);
251:     for (i=0; i<n; i++) {
252:       r[i] = realpart[perm[i]];
253:       c[i] = imagpart[perm[i]];
254:     }
255:     PetscFree(perm);
256:     PetscFree(realpart);
257:   }
258: #else
259:   if (!rank) {
260:     PetscScalar  *work,*eigs;
261:     PetscReal    *rwork;
262:     PetscBLASInt idummy,lwork;
263:     PetscInt     *perm;

265:     idummy   = n;
266:     lwork    = 5*n;
267:     PetscMalloc(5*n*sizeof(PetscScalar),&work);
268:     PetscMalloc(2*n*sizeof(PetscReal),&rwork);
269:     PetscMalloc(n*sizeof(PetscScalar),&eigs);
270: #if defined(PETSC_MISSING_LAPACK_GEEV) 
271:     SETERRQ(PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
272: #else
273:     {
274:       PetscBLASInt lierr;
275:       PetscScalar  sdummy;
276:       PetscBLASInt nb = PetscBLASIntCast(n);
277:       LAPACKgeev_("N","N",&nb,array,&nb,eigs,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,rwork,&lierr);
278:       if (lierr) SETERRQ1(PETSC_ERR_LIB,"Error in LAPACK routine %d",(int)lierr);
279:     }
280: #endif
281:     PetscFree(work);
282:     PetscFree(rwork);
283:     PetscMalloc(n*sizeof(PetscInt),&perm);
284:     for (i=0; i<n; i++) { perm[i] = i;}
285:     for (i=0; i<n; i++) { r[i]    = PetscRealPart(eigs[i]);}
286:     PetscSortRealWithPermutation(n,r,perm);
287:     for (i=0; i<n; i++) {
288:       r[i] = PetscRealPart(eigs[perm[i]]);
289:       c[i] = PetscImaginaryPart(eigs[perm[i]]);
290:     }
291:     PetscFree(perm);
292:     PetscFree(eigs);
293:   }
294: #endif  
295:   if (size > 1) {
296:     MatRestoreArray(A,&array);
297:     MatDestroy(A);
298:   } else {
299:     MatRestoreArray(BA,&array);
300:   }
301:   MatDestroy(BA);
302:   return(0);
303: }