Actual source code: pepdefault.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    Simple default routines for common PEP operations
 12: */

 14: #include <slepc/private/pepimpl.h>     /*I "slepcpep.h" I*/

 16: /*@
 17:    PEPSetWorkVecs - Sets a number of work vectors into a PEP object.

 19:    Collective on PEP

 21:    Input Parameters:
 22: +  pep - polynomial eigensolver context
 23: -  nw  - number of work vectors to allocate

 25:    Developers Note:
 26:    This is SLEPC_EXTERN because it may be required by user plugin PEP
 27:    implementations.

 29:    Level: developer
 30: @*/
 31: PetscErrorCode PEPSetWorkVecs(PEP pep,PetscInt nw)
 32: {
 34:   Vec            t;

 37:   if (pep->nwork < nw) {
 38:     VecDestroyVecs(pep->nwork,&pep->work);
 39:     pep->nwork = nw;
 40:     BVGetColumn(pep->V,0,&t);
 41:     VecDuplicateVecs(t,nw,&pep->work);
 42:     BVRestoreColumn(pep->V,0,&t);
 43:     PetscLogObjectParents(pep,nw,pep->work);
 44:   }
 45:   return(0);
 46: }

 48: /*
 49:   PEPConvergedRelative - Checks convergence relative to the eigenvalue.
 50: */
 51: PetscErrorCode PEPConvergedRelative(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 52: {
 53:   PetscReal w;

 56:   w = SlepcAbsEigenvalue(eigr,eigi);
 57:   *errest = res/w;
 58:   return(0);
 59: }

 61: /*
 62:   PEPConvergedNorm - Checks convergence relative to the matrix norms.
 63: */
 64: PetscErrorCode PEPConvergedNorm(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 65: {
 66:   PetscReal      w=0.0,t;
 67:   PetscInt       j;
 68:   PetscBool      flg;

 72:   /* initialization of matrix norms */
 73:   if (!pep->nrma[pep->nmat-1]) {
 74:     for (j=0;j<pep->nmat;j++) {
 75:       MatHasOperation(pep->A[j],MATOP_NORM,&flg);
 76:       if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"The convergence test related to the matrix norms requires a matrix norm operation");
 77:       MatNorm(pep->A[j],NORM_INFINITY,&pep->nrma[j]);
 78:     }
 79:   }
 80:   t = SlepcAbsEigenvalue(eigr,eigi);
 81:   for (j=pep->nmat-1;j>=0;j--) {
 82:     w = w*t+pep->nrma[j];
 83:   }
 84:   *errest = res/w;
 85:   return(0);
 86: }

 88: /*
 89:   PEPConvergedAbsolute - Checks convergence absolutely.
 90: */
 91: PetscErrorCode PEPConvergedAbsolute(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 92: {
 94:   *errest = res;
 95:   return(0);
 96: }

 98: /*@C
 99:    PEPStoppingBasic - Default routine to determine whether the outer eigensolver
100:    iteration must be stopped.

102:    Collective on PEP

104:    Input Parameters:
105: +  pep    - eigensolver context obtained from PEPCreate()
106: .  its    - current number of iterations
107: .  max_it - maximum number of iterations
108: .  nconv  - number of currently converged eigenpairs
109: .  nev    - number of requested eigenpairs
110: -  ctx    - context (not used here)

112:    Output Parameter:
113: .  reason - result of the stopping test

115:    Notes:
116:    A positive value of reason indicates that the iteration has finished successfully
117:    (converged), and a negative value indicates an error condition (diverged). If
118:    the iteration needs to be continued, reason must be set to PEP_CONVERGED_ITERATING
119:    (zero).

121:    PEPStoppingBasic() will stop if all requested eigenvalues are converged, or if
122:    the maximum number of iterations has been reached.

124:    Use PEPSetStoppingTest() to provide your own test instead of using this one.

126:    Level: advanced

128: .seealso: PEPSetStoppingTest(), PEPConvergedReason, PEPGetConvergedReason()
129: @*/
130: PetscErrorCode PEPStoppingBasic(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ctx)
131: {

135:   *reason = PEP_CONVERGED_ITERATING;
136:   if (nconv >= nev) {
137:     PetscInfo2(pep,"Polynomial eigensolver finished successfully: %D eigenpairs converged at iteration %D\n",nconv,its);
138:     *reason = PEP_CONVERGED_TOL;
139:   } else if (its >= max_it) {
140:     *reason = PEP_DIVERGED_ITS;
141:     PetscInfo1(pep,"Polynomial eigensolver iteration reached maximum number of iterations (%D)\n",its);
142:   }
143:   return(0);
144: }

146: PetscErrorCode PEPBackTransform_Default(PEP pep)
147: {

151:   STBackTransform(pep->st,pep->nconv,pep->eigr,pep->eigi);
152:   return(0);
153: }

155: PetscErrorCode PEPComputeVectors_Default(PEP pep)
156: {
158:   PetscInt       i;
159:   Vec            v;
160: #if !defined(PETSC_USE_COMPLEX)
161:   Vec            v1;
162: #endif

165:   PEPExtractVectors(pep);

167:   /* Fix eigenvectors if balancing was used */
168:   if ((pep->scale==PEP_SCALE_DIAGONAL || pep->scale==PEP_SCALE_BOTH) && pep->Dr && (pep->refine!=PEP_REFINE_MULTIPLE)) {
169:     for (i=0;i<pep->nconv;i++) {
170:       BVGetColumn(pep->V,i,&v);
171:       VecPointwiseMult(v,v,pep->Dr);
172:       BVRestoreColumn(pep->V,i,&v);
173:     }
174:   }

176:   /* normalization */
177:   for (i=0;i<pep->nconv;i++) {
178: #if !defined(PETSC_USE_COMPLEX)
179:     if (pep->eigi[i]!=0.0) {   /* first eigenvalue of a complex conjugate pair */
180:       BVGetColumn(pep->V,i,&v);
181:       BVGetColumn(pep->V,i+1,&v1);
182:       VecNormalizeComplex(v,v1,PETSC_TRUE,NULL);
183:       BVRestoreColumn(pep->V,i,&v);
184:       BVRestoreColumn(pep->V,i+1,&v1);
185:       i++;
186:     } else   /* real eigenvalue */
187: #endif
188:     {
189:       BVGetColumn(pep->V,i,&v);
190:       VecNormalizeComplex(v,NULL,PETSC_FALSE,NULL);
191:       BVRestoreColumn(pep->V,i,&v);
192:     }
193:   }
194:   return(0);
195: }

197: /*
198:   PEPBuildDiagonalScaling - compute two diagonal matrices to be applied for balancing
199:   in polynomial eigenproblems.
200: */
201: PetscErrorCode PEPBuildDiagonalScaling(PEP pep)
202: {
204:   PetscInt       it,i,j,k,nmat,nr,e,nz,lst,lend,nc=0,*cols,emax,emin,emaxl,eminl;
205:   const PetscInt *cidx,*ridx;
206:   Mat            M,*T,A;
207:   PetscMPIInt    n;
208:   PetscBool      cont=PETSC_TRUE,flg=PETSC_FALSE;
209:   PetscScalar    *array,*Dr,*Dl,t;
210:   PetscReal      l2,d,*rsum,*aux,*csum,w=1.0;
211:   MatStructure   str;
212:   MatInfo        info;

215:   l2 = 2*PetscLogReal(2.0);
216:   nmat = pep->nmat;
217:   PetscMPIIntCast(pep->n,&n);
218:   STGetMatStructure(pep->st,&str);
219:   PetscMalloc1(nmat,&T);
220:   for (k=0;k<nmat;k++) {
221:     STGetMatrixTransformed(pep->st,k,&T[k]);
222:   }
223:   /* Form local auxiliar matrix M */
224:   PetscObjectTypeCompareAny((PetscObject)T[0],&cont,MATMPIAIJ,MATSEQAIJ,"");
225:   if (!cont) SETERRQ(PetscObjectComm((PetscObject)T[0]),PETSC_ERR_SUP,"Only for MPIAIJ or SEQAIJ matrix types");
226:   PetscObjectTypeCompare((PetscObject)T[0],MATMPIAIJ,&cont);
227:   if (cont) {
228:     MatMPIAIJGetLocalMat(T[0],MAT_INITIAL_MATRIX,&M);
229:     flg = PETSC_TRUE;
230:   } else {
231:     MatDuplicate(T[0],MAT_COPY_VALUES,&M);
232:   }
233:   MatGetInfo(M,MAT_LOCAL,&info);
234:   nz = (PetscInt)info.nz_used;
235:   MatSeqAIJGetArray(M,&array);
236:   for (i=0;i<nz;i++) {
237:     t = PetscAbsScalar(array[i]);
238:     array[i] = t*t;
239:   }
240:   MatSeqAIJRestoreArray(M,&array);
241:   for (k=1;k<nmat;k++) {
242:     if (flg) {
243:       MatMPIAIJGetLocalMat(T[k],MAT_INITIAL_MATRIX,&A);
244:     } else {
245:       if (str==SAME_NONZERO_PATTERN) {
246:         MatCopy(T[k],A,SAME_NONZERO_PATTERN);
247:       } else {
248:         MatDuplicate(T[k],MAT_COPY_VALUES,&A);
249:       }
250:     }
251:     MatGetInfo(A,MAT_LOCAL,&info);
252:     nz = (PetscInt)info.nz_used;
253:     MatSeqAIJGetArray(A,&array);
254:     for (i=0;i<nz;i++) {
255:       t = PetscAbsScalar(array[i]);
256:       array[i] = t*t;
257:     }
258:     MatSeqAIJRestoreArray(A,&array);
259:     w *= pep->slambda*pep->slambda*pep->sfactor;
260:     MatAXPY(M,w,A,str);
261:     if (flg || str!=SAME_NONZERO_PATTERN || k==nmat-2) {
262:       MatDestroy(&A);
263:     }
264:   }
265:   MatGetRowIJ(M,0,PETSC_FALSE,PETSC_FALSE,&nr,&ridx,&cidx,&cont);
266:   if (!cont) SETERRQ(PetscObjectComm((PetscObject)T[0]),PETSC_ERR_SUP,"It is not possible to compute scaling diagonals for these PEP matrices");
267:   MatGetInfo(M,MAT_LOCAL,&info);
268:   nz = (PetscInt)info.nz_used;
269:   VecGetOwnershipRange(pep->Dl,&lst,&lend);
270:   PetscMalloc4(nr,&rsum,pep->n,&csum,pep->n,&aux,PetscMin(pep->n-lend+lst,nz),&cols);
271:   VecSet(pep->Dr,1.0);
272:   VecSet(pep->Dl,1.0);
273:   VecGetArray(pep->Dl,&Dl);
274:   VecGetArray(pep->Dr,&Dr);
275:   MatSeqAIJGetArray(M,&array);
276:   PetscMemzero(aux,pep->n*sizeof(PetscReal));
277:   for (j=0;j<nz;j++) {
278:     /* Search non-zero columns outsize lst-lend */
279:     if (aux[cidx[j]]==0 && (cidx[j]<lst || lend<=cidx[j])) cols[nc++] = cidx[j];
280:     /* Local column sums */
281:     aux[cidx[j]] += PetscAbsScalar(array[j]);
282:   }
283:   for (it=0;it<pep->sits && cont;it++) {
284:     emaxl = 0; eminl = 0;
285:     /* Column sum  */
286:     if (it>0) { /* it=0 has been already done*/
287:       MatSeqAIJGetArray(M,&array);
288:       PetscMemzero(aux,pep->n*sizeof(PetscReal));
289:       for (j=0;j<nz;j++) aux[cidx[j]] += PetscAbsScalar(array[j]);
290:       MatSeqAIJRestoreArray(M,&array);
291:     }
292:     MPI_Allreduce(aux,csum,n,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)pep->Dr));
293:     /* Update Dr */
294:     for (j=lst;j<lend;j++) {
295:       d = PetscLogReal(csum[j])/l2;
296:       e = -(PetscInt)((d < 0)?(d-0.5):(d+0.5));
297:       d = PetscPowReal(2.0,e);
298:       Dr[j-lst] *= d;
299:       aux[j] = d*d;
300:       emaxl = PetscMax(emaxl,e);
301:       eminl = PetscMin(eminl,e);
302:     }
303:     for (j=0;j<nc;j++) {
304:       d = PetscLogReal(csum[cols[j]])/l2;
305:       e = -(PetscInt)((d < 0)?(d-0.5):(d+0.5));
306:       d = PetscPowReal(2.0,e);
307:       aux[cols[j]] = d*d;
308:       emaxl = PetscMax(emaxl,e);
309:       eminl = PetscMin(eminl,e);
310:     }
311:     /* Scale M */
312:     MatSeqAIJGetArray(M,&array);
313:     for (j=0;j<nz;j++) {
314:       array[j] *= aux[cidx[j]];
315:     }
316:     MatSeqAIJRestoreArray(M,&array);
317:     /* Row sum */
318:     PetscMemzero(rsum,nr*sizeof(PetscReal));
319:     MatSeqAIJGetArray(M,&array);
320:     for (i=0;i<nr;i++) {
321:       for (j=ridx[i];j<ridx[i+1];j++) rsum[i] += PetscAbsScalar(array[j]);
322:       /* Update Dl */
323:       d = PetscLogReal(rsum[i])/l2;
324:       e = -(PetscInt)((d < 0)?(d-0.5):(d+0.5));
325:       d = PetscPowReal(2.0,e);
326:       Dl[i] *= d;
327:       /* Scale M */
328:       for (j=ridx[i];j<ridx[i+1];j++) array[j] *= d*d;
329:       emaxl = PetscMax(emaxl,e);
330:       eminl = PetscMin(eminl,e);
331:     }
332:     MatSeqAIJRestoreArray(M,&array);
333:     /* Compute global max and min */
334:     MPI_Allreduce(&emaxl,&emax,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)pep->Dl));
335:     MPI_Allreduce(&eminl,&emin,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)pep->Dl));
336:     if (emax<=emin+2) cont = PETSC_FALSE;
337:   }
338:   VecRestoreArray(pep->Dr,&Dr);
339:   VecRestoreArray(pep->Dl,&Dl);
340:   /* Free memory*/
341:   MatDestroy(&M);
342:   PetscFree4(rsum,csum,aux,cols);
343:   PetscFree(T);
344:   return(0);
345: }

