Actual source code: eige.c


  2: #include <petsc/private/kspimpl.h>
  3: #include <petscdm.h>
  4: #include <petscblaslapack.h>

  6: typedef struct {
  7:   KSP ksp;
  8:   Vec work;
  9: } Mat_KSP;

 11: static PetscErrorCode MatCreateVecs_KSP(Mat A,Vec *X,Vec *Y)
 12: {
 13:   Mat_KSP        *ctx;
 14:   Mat            M;

 18:   MatShellGetContext(A,&ctx);
 19:   KSPGetOperators(ctx->ksp,&M,NULL);
 20:   MatCreateVecs(M,X,Y);
 21:   return(0);
 22: }

 24: static PetscErrorCode MatMult_KSP(Mat A,Vec X,Vec Y)
 25: {
 26:   Mat_KSP        *ctx;

 30:   MatShellGetContext(A,&ctx);
 31:   KSP_PCApplyBAorAB(ctx->ksp,X,Y,ctx->work);
 32:   return(0);
 33: }

 35: /*@
 36:     KSPComputeOperator - Computes the explicit preconditioned operator, including diagonal scaling and null
 37:     space removal if applicable.

 39:     Collective on ksp

 41:     Input Parameters:
 42: +   ksp - the Krylov subspace context
 43: -   mattype - the matrix type to be used

 45:     Output Parameter:
 46: .   mat - the explicit preconditioned operator

 48:     Notes:
 49:     This computation is done by applying the operators to columns of the
 50:     identity matrix.

 52:     Currently, this routine uses a dense matrix format for the output operator if mattype == NULL.
 53:     This routine is costly in general, and is recommended for use only with relatively small systems.

 55:     Level: advanced

 57: .seealso: KSPComputeEigenvaluesExplicitly(), PCComputeOperator(), KSPSetDiagonalScale(), KSPSetNullSpace(), MatType
 58: @*/
 59: PetscErrorCode  KSPComputeOperator(KSP ksp, MatType mattype, Mat *mat)
 60: {
 62:   PetscInt       N,M,m,n;
 63:   Mat_KSP        ctx;
 64:   Mat            A,Aksp;

 69:   KSPGetOperators(ksp,&A,NULL);
 70:   MatGetLocalSize(A,&m,&n);
 71:   MatGetSize(A,&M,&N);
 72:   MatCreateShell(PetscObjectComm((PetscObject)ksp),m,n,M,N,&ctx,&Aksp);
 73:   MatShellSetOperation(Aksp,MATOP_MULT,(void (*)(void))MatMult_KSP);
 74:   MatShellSetOperation(Aksp,MATOP_CREATE_VECS,(void (*)(void))MatCreateVecs_KSP);
 75:   ctx.ksp = ksp;
 76:   MatCreateVecs(A,&ctx.work,NULL);
 77:   MatComputeOperator(Aksp,mattype,mat);
 78:   VecDestroy(&ctx.work);
 79:   MatDestroy(&Aksp);
 80:   return(0);
 81: }

 83: /*@
 84:    KSPComputeEigenvaluesExplicitly - Computes all of the eigenvalues of the
 85:    preconditioned operator using LAPACK.

 87:    Collective on ksp

 89:    Input Parameters:
 90: +  ksp - iterative context obtained from KSPCreate()
 91: -  n - size of arrays r and c

 93:    Output Parameters:
 94: +  r - real part of computed eigenvalues, provided by user with a dimension at least of n
 95: -  c - complex part of computed eigenvalues, provided by user with a dimension at least of n

 97:    Notes:
 98:    This approach is very slow but will generally provide accurate eigenvalue
 99:    estimates.  This routine explicitly forms a dense matrix representing
100:    the preconditioned operator, and thus will run only for relatively small
101:    problems, say n < 500.

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

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

111:    Level: advanced

113: .seealso: KSPComputeEigenvalues(), KSPMonitorSingularValue(), KSPComputeExtremeSingularValues(), KSPSetOperators(), KSPSolve()
114: @*/
115: PetscErrorCode  KSPComputeEigenvaluesExplicitly(KSP ksp,PetscInt nmax,PetscReal r[],PetscReal c[])
116: {
117:   Mat               BA;
118:   PetscErrorCode    ierr;
119:   PetscMPIInt       size,rank;
120:   MPI_Comm          comm;
121:   PetscScalar       *array;
122:   Mat               A;
123:   PetscInt          m,row,nz,i,n,dummy;
124:   const PetscInt    *cols;
125:   const PetscScalar *vals;

128:   PetscObjectGetComm((PetscObject)ksp,&comm);
129:   KSPComputeOperator(ksp,MATDENSE,&BA);
130:   MPI_Comm_size(comm,&size);
131:   MPI_Comm_rank(comm,&rank);

133:   MatGetSize(BA,&n,&n);
134:   if (size > 1) { /* assemble matrix on first processor */
135:     MatCreate(PetscObjectComm((PetscObject)ksp),&A);
136:     if (rank == 0) {
137:       MatSetSizes(A,n,n,n,n);
138:     } else {
139:       MatSetSizes(A,0,0,n,n);
140:     }
141:     MatSetType(A,MATMPIDENSE);
142:     MatMPIDenseSetPreallocation(A,NULL);
143:     PetscLogObjectParent((PetscObject)BA,(PetscObject)A);

145:     MatGetOwnershipRange(BA,&row,&dummy);
146:     MatGetLocalSize(BA,&m,&dummy);
147:     for (i=0; i<m; i++) {
148:       MatGetRow(BA,row,&nz,&cols,&vals);
149:       MatSetValues(A,1,&row,nz,cols,vals,INSERT_VALUES);
150:       MatRestoreRow(BA,row,&nz,&cols,&vals);
151:       row++;
152:     }

154:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
155:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
156:     MatDenseGetArray(A,&array);
157:   } else {
158:     MatDenseGetArray(BA,&array);
159:   }

161: #if !defined(PETSC_USE_COMPLEX)
162:   if (rank == 0) {
163:     PetscScalar  *work;
164:     PetscReal    *realpart,*imagpart;
165:     PetscBLASInt idummy,lwork;
166:     PetscInt     *perm;

168:     idummy   = n;
169:     lwork    = 5*n;
170:     PetscMalloc2(n,&realpart,n,&imagpart);
171:     PetscMalloc1(5*n,&work);
172:     {
173:       PetscBLASInt lierr;
174:       PetscScalar  sdummy;
175:       PetscBLASInt bn;

177:       PetscBLASIntCast(n,&bn);
178:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
179:       PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&bn,array,&bn,realpart,imagpart,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,&lierr));
180:       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in LAPACK routine %d",(int)lierr);
181:       PetscFPTrapPop();
182:     }
183:     PetscFree(work);
184:     PetscMalloc1(n,&perm);

