Actual source code: pepdefault.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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>

 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:    Developer Notes:
 26:    This is SLEPC_EXTERN because it may be required by user plugin PEP
 27:    implementations.

 29:    Level: developer

 31: .seealso: PEPSetUp()
 32: @*/
 33: PetscErrorCode PEPSetWorkVecs(PEP pep,PetscInt nw)
 34: {
 35:   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:   PetscFunctionReturn(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;

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

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

 69:   /* initialization of matrix norms */
 70:   if (!pep->nrma[pep->nmat-1]) {
 71:     for (j=0;j<pep->nmat;j++) {
 72:       MatHasOperation(pep->A[j],MATOP_NORM,&flg);
 74:       MatNorm(pep->A[j],NORM_INFINITY,&pep->nrma[j]);
 75:     }
 76:   }
 77:   t = SlepcAbsEigenvalue(eigr,eigi);
 78:   for (j=pep->nmat-1;j>=0;j--) {
 79:     w = w*t+pep->nrma[j];
 80:   }
 81:   *errest = res/w;
 82:   PetscFunctionReturn(0);
 83: }

 85: /*
 86:   PEPSetWhichEigenpairs_Default - Sets the default value for which,
 87:   depending on the ST.
 88:  */
 89: PetscErrorCode PEPSetWhichEigenpairs_Default(PEP pep)
 90: {
 91:   PetscBool      target;

 93:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&target);
 94:   if (target) pep->which = PEP_TARGET_MAGNITUDE;
 95:   else pep->which = PEP_LARGEST_MAGNITUDE;
 96:   PetscFunctionReturn(0);
 97: }

 99: /*
100:   PEPConvergedAbsolute - Checks convergence absolutely.
101: */
102: PetscErrorCode PEPConvergedAbsolute(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
103: {
104:   *errest = res;
105:   PetscFunctionReturn(0);
106: }

108: /*@C
109:    PEPStoppingBasic - Default routine to determine whether the outer eigensolver
110:    iteration must be stopped.

112:    Collective on pep

114:    Input Parameters:
115: +  pep    - eigensolver context obtained from PEPCreate()
116: .  its    - current number of iterations
117: .  max_it - maximum number of iterations
118: .  nconv  - number of currently converged eigenpairs
119: .  nev    - number of requested eigenpairs
120: -  ctx    - context (not used here)

122:    Output Parameter:
123: .  reason - result of the stopping test

125:    Notes:
126:    A positive value of reason indicates that the iteration has finished successfully
127:    (converged), and a negative value indicates an error condition (diverged). If
128:    the iteration needs to be continued, reason must be set to PEP_CONVERGED_ITERATING
129:    (zero).

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

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

136:    Level: advanced

138: .seealso: PEPSetStoppingTest(), PEPConvergedReason, PEPGetConvergedReason()
139: @*/
140: PetscErrorCode PEPStoppingBasic(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ctx)
141: {
142:   *reason = PEP_CONVERGED_ITERATING;
143:   if (nconv >= nev) {
144:     PetscInfo(pep,"Polynomial eigensolver finished successfully: %" PetscInt_FMT " eigenpairs converged at iteration %" PetscInt_FMT "\n",nconv,its);
145:     *reason = PEP_CONVERGED_TOL;
146:   } else if (its >= max_it) {
147:     *reason = PEP_DIVERGED_ITS;
148:     PetscInfo(pep,"Polynomial eigensolver iteration reached maximum number of iterations (%" PetscInt_FMT ")\n",its);
149:   }
150:   PetscFunctionReturn(0);
151: }

153: PetscErrorCode PEPBackTransform_Default(PEP pep)
154: {
155:   STBackTransform(pep->st,pep->nconv,pep->eigr,pep->eigi);
156:   PetscFunctionReturn(0);
157: }

159: PetscErrorCode PEPComputeVectors_Default(PEP pep)
160: {
161:   PetscInt       i;
162:   Vec            v;

164:   PEPExtractVectors(pep);

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

175:   /* normalization */
176:   BVNormalize(pep->V,pep->eigi);
177:   PetscFunctionReturn(0);
178: }

180: /*
181:   PEPBuildDiagonalScaling - compute two diagonal matrices to be applied for balancing
182:   in polynomial eigenproblems.
183: */
184: PetscErrorCode PEPBuildDiagonalScaling(PEP pep)
185: {
186:   PetscInt       it,i,j,k,nmat,nr,e,nz,lst,lend,nc=0,*cols,emax,emin,emaxl,eminl;
187:   const PetscInt *cidx,*ridx;
188:   Mat            M,*T,A;
189:   PetscMPIInt    n;
190:   PetscBool      cont=PETSC_TRUE,flg=PETSC_FALSE;
191:   PetscScalar    *array,*Dr,*Dl,t;
192:   PetscReal      l2,d,*rsum,*aux,*csum,w=1.0;
193:   MatStructure   str;
194:   MatInfo        info;

196:   l2 = 2*PetscLogReal(2.0);
197:   nmat = pep->nmat;
198:   PetscMPIIntCast(pep->n,&n);
199:   STGetMatStructure(pep->st,&str);
200:   PetscMalloc1(nmat,&T);
201:   for (k=0;k<nmat;k++) STGetMatrixTransformed(pep->st,k,&T[k]);
202:   /* Form local auxiliary matrix M */
203:   PetscObjectBaseTypeCompareAny((PetscObject)T[0],&cont,MATMPIAIJ,MATSEQAIJ,"");
205:   PetscObjectBaseTypeCompare((PetscObject)T[0],MATMPIAIJ,&cont);
206:   if (cont) {
207:     MatMPIAIJGetLocalMat(T[0],MAT_INITIAL_MATRIX,&M);
208:     flg = PETSC_TRUE;
209:   } else MatDuplicate(T[0],MAT_COPY_VALUES,&M);
210:   MatGetInfo(M,MAT_LOCAL,&info);
211:   nz = (PetscInt)info.nz_used;
212:   MatSeqAIJGetArray(M,&array);
213:   for (i=0;i<nz;i++) {
214:     t = PetscAbsScalar(array[i]);
215:     array[i] = t*t;
216:   }
217:   MatSeqAIJRestoreArray(M,&array);
218:   for (k=1;k<nmat;k++) {
219:     if (flg) MatMPIAIJGetLocalMat(T[k],MAT_INITIAL_MATRIX,&A);
220:     else {
221:       if (str==SAME_NONZERO_PATTERN) MatCopy(T[k],A,SAME_NONZERO_PATTERN);
222:       else MatDuplicate(T[k],MAT_COPY_VALUES,&A);
223:     }
224:     MatGetInfo(A,MAT_LOCAL,&info);
225:     nz = (PetscInt)info.nz_used;
226:     MatSeqAIJGetArray(A,&array);
227:     for (i=0;i<nz;i++) {
228:       t = PetscAbsScalar(array[i]);
229:       array[i] = t*t;
230:     }
231:     MatSeqAIJRestoreArray(A,&array);
232:     w *= pep->slambda*pep->slambda*pep->sfactor;
233:     MatAXPY(M,w,A,str);
234:     if (flg || str!=SAME_NONZERO_PATTERN || k==nmat-2) MatDestroy(&A);
235:   }
236:   MatGetRowIJ(M,0,PETSC_FALSE,PETSC_FALSE,&nr,&ridx,&cidx,&cont);
238:   MatGetInfo(M,MAT_LOCAL,&info);
239:   nz = (PetscInt)info.nz_used;
240:   VecGetOwnershipRange(pep->Dl,&lst,&lend);
241:   PetscMalloc4(nr,&rsum,pep->n,&csum,pep->n,&aux,PetscMin(pep->n-lend+lst,nz),&cols);
242:   VecSet(pep->Dr,1.0);
243:   VecSet(pep->Dl,1.0);
244:   VecGetArray(pep->Dl,&Dl);
245:   VecGetArray(pep->Dr,&Dr);
246:   MatSeqAIJGetArray(M,&array);
247:   PetscArrayzero(aux,pep->n);
248:   for (j=0;j<nz;j++) {
249:     /* Search non-zero columns outsize lst-lend */
250:     if (aux[cidx[j]]==0 && (cidx[j]<lst || lend<=cidx[j])) cols[nc++] = cidx[j];
251:     /* Local column sums */
252:     aux[cidx[j]] += PetscAbsScalar(array[j]);
253:   }
254:   for (it=0;it<pep->sits && cont;it++) {
255:     emaxl = 0; eminl = 0;
256:     /* Column sum  */
257:     if (it>0) { /* it=0 has been already done*/
258:       MatSeqAIJGetArray(M,&array);
259:       PetscArrayzero(aux,pep->n);
260:       for (j=0;j<nz;j++) aux[cidx[j]] += PetscAbsScalar(array[j]);
261:       MatSeqAIJRestoreArray(M,&array);
262:     }
263:     MPIU_Allreduce(aux,csum,n,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)pep->Dr));
264:     /* Update Dr */
265:     for (j=lst;j<lend;j++) {
266:       d = PetscLogReal(csum[j])/l2;
267:       e = -(PetscInt)((d < 0)?(d-0.5):(d+0.5));
268:       d = PetscPowReal(2.0,e);
269:       Dr[j-lst] *= d;
270:       aux[j] = d*d;
271:       emaxl = PetscMax(emaxl,e);
272:       eminl = PetscMin(eminl,e);
273:     }
274:     for (j=0;j<nc;j++) {
275:       d = PetscLogReal(csum[cols[j]])/l2;
276:       e = -(PetscInt)((d < 0)?(d-0.5):(d+0.5));
277:       d = PetscPowReal(2.0,e);
278:       aux[cols[j]] = d*d;
279:       emaxl = PetscMax(emaxl,e);
280:       eminl = PetscMin(eminl,e);
281:     }
282:     /* Scale M */
283:     MatSeqAIJGetArray(M,&array);
284:     for (j=0;j<nz;j++) {
285:       array[j] *= aux[cidx[j]];
286:     }
287:     MatSeqAIJRestoreArray(M,&array);
288:     /* Row sum */
289:     PetscArrayzero(rsum,nr);
290:     MatSeqAIJGetArray(M,&array);
291:     for (i=0;i<nr;i++) {
292:       for (j=ridx[i];j<ridx[i+1];j++) rsum[i] += PetscAbsScalar(array[j]);
293:       /* Update Dl */
294:       d = PetscLogReal(rsum[i])/l2;
295:       e = -(PetscInt)((d < 0)?(d-0.5):(d+0.5));
296:       d = PetscPowReal(2.0,e);
297:       Dl[i] *= d;
298:       /* Scale M */
299:       for (j=ridx[i];j<ridx[i+1];j++) array[j] *= d*d;
300:       emaxl = PetscMax(emaxl,e);
301:       eminl = PetscMin(eminl,e);
302:     }
303:     MatSeqAIJRestoreArray(M,&array);
304:     /* Compute global max and min */
305:     MPIU_Allreduce(&emaxl,&emax,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)pep->Dl));
306:     MPIU_Allreduce(&eminl,&emin,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)pep->Dl));
307:     if (emax<=emin+2) cont = PETSC_FALSE;
308:   }
309:   VecRestoreArray(pep->Dr,&Dr);
310:   VecRestoreArray(pep->Dl,&Dl);
311:   /* Free memory*/
312:   MatDestroy(&M);
313:   PetscFree4(rsum,csum,aux,cols);
314:   PetscFree(T);
315:   PetscFunctionReturn(0);
316: }