347: /*
348:    PEPComputeScaleFactor - compute sfactor as described in [Betcke 2008].
349: */
350: PetscErrorCode PEPComputeScaleFactor(PEP pep)
351: {
353:   PetscBool      has0,has1,flg;
354:   PetscReal      norm0,norm1;
355:   Mat            T[2];
356:   PEPBasis       basis;
357:   PetscInt       i;

360:   if (pep->scale==PEP_SCALE_NONE || pep->scale==PEP_SCALE_DIAGONAL) {  /* no scalar scaling */
361:     pep->sfactor = 1.0;
362:     pep->dsfactor = 1.0;
363:     return(0);
364:   }
365:   if (pep->sfactor_set) return(0);  /* user provided value */
366:   pep->sfactor = 1.0;
367:   pep->dsfactor = 1.0;
368:   PEPGetBasis(pep,&basis);
369:   if (basis==PEP_BASIS_MONOMIAL) {
370:     STGetTransform(pep->st,&flg);
371:     if (flg) {
372:       STGetMatrixTransformed(pep->st,0,&T[0]);
373:       STGetMatrixTransformed(pep->st,pep->nmat-1,&T[1]);
374:     } else {
375:       T[0] = pep->A[0];
376:       T[1] = pep->A[pep->nmat-1];
377:     }
378:     if (pep->nmat>2) {
379:       MatHasOperation(T[0],MATOP_NORM,&has0);
380:       MatHasOperation(T[1],MATOP_NORM,&has1);
381:       if (has0 && has1) {
382:         MatNorm(T[0],NORM_INFINITY,&norm0);
383:         MatNorm(T[1],NORM_INFINITY,&norm1);
384:         pep->sfactor = PetscPowReal(norm0/norm1,1.0/(pep->nmat-1));
385:         pep->dsfactor = norm1;
386:         for (i=pep->nmat-2;i>0;i--) {
387:           STGetMatrixTransformed(pep->st,i,&T[1]);
388:           MatHasOperation(T[1],MATOP_NORM,&has1);
389:           if (has1) {
390:             MatNorm(T[1],NORM_INFINITY,&norm1);
391:             pep->dsfactor = pep->dsfactor*pep->sfactor+norm1;
392:           } else break;
393:         }
394:         if (has1) {
395:           pep->dsfactor = pep->dsfactor*pep->sfactor+norm0;
396:           pep->dsfactor = pep->nmat/pep->dsfactor;
397:         } else pep->dsfactor = 1.0;
398:       }
399:     }
400:   }
401:   return(0);
402: }