186:     for (i=0; i<n; i++)  perm[i] = i;
187:     PetscSortRealWithPermutation(n,realpart,perm);
188:     for (i=0; i<n; i++) {
189:       r[i] = realpart[perm[i]];
190:       c[i] = imagpart[perm[i]];
191:     }
192:     PetscFree(perm);
193:     PetscFree2(realpart,imagpart);
194:   }
195: #else
196:   if (rank == 0) {
197:     PetscScalar  *work,*eigs;
198:     PetscReal    *rwork;
199:     PetscBLASInt idummy,lwork;
200:     PetscInt     *perm;

202:     idummy = n;
203:     lwork  = 5*n;
204:     PetscMalloc1(5*n,&work);
205:     PetscMalloc1(2*n,&rwork);
206:     PetscMalloc1(n,&eigs);
207:     {
208:       PetscBLASInt lierr;
209:       PetscScalar  sdummy;
210:       PetscBLASInt nb;
211:       PetscBLASIntCast(n,&nb);
212:       PetscFPTrapPush(PETSC_FP_TRAP_OFF);
213:       PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&nb,array,&nb,eigs,&sdummy,&idummy,&sdummy,&idummy,work,&lwork,rwork,&lierr));
214:       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in LAPACK routine %d",(int)lierr);
215:       PetscFPTrapPop();
216:     }
217:     PetscFree(work);
218:     PetscFree(rwork);
219:     PetscMalloc1(n,&perm);
220:     for (i=0; i<n; i++) perm[i] = i;
221:     for (i=0; i<n; i++) r[i]    = PetscRealPart(eigs[i]);
222:     PetscSortRealWithPermutation(n,r,perm);
223:     for (i=0; i<n; i++) {
224:       r[i] = PetscRealPart(eigs[perm[i]]);
225:       c[i] = PetscImaginaryPart(eigs[perm[i]]);
226:     }
227:     PetscFree(perm);
228:     PetscFree(eigs);
229:   }
230: #endif
231:   if (size > 1) {
232:     MatDenseRestoreArray(A,&array);
233:     MatDestroy(&A);
234:   } else {
235:     MatDenseRestoreArray(BA,&array);
236:   }
237:   MatDestroy(&BA);
238:   return(0);
239: }