318: /*
319:    PEPComputeScaleFactor - compute sfactor as described in [Betcke 2008].
320: */
321: PetscErrorCode PEPComputeScaleFactor(PEP pep)
322: {
323:   PetscBool      has0,has1,flg;
324:   PetscReal      norm0,norm1;
325:   Mat            T[2];
326:   PEPBasis       basis;
327:   PetscInt       i;

329:   if (pep->scale==PEP_SCALE_NONE || pep->scale==PEP_SCALE_DIAGONAL) {  /* no scalar scaling */
330:     pep->sfactor = 1.0;
331:     pep->dsfactor = 1.0;
332:     PetscFunctionReturn(0);
333:   }
334:   if (pep->sfactor_set) PetscFunctionReturn(0);  /* user provided value */
335:   pep->sfactor = 1.0;
336:   pep->dsfactor = 1.0;
337:   PEPGetBasis(pep,&basis);
338:   if (basis==PEP_BASIS_MONOMIAL) {
339:     STGetTransform(pep->st,&flg);
340:     if (flg) {
341:       STGetMatrixTransformed(pep->st,0,&T[0]);
342:       STGetMatrixTransformed(pep->st,pep->nmat-1,&T[1]);
343:     } else {
344:       T[0] = pep->A[0];
345:       T[1] = pep->A[pep->nmat-1];
346:     }
347:     if (pep->nmat>2) {
348:       MatHasOperation(T[0],MATOP_NORM,&has0);
349:       MatHasOperation(T[1],MATOP_NORM,&has1);
350:       if (has0 && has1) {
351:         MatNorm(T[0],NORM_INFINITY,&norm0);
352:         MatNorm(T[1],NORM_INFINITY,&norm1);
353:         pep->sfactor = PetscPowReal(norm0/norm1,1.0/(pep->nmat-1));
354:         pep->dsfactor = norm1;
355:         for (i=pep->nmat-2;i>0;i--) {
356:           STGetMatrixTransformed(pep->st,i,&T[1]);
357:           MatHasOperation(T[1],MATOP_NORM,&has1);
358:           if (has1) {
359:             MatNorm(T[1],NORM_INFINITY,&norm1);
360:             pep->dsfactor = pep->dsfactor*pep->sfactor+norm1;
361:           } else break;
362:         }
363:         if (has1) {
364:           pep->dsfactor = pep->dsfactor*pep->sfactor+norm0;
365:           pep->dsfactor = pep->nmat/pep->dsfactor;
366:         } else pep->dsfactor = 1.0;
367:       }
368:     }
369:   }
370:   PetscFunctionReturn(0);
371: }