404: /*
405:    PEPBasisCoefficients - compute polynomial basis coefficients
406: */
407: PetscErrorCode PEPBasisCoefficients(PEP pep,PetscReal *pbc)
408: {
409:   PetscReal *ca,*cb,*cg;
410:   PetscInt  k,nmat=pep->nmat;

413:   ca = pbc;
414:   cb = pbc+nmat;
415:   cg = pbc+2*nmat;
416:   switch (pep->basis) {
417:   case PEP_BASIS_MONOMIAL:
418:     for (k=0;k<nmat;k++) {
419:       ca[k] = 1.0; cb[k] = 0.0; cg[k] = 0.0;
420:     }
421:     break;
422:   case PEP_BASIS_CHEBYSHEV1:
423:     ca[0] = 1.0; cb[0] = 0.0; cg[0] = 0.0;
424:     for (k=1;k<nmat;k++) {
425:       ca[k] = .5; cb[k] = 0.0; cg[k] = .5;
426:     }
427:     break;
428:   case PEP_BASIS_CHEBYSHEV2:
429:     ca[0] = .5; cb[0] = 0.0; cg[0] = 0.0;
430:     for (k=1;k<nmat;k++) {
431:       ca[k] = .5; cb[k] = 0.0; cg[k] = .5;
432:     }
433:     break;
434:   case PEP_BASIS_LEGENDRE:
435:     ca[0] = 1.0; cb[0] = 0.0; cg[0] = 0.0;
436:     for (k=1;k<nmat;k++) {
437:       ca[k] = k+1; cb[k] = -2*k; cg[k] = k;
438:     }
439:     break;
440:   case PEP_BASIS_LAGUERRE:
441:     ca[0] = -1.0; cb[0] = 0.0; cg[0] = 0.0;
442:     for (k=1;k<nmat;k++) {
443:       ca[k] = -(k+1); cb[k] = 2*k+1; cg[k] = -k;
444:     }
445:     break;
446:   case PEP_BASIS_HERMITE:
447:     ca[0] = .5; cb[0] = 0.0; cg[0] = 0.0;
448:     for (k=1;k<nmat;k++) {
449:       ca[k] = .5; cb[k] = 0.0; cg[k] = -k;
450:     }
451:     break;
452:   }
453:   return(0);
454: }