241: static PetscErrorCode PolyEval(PetscInt nroots,const PetscReal *r,const PetscReal *c,PetscReal x,PetscReal y,PetscReal *px,PetscReal *py)
242: {
243:   PetscInt  i;
244:   PetscReal rprod = 1,iprod = 0;

247:   for (i=0; i<nroots; i++) {
248:     PetscReal rnew = rprod*(x - r[i]) - iprod*(y - c[i]);
249:     PetscReal inew = rprod*(y - c[i]) + iprod*(x - r[i]);
250:     rprod = rnew;
251:     iprod = inew;
252:   }
253:   *px = rprod;
254:   *py = iprod;
255:   return(0);
256: }

258: #include <petscdraw.h>
259: /* collective on ksp */
260: PetscErrorCode KSPPlotEigenContours_Private(KSP ksp,PetscInt neig,const PetscReal *r,const PetscReal *c)
261: {
263:   PetscReal      xmin,xmax,ymin,ymax,*xloc,*yloc,*value,px0,py0,rscale,iscale;
264:   PetscInt       M,N,i,j;
265:   PetscMPIInt    rank;
266:   PetscViewer    viewer;
267:   PetscDraw      draw;
268:   PetscDrawAxis  drawaxis;

271:   MPI_Comm_rank(PetscObjectComm((PetscObject)ksp),&rank);
272:   if (rank) return(0);
273:   M    = 80;
274:   N    = 80;
275:   xmin = r[0]; xmax = r[0];
276:   ymin = c[0]; ymax = c[0];
277:   for (i=1; i<neig; i++) {
278:     xmin = PetscMin(xmin,r[i]);
279:     xmax = PetscMax(xmax,r[i]);
280:     ymin = PetscMin(ymin,c[i]);
281:     ymax = PetscMax(ymax,c[i]);
282:   }
283:   PetscMalloc3(M,&xloc,N,&yloc,M*N,&value);
284:   for (i=0; i<M; i++) xloc[i] = xmin - 0.1*(xmax-xmin) + 1.2*(xmax-xmin)*i/(M-1);
285:   for (i=0; i<N; i++) yloc[i] = ymin - 0.1*(ymax-ymin) + 1.2*(ymax-ymin)*i/(N-1);
286:   PolyEval(neig,r,c,0,0,&px0,&py0);
287:   rscale = px0/(PetscSqr(px0)+PetscSqr(py0));
288:   iscale = -py0/(PetscSqr(px0)+PetscSqr(py0));
289:   for (j=0; j<N; j++) {
290:     for (i=0; i<M; i++) {
291:       PetscReal px,py,tx,ty,tmod;
292:       PolyEval(neig,r,c,xloc[i],yloc[j],&px,&py);
293:       tx   = px*rscale - py*iscale;
294:       ty   = py*rscale + px*iscale;
295:       tmod = PetscSqr(tx) + PetscSqr(ty); /* modulus of the complex polynomial */
296:       if (tmod > 1) tmod = 1.0;
297:       if (tmod > 0.5 && tmod < 1) tmod = 0.5;
298:       if (tmod > 0.2 && tmod < 0.5) tmod = 0.2;
299:       if (tmod > 0.05 && tmod < 0.2) tmod = 0.05;
300:       if (tmod < 1e-3) tmod = 1e-3;
301:       value[i+j*M] = PetscLogReal(tmod) / PetscLogReal(10.0);
302:     }
303:   }
304:   PetscViewerDrawOpen(PETSC_COMM_SELF,NULL,"Iteratively Computed Eigen-contours",PETSC_DECIDE,PETSC_DECIDE,450,450,&viewer);
305:   PetscViewerDrawGetDraw(viewer,0,&draw);
306:   PetscDrawTensorContour(draw,M,N,NULL,NULL,value);
307:   if (0) {
308:     PetscDrawAxisCreate(draw,&drawaxis);
309:     PetscDrawAxisSetLimits(drawaxis,xmin,xmax,ymin,ymax);
310:     PetscDrawAxisSetLabels(drawaxis,"Eigen-counters","real","imag");
311:     PetscDrawAxisDraw(drawaxis);
312:     PetscDrawAxisDestroy(&drawaxis);
313:   }
314:   PetscViewerDestroy(&viewer);
315:   PetscFree3(xloc,yloc,value);
316:   return(0);
317: }