373: /*
374:    PEPBasisCoefficients - compute polynomial basis coefficients
375: */
376: PetscErrorCode PEPBasisCoefficients(PEP pep,PetscReal *pbc)
377: {
378:   PetscReal *ca,*cb,*cg;
379:   PetscInt  k,nmat=pep->nmat;

381:   ca = pbc;
382:   cb = pbc+nmat;
383:   cg = pbc+2*nmat;
384:   switch (pep->basis) {
385:   case PEP_BASIS_MONOMIAL:
386:     for (k=0;k<nmat;k++) {
387:       ca[k] = 1.0; cb[k] = 0.0; cg[k] = 0.0;
388:     }
389:     break;
390:   case PEP_BASIS_CHEBYSHEV1:
391:     ca[0] = 1.0; cb[0] = 0.0; cg[0] = 0.0;
392:     for (k=1;k<nmat;k++) {
393:       ca[k] = .5; cb[k] = 0.0; cg[k] = .5;
394:     }
395:     break;
396:   case PEP_BASIS_CHEBYSHEV2:
397:     ca[0] = .5; cb[0] = 0.0; cg[0] = 0.0;
398:     for (k=1;k<nmat;k++) {
399:       ca[k] = .5; cb[k] = 0.0; cg[k] = .5;
400:     }
401:     break;
402:   case PEP_BASIS_LEGENDRE:
403:     ca[0] = 1.0; cb[0] = 0.0; cg[0] = 0.0;
404:     for (k=1;k<nmat;k++) {
405:       ca[k] = k+1; cb[k] = -2*k; cg[k] = k;
406:     }
407:     break;
408:   case PEP_BASIS_LAGUERRE:
409:     ca[0] = -1.0; cb[0] = 0.0; cg[0] = 0.0;
410:     for (k=1;k<nmat;k++) {
411:       ca[k] = -(k+1); cb[k] = 2*k+1; cg[k] = -k;
412:     }
413:     break;
414:   case PEP_BASIS_HERMITE:
415:     ca[0] = .5; cb[0] = 0.0; cg[0] = 0.0;
416:     for (k=1;k<nmat;k++) {
417:       ca[k] = .5; cb[k] = 0.0; cg[k] = -k;
418:     }
419:     break;
420:   }
421:   PetscFunctionReturn(0);
422: }