Actual source code: pepdefault.c

  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

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

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

 29:    Level: developer

 31: .seealso: [](ch:pep), `PEPSetUp()`
 32: @*/
 33: PetscErrorCode PEPSetWorkVecs(PEP pep,PetscInt nw)
 34: {
 35:   Vec            t;

 37:   PetscFunctionBegin;
 40:   PetscCheck(nw>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"nw must be > 0: nw = %" PetscInt_FMT,nw);
 41:   if (pep->nwork < nw) {
 42:     PetscCall(VecDestroyVecs(pep->nwork,&pep->work));
 43:     pep->nwork = nw;
 44:     PetscCall(BVGetColumn(pep->V,0,&t));
 45:     PetscCall(VecDuplicateVecs(t,nw,&pep->work));
 46:     PetscCall(BVRestoreColumn(pep->V,0,&t));
 47:   }
 48:   PetscFunctionReturn(PETSC_SUCCESS);
 49: }

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

 58:   PetscFunctionBegin;
 59:   w = SlepcAbsEigenvalue(eigr,eigi);
 60:   *errest = (w!=0.0)? res/w: PETSC_MAX_REAL;
 61:   PetscFunctionReturn(PETSC_SUCCESS);
 62: }

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

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

 90: /*
 91:   PEPSetWhichEigenpairs_Default - Sets the default value for which,
 92:   depending on the ST.
 93:  */
 94: PetscErrorCode PEPSetWhichEigenpairs_Default(PEP pep)
 95: {
 96:   PetscBool      target;

 98:   PetscFunctionBegin;
 99:   PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&target));
100:   if (target) pep->which = PEP_TARGET_MAGNITUDE;
101:   else pep->which = PEP_LARGEST_MAGNITUDE;
102:   PetscFunctionReturn(PETSC_SUCCESS);
103: }

105: /*
106:   PEPConvergedAbsolute - Checks convergence absolutely.
107: */
108: PetscErrorCode PEPConvergedAbsolute(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
109: {
110:   PetscFunctionBegin;
111:   *errest = res;
112:   PetscFunctionReturn(PETSC_SUCCESS);
113: }

115: /*@C
116:    PEPStoppingBasic - Default routine to determine whether the outer eigensolver
117:    iteration must be stopped.

119:    Collective

121:    Input Parameters:
122: +  pep    - the polynomial eigensolver context
123: .  its    - current number of iterations
124: .  max_it - maximum number of iterations
125: .  nconv  - number of currently converged eigenpairs
126: .  nev    - number of requested eigenpairs
127: -  ctx    - context (not used here)

129:    Output Parameter:
130: .  reason - result of the stopping test

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

136:    This is the default stopping test.
137:    Use `PEPSetStoppingTest()` to provide your own test instead of using this one.

139:    Level: advanced

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

157: PetscErrorCode PEPBackTransform_Default(PEP pep)
158: {
159:   PetscFunctionBegin;
160:   PetscCall(STBackTransform(pep->st,pep->nconv,pep->eigr,pep->eigi));
161:   PetscFunctionReturn(PETSC_SUCCESS);
162: }

164: PetscErrorCode PEPComputeVectors_Default(PEP pep)
165: {
166:   PetscInt       i;
167:   Vec            v;

169:   PetscFunctionBegin;
170:   PetscCall(PEPExtractVectors(pep));

172:   /* Fix eigenvectors if balancing was used */
173:   if ((pep->scale==PEP_SCALE_DIAGONAL || pep->scale==PEP_SCALE_BOTH) && pep->Dr && (pep->refine!=PEP_REFINE_MULTIPLE)) {
174:     for (i=0;i<pep->nconv;i++) {
175:       PetscCall(BVGetColumn(pep->V,i,&v));
176:       PetscCall(VecPointwiseMult(v,v,pep->Dr));
177:       PetscCall(BVRestoreColumn(pep->V,i,&v));
178:     }
179:   }

181:   /* normalization */
182:   PetscCall(BVNormalize(pep->V,pep->eigi));
183:   PetscFunctionReturn(PETSC_SUCCESS);
184: }

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

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

325: /*
326:    PEPComputeScaleFactor - compute sfactor as described in [Betcke 2008].
327: */
328: PetscErrorCode PEPComputeScaleFactor(PEP pep)
329: {
330:   PetscBool      has0,has1,flg;
331:   PetscReal      norm0,norm1;
332:   Mat            T[2];
333:   PEPBasis       basis;
334:   PetscInt       i;

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

381: /*
382:    PEPBasisCoefficients - compute polynomial basis coefficients
383: */
384: PetscErrorCode PEPBasisCoefficients(PEP pep,PetscReal *pbc)
385: {
386:   PetscReal *ca,*cb,*cg;
387:   PetscInt  k,nmat=pep->nmat;

389:   PetscFunctionBegin;
390:   ca = pbc;
391:   cb = pbc+nmat;
392:   cg = pbc+2*nmat;
393:   switch (pep->basis) {
394:   case PEP_BASIS_MONOMIAL:
395:     for (k=0;k<nmat;k++) {
396:       ca[k] = 1.0; cb[k] = 0.0; cg[k] = 0.0;
397:     }
398:     break;
399:   case PEP_BASIS_CHEBYSHEV1:
400:     ca[0] = 1.0; cb[0] = 0.0; cg[0] = 0.0;
401:     for (k=1;k<nmat;k++) {
402:       ca[k] = .5; cb[k] = 0.0; cg[k] = .5;
403:     }
404:     break;
405:   case PEP_BASIS_CHEBYSHEV2:
406:     ca[0] = .5; cb[0] = 0.0; cg[0] = 0.0;
407:     for (k=1;k<nmat;k++) {
408:       ca[k] = .5; cb[k] = 0.0; cg[k] = .5;
409:     }
410:     break;
411:   case PEP_BASIS_LEGENDRE:
412:     ca[0] = 1.0; cb[0] = 0.0; cg[0] = 0.0;
413:     for (k=1;k<nmat;k++) {
414:       ca[k] = k+1; cb[k] = -2*k; cg[k] = k;
415:     }
416:     break;
417:   case PEP_BASIS_LAGUERRE:
418:     ca[0] = -1.0; cb[0] = 0.0; cg[0] = 0.0;
419:     for (k=1;k<nmat;k++) {
420:       ca[k] = -(k+1); cb[k] = 2*k+1; cg[k] = -k;
421:     }
422:     break;
423:   case PEP_BASIS_HERMITE:
424:     ca[0] = .5; cb[0] = 0.0; cg[0] = 0.0;
425:     for (k=1;k<nmat;k++) {
426:       ca[k] = .5; cb[k] = 0.0; cg[k] = -k;
427:     }
428:     break;
429:   }
430:   PetscFunctionReturn(PETSC_SUCCESS);
431: }