Actual source code: krylovschur.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*

  3:    SLEPc eigensolver: "krylovschur"

  5:    Method: Krylov-Schur

  7:    Algorithm:

  9:        Single-vector Krylov-Schur method for non-symmetric problems,
 10:        including harmonic extraction.

 12:    References:

 14:        [1] "Krylov-Schur Methods in SLEPc", SLEPc Technical Report STR-7,
 15:            available at http://slepc.upv.es.

 17:        [2] G.W. Stewart, "A Krylov-Schur Algorithm for Large Eigenproblems",
 18:            SIAM J. Matrix Anal. App. 23(3):601-614, 2001.

 20:        [3] "Practical Implementation of Harmonic Krylov-Schur", SLEPc Technical
 21:             Report STR-9, available at http://slepc.upv.es.

 23:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 24:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 25:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

 27:    This file is part of SLEPc.

 29:    SLEPc is free software: you can redistribute it and/or modify it under  the
 30:    terms of version 3 of the GNU Lesser General Public License as published by
 31:    the Free Software Foundation.

 33:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 34:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 35:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 36:    more details.

 38:    You  should have received a copy of the GNU Lesser General  Public  License
 39:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 40:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 41: */

 43: #include <slepc/private/epsimpl.h>                /*I "slepceps.h" I*/
 44:  #include krylovschur.h

 48: PetscErrorCode EPSGetArbitraryValues(EPS eps,PetscScalar *rr,PetscScalar *ri)
 49: {
 51:   PetscInt       i,newi,ld,n,l;
 52:   Vec            xr=eps->work[0],xi=eps->work[1];
 53:   PetscScalar    re,im,*Zr,*Zi,*X;

 56:   DSGetLeadingDimension(eps->ds,&ld);
 57:   DSGetDimensions(eps->ds,&n,NULL,&l,NULL,NULL);
 58:   for (i=l;i<n;i++) {
 59:     re = eps->eigr[i];
 60:     im = eps->eigi[i];
 61:     STBackTransform(eps->st,1,&re,&im);
 62:     newi = i;
 63:     DSVectors(eps->ds,DS_MAT_X,&newi,NULL);
 64:     DSGetArray(eps->ds,DS_MAT_X,&X);
 65:     Zr = X+i*ld;
 66:     if (newi==i+1) Zi = X+newi*ld;
 67:     else Zi = NULL;
 68:     EPSComputeRitzVector(eps,Zr,Zi,eps->V,xr,xi);
 69:     DSRestoreArray(eps->ds,DS_MAT_X,&X);
 70:     (*eps->arbitrary)(re,im,xr,xi,rr+i,ri+i,eps->arbitraryctx);
 71:   }
 72:   return(0);
 73: }

 77: PetscErrorCode EPSSetUp_KrylovSchur(EPS eps)
 78: {
 79:   PetscErrorCode    ierr;
 80:   PetscReal         eta;
 81:   BVOrthogType      otype;
 82:   BVOrthogBlockType obtype;
 83:   EPS_KRYLOVSCHUR   *ctx = (EPS_KRYLOVSCHUR*)eps->data;
 84:   enum { EPS_KS_DEFAULT,EPS_KS_SYMM,EPS_KS_SLICE,EPS_KS_INDEF } variant;

 87:   /* spectrum slicing requires special treatment of default values */
 88:   if (eps->which==EPS_ALL) {
 89:     EPSSetUp_KrylovSchur_Slice(eps);
 90:   } else {
 91:     EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
 92:     if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
 93:     if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
 94:     if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
 95:   }
 96:   if (!ctx->lock && eps->mpd<eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");

 98:   if (eps->isgeneralized && eps->ishermitian && !eps->ispositive && eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not implemented for indefinite problems");
 99:   if (eps->ishermitian && eps->ispositive && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");

101:   if (!eps->extraction) {
102:     EPSSetExtraction(eps,EPS_RITZ);
103:   } else if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC)
104:     SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
105:   if (eps->extraction==EPS_HARMONIC && ctx->lock) { PetscInfo(eps,"Locking was requested but will be deactivated since is not supported with harmonic extraction\n"); }

107:   if (!ctx->keep) ctx->keep = 0.5;

109:   EPSAllocateSolution(eps,1);
110:   EPS_SetInnerProduct(eps);
111:   if (eps->arbitrary) {
112:     EPSSetWorkVecs(eps,2);
113:   } else if (eps->ishermitian && !eps->ispositive){
114:     EPSSetWorkVecs(eps,1);
115:   }

117:   /* dispatch solve method */
118:   if (eps->ishermitian) {
119:     if (eps->which==EPS_ALL) {
120:       if (eps->isgeneralized && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing not implemented for indefinite problems");
121:       else variant = EPS_KS_SLICE;
122:     } else if (eps->isgeneralized && !eps->ispositive) {
123:       variant = EPS_KS_INDEF;
124:     } else {
125:       switch (eps->extraction) {
126:         case EPS_RITZ:     variant = EPS_KS_SYMM; break;
127:         case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
128:         default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
129:       }
130:     }
131:   } else {
132:     switch (eps->extraction) {
133:       case EPS_RITZ:     variant = EPS_KS_DEFAULT; break;
134:       case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
135:       default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
136:     }
137:   }
138:   switch (variant) {
139:     case EPS_KS_DEFAULT:
140:       eps->ops->solve = EPSSolve_KrylovSchur_Default;
141:       eps->ops->computevectors = EPSComputeVectors_Schur;
142:       DSSetType(eps->ds,DSNHEP);
143:       DSAllocate(eps->ds,eps->ncv+1);
144:       break;
145:     case EPS_KS_SYMM:
146:       eps->ops->solve = EPSSolve_KrylovSchur_Symm;
147:       eps->ops->computevectors = EPSComputeVectors_Hermitian;
148:       DSSetType(eps->ds,DSHEP);
149:       DSSetCompact(eps->ds,PETSC_TRUE);
150:       DSSetExtraRow(eps->ds,PETSC_TRUE);
151:       DSAllocate(eps->ds,eps->ncv+1);
152:       break;
153:     case EPS_KS_SLICE:
154:       if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing does not support user-defined stopping test");
155:       eps->ops->solve = EPSSolve_KrylovSchur_Slice;
156:       eps->ops->computevectors = EPSComputeVectors_Slice;
157:       break;
158:     case EPS_KS_INDEF:
159:       eps->ops->solve = EPSSolve_KrylovSchur_Indefinite;
160:       eps->ops->computevectors = EPSComputeVectors_Indefinite;
161:       DSSetType(eps->ds,DSGHIEP);
162:       DSSetCompact(eps->ds,PETSC_TRUE);
163:       DSAllocate(eps->ds,eps->ncv+1);
164:       /* force reorthogonalization for pseudo-Lanczos */
165:       BVGetOrthogonalization(eps->V,&otype,NULL,&eta,&obtype);
166:       BVSetOrthogonalization(eps->V,otype,BV_ORTHOG_REFINE_ALWAYS,eta,obtype);
167:       break;
168:     default: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error");
169:   }
170:   return(0);
171: }

175: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
176: {
177:   PetscErrorCode  ierr;
178:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
179:   PetscInt        i,j,*pj,k,l,nv,ld,nconv;
180:   Mat             U;
181:   PetscScalar     *S,*Q,*g;
182:   PetscReal       beta,gamma=1.0;
183:   PetscBool       breakdown,harmonic;

186:   DSGetLeadingDimension(eps->ds,&ld);
187:   harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
188:   if (harmonic) { PetscMalloc1(ld,&g); }
189:   if (eps->arbitrary) pj = &j;
190:   else pj = NULL;

192:   /* Get the starting Arnoldi vector */
193:   EPSGetStartVector(eps,0,NULL);
194:   l = 0;

196:   /* Restart loop */
197:   while (eps->reason == EPS_CONVERGED_ITERATING) {
198:     eps->its++;

200:     /* Compute an nv-step Arnoldi factorization */
201:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
202:     DSGetArray(eps->ds,DS_MAT_A,&S);
203:     EPSBasicArnoldi(eps,PETSC_FALSE,S,ld,eps->nconv+l,&nv,&beta,&breakdown);
204:     DSRestoreArray(eps->ds,DS_MAT_A,&S);
205:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
206:     if (l==0) {
207:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
208:     } else {
209:       DSSetState(eps->ds,DS_STATE_RAW);
210:     }
211:     BVSetActiveColumns(eps->V,eps->nconv,nv);

213:     /* Compute translation of Krylov decomposition if harmonic extraction used */
214:     if (harmonic) {
215:       DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,g,&gamma);
216:     }

218:     /* Solve projected problem */
219:     DSSolve(eps->ds,eps->eigr,eps->eigi);
220:     if (eps->arbitrary) {
221:       EPSGetArbitraryValues(eps,eps->rr,eps->ri);
222:       j=1;
223:     }
224:     DSSort(eps->ds,eps->eigr,eps->eigi,eps->rr,eps->ri,pj);

226:     /* Check convergence */
227:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,gamma,&k);
228:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
229:     nconv = k;

231:     /* Update l */
232:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
233:     else {
234:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
235: #if !defined(PETSC_USE_COMPLEX)
236:       DSGetArray(eps->ds,DS_MAT_A,&S);
237:       if (S[k+l+(k+l-1)*ld] != 0.0) {
238:         if (k+l<nv-1) l = l+1;
239:         else l = l-1;
240:       }
241:       DSRestoreArray(eps->ds,DS_MAT_A,&S);
242: #endif
243:     }
244:     if ((!ctx->lock || harmonic) && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */

246:     if (eps->reason == EPS_CONVERGED_ITERATING) {
247:       if (breakdown) {
248:         /* Start a new Arnoldi factorization */
249:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
250:         if (k<eps->nev) {
251:           EPSGetStartVector(eps,k,&breakdown);
252:           if (breakdown) {
253:             eps->reason = EPS_DIVERGED_BREAKDOWN;
254:             PetscInfo(eps,"Unable to generate more start vectors\n");
255:           }
256:         }
257:       } else {
258:         /* Undo translation of Krylov decomposition */
259:         if (harmonic) {
260:           DSSetDimensions(eps->ds,nv,0,k,l);
261:           DSTranslateHarmonic(eps->ds,0.0,beta,PETSC_TRUE,g,&gamma);
262:           /* gamma u^ = u - U*g~ */
263:           BVMultColumn(eps->V,-1.0,1.0,nv,g);
264:           BVScaleColumn(eps->V,nv,1.0/gamma);
265:         }
266:         /* Prepare the Rayleigh quotient for restart */
267:         DSGetArray(eps->ds,DS_MAT_A,&S);
268:         DSGetArray(eps->ds,DS_MAT_Q,&Q);
269:         for (i=k;i<k+l;i++) {
270:           S[k+l+i*ld] = Q[nv-1+i*ld]*beta*gamma;
271:         }
272:         DSRestoreArray(eps->ds,DS_MAT_A,&S);
273:         DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
274:       }
275:     }
276:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
277:     DSGetMat(eps->ds,DS_MAT_Q,&U);
278:     BVMultInPlace(eps->V,U,eps->nconv,k+l);
279:     MatDestroy(&U);

281:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
282:       BVCopyColumn(eps->V,nv,k+l);
283:     }
284:     eps->nconv = k;
285:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv);
286:   }

288:   if (harmonic) { PetscFree(g); }
289:   /* truncate Schur decomposition and change the state to raw so that
290:      DSVectors() computes eigenvectors from scratch */
291:   DSSetDimensions(eps->ds,eps->nconv,0,0,0);
292:   DSSetState(eps->ds,DS_STATE_RAW);
293:   return(0);
294: }

298: static PetscErrorCode EPSKrylovSchurSetRestart_KrylovSchur(EPS eps,PetscReal keep)
299: {
300:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

303:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
304:   else {
305:     if (keep<0.1 || keep>0.9) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
306:     ctx->keep = keep;
307:   }
308:   return(0);
309: }

313: /*@
314:    EPSKrylovSchurSetRestart - Sets the restart parameter for the Krylov-Schur
315:    method, in particular the proportion of basis vectors that must be kept
316:    after restart.

318:    Logically Collective on EPS

320:    Input Parameters:
321: +  eps - the eigenproblem solver context
322: -  keep - the number of vectors to be kept at restart

324:    Options Database Key:
325: .  -eps_krylovschur_restart - Sets the restart parameter

327:    Notes:
328:    Allowed values are in the range [0.1,0.9]. The default is 0.5.

330:    Level: advanced

332: .seealso: EPSKrylovSchurGetRestart()
333: @*/
334: PetscErrorCode EPSKrylovSchurSetRestart(EPS eps,PetscReal keep)
335: {

341:   PetscTryMethod(eps,"EPSKrylovSchurSetRestart_C",(EPS,PetscReal),(eps,keep));
342:   return(0);
343: }

347: static PetscErrorCode EPSKrylovSchurGetRestart_KrylovSchur(EPS eps,PetscReal *keep)
348: {
349:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

352:   *keep = ctx->keep;
353:   return(0);
354: }

358: /*@
359:    EPSKrylovSchurGetRestart - Gets the restart parameter used in the
360:    Krylov-Schur method.

362:    Not Collective

364:    Input Parameter:
365: .  eps - the eigenproblem solver context

367:    Output Parameter:
368: .  keep - the restart parameter

370:    Level: advanced

372: .seealso: EPSKrylovSchurSetRestart()
373: @*/
374: PetscErrorCode EPSKrylovSchurGetRestart(EPS eps,PetscReal *keep)
375: {

381:   PetscUseMethod(eps,"EPSKrylovSchurGetRestart_C",(EPS,PetscReal*),(eps,keep));
382:   return(0);
383: }

387: static PetscErrorCode EPSKrylovSchurSetLocking_KrylovSchur(EPS eps,PetscBool lock)
388: {
389:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

392:   ctx->lock = lock;
393:   return(0);
394: }

398: /*@
399:    EPSKrylovSchurSetLocking - Choose between locking and non-locking variants of
400:    the Krylov-Schur method.

402:    Logically Collective on EPS

404:    Input Parameters:
405: +  eps  - the eigenproblem solver context
406: -  lock - true if the locking variant must be selected

408:    Options Database Key:
409: .  -eps_krylovschur_locking - Sets the locking flag

411:    Notes:
412:    The default is to lock converged eigenpairs when the method restarts.
413:    This behaviour can be changed so that all directions are kept in the
414:    working subspace even if already converged to working accuracy (the
415:    non-locking variant).

417:    Level: advanced

419: .seealso: EPSKrylovSchurGetLocking()
420: @*/
421: PetscErrorCode EPSKrylovSchurSetLocking(EPS eps,PetscBool lock)
422: {

428:   PetscTryMethod(eps,"EPSKrylovSchurSetLocking_C",(EPS,PetscBool),(eps,lock));
429:   return(0);
430: }

434: static PetscErrorCode EPSKrylovSchurGetLocking_KrylovSchur(EPS eps,PetscBool *lock)
435: {
436:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

439:   *lock = ctx->lock;
440:   return(0);
441: }

445: /*@
446:    EPSKrylovSchurGetLocking - Gets the locking flag used in the Krylov-Schur
447:    method.

449:    Not Collective

451:    Input Parameter:
452: .  eps - the eigenproblem solver context

454:    Output Parameter:
455: .  lock - the locking flag

457:    Level: advanced

459: .seealso: EPSKrylovSchurSetLocking()
460: @*/
461: PetscErrorCode EPSKrylovSchurGetLocking(EPS eps,PetscBool *lock)
462: {

468:   PetscUseMethod(eps,"EPSKrylovSchurGetLocking_C",(EPS,PetscBool*),(eps,lock));
469:   return(0);
470: }

474: static PetscErrorCode EPSKrylovSchurSetPartitions_KrylovSchur(EPS eps,PetscInt npart)
475: {
476:   PetscErrorCode  ierr;
477:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
478:   PetscMPIInt     size;

481:   if (ctx->npart!=npart) {
482:     if (ctx->commset) { PetscSubcommDestroy(&ctx->subc); }
483:     if (ctx->eps) { EPSDestroy(&ctx->eps); }
484:   } 
485:   if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
486:     ctx->npart = 1;
487:   } else {
488:     MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size);
489:     if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
490:     ctx->npart = npart;
491:   }
492:   eps->state = EPS_STATE_INITIAL;
493:   return(0);
494: }

498: /*@
499:    EPSKrylovSchurSetPartitions - Sets the number of partitions for the
500:    case of doing spectrum slicing for a computational interval with the
501:    communicator split in several sub-communicators.

503:    Logically Collective on EPS

505:    Input Parameters:
506: +  eps   - the eigenproblem solver context
507: -  npart - number of partitions

509:    Options Database Key:
510: .  -eps_krylovschur_partitions <npart> - Sets the number of partitions

512:    Notes:
513:    By default, npart=1 so all processes in the communicator participate in
514:    the processing of the whole interval. If npart>1 then the interval is
515:    divided into npart subintervals, each of them being processed by a
516:    subset of processes.

518:    The interval is split proportionally unless the separation points are
519:    specified with EPSKrylovSchurSetSubintervals().

521:    Level: advanced

523: .seealso: EPSKrylovSchurSetSubintervals(), EPSSetInterval()
524: @*/
525: PetscErrorCode EPSKrylovSchurSetPartitions(EPS eps,PetscInt npart)
526: {

532:   PetscTryMethod(eps,"EPSKrylovSchurSetPartitions_C",(EPS,PetscInt),(eps,npart));
533:   return(0);
534: }

538: static PetscErrorCode EPSKrylovSchurGetPartitions_KrylovSchur(EPS eps,PetscInt *npart)
539: {
540:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

543:   *npart  = ctx->npart;
544:   return(0);
545: }

549: /*@
550:    EPSKrylovSchurGetPartitions - Gets the number of partitions of the
551:    communicator in case of spectrum slicing.

553:    Not Collective

555:    Input Parameter:
556: .  eps - the eigenproblem solver context

558:    Output Parameter:
559: .  npart - number of partitions

561:    Level: advanced

563: .seealso: EPSKrylovSchurSetPartitions()
564: @*/
565: PetscErrorCode EPSKrylovSchurGetPartitions(EPS eps,PetscInt *npart)
566: {

572:   PetscUseMethod(eps,"EPSKrylovSchurGetPartitions_C",(EPS,PetscInt*),(eps,npart));
573:   return(0);
574: }

578: static PetscErrorCode EPSKrylovSchurSetDetectZeros_KrylovSchur(EPS eps,PetscBool detect)
579: {
580:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

583:   ctx->detect = detect;
584:   eps->state  = EPS_STATE_INITIAL;
585:   return(0);
586: }

590: /*@
591:    EPSKrylovSchurSetDetectZeros - Sets a flag to enforce detection of
592:    zeros during the factorizations throughout the spectrum slicing computation.

594:    Logically Collective on EPS

596:    Input Parameters:
597: +  eps    - the eigenproblem solver context
598: -  detect - check for zeros

600:    Options Database Key:
601: .  -eps_krylovschur_detect_zeros - Check for zeros; this takes an optional
602:    bool value (0/1/no/yes/true/false)

604:    Notes:
605:    A zero in the factorization indicates that a shift coincides with an eigenvalue.

607:    This flag is turned off by default, and may be necessary in some cases,
608:    especially when several partitions are being used. This feature currently
609:    requires an external package for factorizations with support for zero
610:    detection, e.g. MUMPS.

612:    Level: advanced

614: .seealso: EPSKrylovSchurSetPartitions(), EPSSetInterval()
615: @*/
616: PetscErrorCode EPSKrylovSchurSetDetectZeros(EPS eps,PetscBool detect)
617: {

623:   PetscTryMethod(eps,"EPSKrylovSchurSetDetectZeros_C",(EPS,PetscBool),(eps,detect));
624:   return(0);
625: }

629: static PetscErrorCode EPSKrylovSchurGetDetectZeros_KrylovSchur(EPS eps,PetscBool *detect)
630: {
631:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

634:   *detect = ctx->detect;
635:   return(0);
636: }

640: /*@
641:    EPSKrylovSchurGetDetectZeros - Gets the flag that enforces zero detection
642:    in spectrum slicing.

644:    Not Collective

646:    Input Parameter:
647: .  eps - the eigenproblem solver context

649:    Output Parameter:
650: .  detect - whether zeros detection is enforced during factorizations

652:    Level: advanced

654: .seealso: EPSKrylovSchurSetDetectZeros()
655: @*/
656: PetscErrorCode EPSKrylovSchurGetDetectZeros(EPS eps,PetscBool *detect)
657: {

663:   PetscUseMethod(eps,"EPSKrylovSchurGetDetectZeros_C",(EPS,PetscBool*),(eps,detect));
664:   return(0);
665: }

669: static PetscErrorCode EPSKrylovSchurSetDimensions_KrylovSchur(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
670: {
671:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

674:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
675:   ctx->nev = nev;
676:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
677:     ctx->ncv = 0;
678:   } else {
679:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
680:     ctx->ncv = ncv;
681:   }
682:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
683:     ctx->mpd = 0;
684:   } else {
685:     if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
686:     ctx->mpd = mpd;
687:   }
688:   eps->state = EPS_STATE_INITIAL;
689:   return(0);
690: }

694: /*@
695:    EPSKrylovSchurSetDimensions - Sets the dimensions used for each subsolve
696:    step in case of doing spectrum slicing for a computational interval.
697:    The meaning of the parameters is the same as in EPSSetDimensions().

699:    Logically Collective on EPS

701:    Input Parameters:
702: +  eps - the eigenproblem solver context
703: .  nev - number of eigenvalues to compute
704: .  ncv - the maximum dimension of the subspace to be used by the subsolve
705: -  mpd - the maximum dimension allowed for the projected problem

707:    Options Database Key:
708: +  -eps_krylovschur_nev <nev> - Sets the number of eigenvalues
709: .  -eps_krylovschur_ncv <ncv> - Sets the dimension of the subspace
710: -  -eps_krylovschur_mpd <mpd> - Sets the maximum projected dimension

712:    Level: advanced

714: .seealso: EPSKrylovSchurGetDimensions(), EPSSetDimensions(), EPSSetInterval()
715: @*/
716: PetscErrorCode EPSKrylovSchurSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
717: {

725:   PetscTryMethod(eps,"EPSKrylovSchurSetDimensions_C",(EPS,PetscInt,PetscInt,PetscInt),(eps,nev,ncv,mpd));
726:   return(0);
727: }

731: static PetscErrorCode EPSKrylovSchurGetDimensions_KrylovSchur(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
732: {
733:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

736:   if (nev) *nev = ctx->nev;
737:   if (ncv) *ncv = ctx->ncv;
738:   if (mpd) *mpd = ctx->mpd;
739:   return(0);
740: }

744: /*@
745:    EPSKrylovSchurGetDimensions - Gets the dimensions used for each subsolve
746:    step in case of doing spectrum slicing for a computational interval.

748:    Not Collective

750:    Input Parameter:
751: .  eps - the eigenproblem solver context

753:    Output Parameters:
754: +  nev - number of eigenvalues to compute
755: .  ncv - the maximum dimension of the subspace to be used by the subsolve
756: -  mpd - the maximum dimension allowed for the projected problem

758:    Level: advanced

760: .seealso: EPSKrylovSchurSetDimensions()
761: @*/
762: PetscErrorCode EPSKrylovSchurGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
763: {

768:   PetscUseMethod(eps,"EPSKrylovSchurGetDimensions_C",(EPS,PetscInt*,PetscInt*,PetscInt*),(eps,nev,ncv,mpd));
769:   return(0);
770: }

774: static PetscErrorCode EPSKrylovSchurSetSubintervals_KrylovSchur(EPS eps,PetscReal* subint)
775: {
776:   PetscErrorCode  ierr;
777:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
778:   PetscInt        i;

781:   if (subint[0]!=eps->inta || subint[ctx->npart]!=eps->intb) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"First and last values must match the endpoints of EPSSetInterval()");
782:   for (i=0;i<ctx->npart;i++) if (subint[i]>=subint[i+1]) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Array must contain values in ascending order");
783:   if (ctx->subintervals) { PetscFree(ctx->subintervals); }
784:   PetscMalloc1(ctx->npart+1,&ctx->subintervals);
785:   for (i=0;i<ctx->npart+1;i++) ctx->subintervals[i] = subint[i];
786:   ctx->subintset = PETSC_TRUE;
787:   eps->state = EPS_STATE_INITIAL;
788:   return(0);
789: }

793: /*@C
794:    EPSKrylovSchurSetSubintervals - Sets the points that delimit the
795:    subintervals to be used in spectrum slicing with several partitions.

797:    Logically Collective on EPS

799:    Input Parameters:
800: +  eps    - the eigenproblem solver context
801: -  subint - array of real values specifying subintervals

803:    Notes:
804:    This function must be called after EPSKrylovSchurSetPartitions(). For npart
805:    partitions, the argument subint must contain npart+1 real values sorted in
806:    ascending order: subint_0, subint_1, ..., subint_npart, where the first
807:    and last values must coincide with the interval endpoints set with
808:    EPSSetInterval().

810:    The subintervals are then defined by two consecutive points: [subint_0,subint_1],
811:    [subint_1,subint_2], and so on.

813:    Level: advanced

815: .seealso: EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubintervals(), EPSSetInterval()
816: @*/
817: PetscErrorCode EPSKrylovSchurSetSubintervals(EPS eps,PetscReal *subint)
818: {

823:   PetscTryMethod(eps,"EPSKrylovSchurSetSubintervals_C",(EPS,PetscReal*),(eps,subint));
824:   return(0);
825: }

829: static PetscErrorCode EPSKrylovSchurGetSubintervals_KrylovSchur(EPS eps,PetscReal **subint)
830: {
831:   PetscErrorCode  ierr;
832:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
833:   PetscInt        i;

836:   if (!ctx->subintset) {
837:     if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
838:     if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
839:   }
840:   PetscMalloc1(ctx->npart+1,subint);
841:   for (i=0;i<=ctx->npart;i++) (*subint)[i] = ctx->subintervals[i];
842:   return(0);
843: }

847: /*@C
848:    EPSKrylovSchurGetSubintervals - Returns the points that delimit the
849:    subintervals used in spectrum slicing with several partitions.

851:    Logically Collective on EPS

853:    Input Parameter:
854: .  eps    - the eigenproblem solver context

856:    Output Parameter:
857: .  subint - array of real values specifying subintervals

859:    Notes:
860:    If the user passed values with EPSKrylovSchurSetSubintervals(), then the
861:    same values are returned. Otherwise, the values computed internally are
862:    obtained.

864:    This function is only available for spectrum slicing runs.

866:    The returned array has length npart+1 (see EPSKrylovSchurGetPartitions())
867:    and should be freed by the user.

869:    Fortran Notes:
870:    The calling sequence from Fortran is
871: .vb
872:    EPSKrylovSchurGetSubintervals(eps,subint,ierr)
873:    double precision subint(npart+1) output
874: .ve

876:    Level: advanced

878: .seealso: EPSKrylovSchurSetSubintervals(), EPSKrylovSchurGetPartitions(), EPSSetInterval()
879: @*/
880: PetscErrorCode EPSKrylovSchurGetSubintervals(EPS eps,PetscReal** subint)
881: {

887:   PetscUseMethod(eps,"EPSKrylovSchurGetSubintervals_C",(EPS,PetscReal**),(eps,subint));
888:   return(0);
889: }

893: static PetscErrorCode EPSKrylovSchurGetInertias_KrylovSchur(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
894: {
895:   PetscErrorCode  ierr;
896:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
897:   PetscInt        i;
898:   EPS_SR          sr = ctx->sr;

901:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
902:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
903:   switch (eps->state) {
904:   case EPS_STATE_INITIAL:
905:     break;
906:   case EPS_STATE_SETUP:
907:     *n = ctx->npart+1;
908:     PetscMalloc1(*n,shifts);
909:     PetscMalloc1(*n,inertias);
910:     (*shifts)[0] = eps->inta;
911:     (*inertias)[0] = (sr->dir==1)?sr->inertia0:sr->inertia1;
912:     if (ctx->npart==1) {
913:       (*shifts)[1] = eps->intb;
914:       (*inertias)[1] = (sr->dir==1)?sr->inertia1:sr->inertia0;
915:     } else {
916:       for (i=1;i<*n;i++) {
917:         (*shifts)[i] = ctx->subintervals[i];
918:         (*inertias)[i] = (*inertias)[i-1]+ctx->nconv_loc[i-1];
919:       }
920:     }
921:     break;
922:   case EPS_STATE_SOLVED:
923:   case EPS_STATE_EIGENVECTORS:
924:     *n = ctx->nshifts;
925:     PetscMalloc1(*n,shifts);
926:     PetscMalloc1(*n,inertias);
927:     for (i=0;i<*n;i++) {
928:       (*shifts)[i] = ctx->shifts[i];
929:       (*inertias)[i] = ctx->inertias[i];
930:     }
931:     break;
932:   }
933:   return(0);
934: }

938: /*@C
939:    EPSKrylovSchurGetInertias - Gets the values of the shifts and their
940:    corresponding inertias in case of doing spectrum slicing for a
941:    computational interval.

943:    Not Collective

945:    Input Parameter:
946: .  eps - the eigenproblem solver context

948:    Output Parameters:
949: +  n        - number of shifts, including the endpoints of the interval
950: .  shifts   - the values of the shifts used internally in the solver
951: -  inertias - the values of the inertia in each shift

953:    Notes:
954:    If called after EPSSolve(), all shifts used internally by the solver are
955:    returned (including both endpoints and any intermediate ones). If called
956:    before EPSSolve() and after EPSSetUp() then only the information of the
957:    endpoints of subintervals is available.

959:    This function is only available for spectrum slicing runs.

961:    The returned arrays should be freed by the user.

963:    Level: advanced

965: .seealso: EPSSetInterval(), EPSKrylovSchurSetSubintervals()
966: @*/
967: PetscErrorCode EPSKrylovSchurGetInertias(EPS eps,PetscInt *n,PetscReal **shifts,PetscInt **inertias)
968: {

974:   PetscUseMethod(eps,"EPSKrylovSchurGetInertias_C",(EPS,PetscInt*,PetscReal**,PetscInt**),(eps,n,shifts,inertias));
975:   return(0);
976: }

980: static PetscErrorCode EPSKrylovSchurGetSubcommInfo_KrylovSchur(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
981: {
982:   PetscErrorCode  ierr;
983:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
984:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

987:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
988:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
989:   if (k) *k = (ctx->npart==1)? 0: ctx->subc->color;
990:   if (n) *n = sr->numEigs;
991:   if (v) {
992:     BVCreateVec(sr->V,v);
993:   }
994:   return(0);
995: }

999: /*@C
1000:    EPSKrylovSchurGetSubcommInfo - Gets information related to the case of
1001:    doing spectrum slicing for a computational interval with multiple
1002:    communicators.

1004:    Collective on the subcommunicator (if v is given)

1006:    Input Parameter:
1007: .  eps - the eigenproblem solver context

1009:    Output Parameters:
1010: +  k - index of the subinterval for the calling process
1011: .  n - number of eigenvalues found in the k-th subinterval
1012: -  v - a vector owned by processes in the subcommunicator with dimensions
1013:        compatible for locally computed eigenvectors (or NULL)

1015:    Notes:
1016:    This function is only available for spectrum slicing runs.

1018:    The returned Vec should be destroyed by the user.

1020:    Level: advanced

1022: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommPairs()
1023: @*/
1024: PetscErrorCode EPSKrylovSchurGetSubcommInfo(EPS eps,PetscInt *k,PetscInt *n,Vec *v)
1025: {

1030:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommInfo_C",(EPS,PetscInt*,PetscInt*,Vec*),(eps,k,n,v));
1031:   return(0);
1032: }

1036: static PetscErrorCode EPSKrylovSchurGetSubcommPairs_KrylovSchur(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1037: {
1038:   PetscErrorCode  ierr;
1039:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1040:   EPS_SR          sr = ((EPS_KRYLOVSCHUR*)ctx->eps->data)->sr;

1043:   EPSCheckSolved(eps,1);
1044:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1045:   if (i<0 || i>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
1046:   if (eig) *eig = sr->eigr[sr->perm[i]];
1047:   if (v) { BVCopyVec(sr->V,sr->perm[i],v); }
1048:   return(0);
1049: }

1053: /*@C
1054:    EPSKrylovSchurGetSubcommPairs - Gets the i-th eigenpair stored
1055:    internally in the subcommunicator to which the calling process belongs.

1057:    Collective on the subcommunicator (if v is given)

1059:    Input Parameter:
1060: +  eps - the eigenproblem solver context
1061: -  i   - index of the solution

1063:    Output Parameters:
1064: +  eig - the eigenvalue
1065: -  v   - the eigenvector

1067:    Notes:
1068:    It is allowed to pass NULL for v if the eigenvector is not required.
1069:    Otherwise, the caller must provide a valid Vec objects, i.e.,
1070:    it must be created by the calling program with EPSKrylovSchurGetSubcommInfo().

1072:    The index i should be a value between 0 and n-1, where n is the number of
1073:    vectors in the local subinterval, see EPSKrylovSchurGetSubcommInfo().

1075:    Level: advanced

1077: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo(), EPSKrylovSchurGetSubcommMats()
1078: @*/
1079: PetscErrorCode EPSKrylovSchurGetSubcommPairs(EPS eps,PetscInt i,PetscScalar *eig,Vec v)
1080: {

1086:   PetscUseMethod(eps,"EPSKrylovSchurGetSubcommPairs_C",(EPS,PetscInt,PetscScalar*,Vec),(eps,i,eig,v));
1087:   return(0);
1088: }

1092: static PetscErrorCode EPSKrylovSchurGetSubcommMats_KrylovSchur(EPS eps,Mat *A,Mat *B)
1093: {
1094:   PetscErrorCode  ierr;
1095:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

1098:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1099:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
1100:   EPSGetOperators(ctx->eps,A,B);
1101:   return(0);
1102: }

1106: /*@C
1107:    EPSKrylovSchurGetSubcommMats - Gets the eigenproblem matrices stored
1108:    internally in the subcommunicator to which the calling process belongs.

1110:    Collective on the subcommunicator

1112:    Input Parameter:
1113: .  eps - the eigenproblem solver context

1115:    Output Parameters:
1116: +  A  - the matrix associated with the eigensystem
1117: -  B  - the second matrix in the case of generalized eigenproblems

1119:    Notes:
1120:    This is the analog of EPSGetOperators(), but returns the matrices distributed
1121:    differently (in the subcommunicator rather than in the parent communicator).

1123:    These matrices should not be modified by the user.

1125:    Level: advanced

1127: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommInfo()
1128: @*/
1129: PetscErrorCode EPSKrylovSchurGetSubcommMats(EPS eps,Mat *A,Mat *B)
1130: {

1135:   PetscTryMethod(eps,"EPSKrylovSchurGetSubcommMats_C",(EPS,Mat*,Mat*),(eps,A,B));
1136:   return(0);
1137: }

1141: static PetscErrorCode EPSKrylovSchurUpdateSubcommMats_KrylovSchur(EPS eps,PetscScalar a,PetscScalar ap,Mat Au,PetscScalar b,PetscScalar bp, Mat Bu,MatStructure str,PetscBool globalup)
1142: {
1143:   PetscErrorCode  ierr;
1144:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data,*subctx;
1145:   Mat             A,B=NULL,Ag,Bg=NULL;
1146:   PetscBool       reuse=PETSC_TRUE;
1147:    
1149:   if (!ctx->sr) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Only available in interval computations, see EPSSetInterval()");
1150:   if (!eps->state) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSetUp() first");
1151:   EPSGetOperators(eps,&Ag,&Bg);
1152:   EPSGetOperators(ctx->eps,&A,&B);
1153:   
1154:   MatScale(A,a);
1155:   if (Au) {
1156:     MatAXPY(A,ap,Au,str);
1157:   }
1158:   if (B) MatScale(B,b);
1159:   if (Bu) {
1160:     MatAXPY(B,bp,Bu,str);
1161:   }
1162:   EPSSetOperators(ctx->eps,A,B);

1164:   /* Update stored matrix state */
1165:   subctx = (EPS_KRYLOVSCHUR*)ctx->eps->data;
1166:   PetscObjectStateGet((PetscObject)A,&subctx->Astate);
1167:   if (B) { PetscObjectStateGet((PetscObject)B,&subctx->Bstate); }

1169:   /* Update matrices in the parent communicator if requested by user */
1170:   if (globalup) {
1171:     if (ctx->npart>1) {
1172:       if (!ctx->isrow) {
1173:         MatGetOwnershipIS(Ag,&ctx->isrow,&ctx->iscol);
1174:         reuse = PETSC_FALSE;
1175:       }
1176:       if (str==DIFFERENT_NONZERO_PATTERN) reuse = PETSC_FALSE;
1177:       if (ctx->submata && !reuse) {
1178:         MatDestroyMatrices(1,&ctx->submata);
1179:       }
1180:       MatGetSubMatrices(A,1,&ctx->isrow,&ctx->iscol,(reuse)?MAT_REUSE_MATRIX:MAT_INITIAL_MATRIX,&ctx->submata);
1181:       MatCreateMPIMatConcatenateSeqMat(((PetscObject)Ag)->comm,ctx->submata[0],PETSC_DECIDE,MAT_REUSE_MATRIX,&Ag);
1182:       if (B) {
1183:         if (ctx->submatb && !reuse) {
1184:           MatDestroyMatrices(1,&ctx->submatb);
1185:         }
1186:         MatGetSubMatrices(B,1,&ctx->isrow,&ctx->iscol,(reuse)?MAT_REUSE_MATRIX:MAT_INITIAL_MATRIX,&ctx->submatb);
1187:         MatCreateMPIMatConcatenateSeqMat(((PetscObject)Bg)->comm,ctx->submatb[0],PETSC_DECIDE,MAT_REUSE_MATRIX,&Bg);
1188:       }
1189:     }
1190:     PetscObjectStateGet((PetscObject)Ag,&ctx->Astate);
1191:     if (Bg) { PetscObjectStateGet((PetscObject)Bg,&ctx->Bstate); }
1192:   }
1193:   EPSSetOperators(eps,Ag,Bg);
1194:   return(0);
1195: }

1199: /*@C
1200:    EPSKrylovSchurUpdateSubcommMats - Update the eigenproblem matrices stored
1201:    internally in the subcommunicator to which the calling process belongs.

1203:    Collective on EPS

1205:    Input Parameters:
1206: +  eps - the eigenproblem solver context
1207: .  s   - scalar that multiplies the existing A matrix
1208: .  a   - scalar used in the axpy operation on A
1209: .  Au  - matrix used in the axpy operation on A
1210: .  t   - scalar that multiplies the existing B matrix
1211: .  b   - scalar used in the axpy operation on B
1212: .  Bu  - matrix used in the axpy operation on B
1213: .  str - structure flag
1214: -  globalup - flag indicating if global matrices must be updated

1216:    Notes:
1217:    This function modifies the eigenproblem matrices at the subcommunicator level,
1218:    and optionally updates the global matrices in the parent communicator. The updates
1219:    are expressed as A <-- s*A + a*Au,  B <-- t*B + b*Bu.

1221:    It is possible to update one of the matrices, or both.

1223:    The matrices Au and Bu must be equal in all subcommunicators.

1225:    The str flag is passed to the MatAXPY() operations to perform the updates.

1227:    If globalup is true, communication is carried out to reconstruct the updated
1228:    matrices in the parent communicator. The user must be warned that if global
1229:    matrices are not in sync with subcommunicator matrices, the errors computed
1230:    by EPSComputeError() will be wrong even if the computed solution is correct
1231:    (the synchronization may be done only once at the end).

1233:    Level: advanced

1235: .seealso: EPSSetInterval(), EPSKrylovSchurSetPartitions(), EPSKrylovSchurGetSubcommMats()
1236: @*/
1237: PetscErrorCode EPSKrylovSchurUpdateSubcommMats(EPS eps,PetscScalar s,PetscScalar a,Mat Au,PetscScalar t,PetscScalar b, Mat Bu,MatStructure str,PetscBool globalup)
1238: {

1251:   PetscTryMethod(eps,"EPSKrylovSchurUpdateSubcommMats_C",(EPS,PetscScalar,PetscScalar,Mat,PetscScalar,PetscScalar,Mat,MatStructure,PetscBool),(eps,s,a,Au,t,b,Bu,str,globalup));
1252:   return(0);
1253: }

1257: PetscErrorCode EPSSetFromOptions_KrylovSchur(PetscOptionItems *PetscOptionsObject,EPS eps)
1258: {
1259:   PetscErrorCode  ierr;
1260:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1261:   PetscBool       flg,lock,b,f1,f2,f3;
1262:   PetscReal       keep;
1263:   PetscInt        i,j,k;

1266:   PetscOptionsHead(PetscOptionsObject,"EPS Krylov-Schur Options");
1267:   PetscOptionsReal("-eps_krylovschur_restart","Proportion of vectors kept after restart","EPSKrylovSchurSetRestart",0.5,&keep,&flg);
1268:   if (flg) {
1269:     EPSKrylovSchurSetRestart(eps,keep);
1270:   }
1271:   PetscOptionsBool("-eps_krylovschur_locking","Choose between locking and non-locking variants","EPSKrylovSchurSetLocking",PETSC_TRUE,&lock,&flg);
1272:   if (flg) {
1273:     EPSKrylovSchurSetLocking(eps,lock);
1274:   }
1275:   i = ctx->npart;
1276:   PetscOptionsInt("-eps_krylovschur_partitions","Number of partitions of the communicator for spectrum slicing","EPSKrylovSchurSetPartitions",ctx->npart,&i,&flg);
1277:   if (flg) {
1278:     EPSKrylovSchurSetPartitions(eps,i);
1279:   }
1280:   b = ctx->detect;
1281:   PetscOptionsBool("-eps_krylovschur_detect_zeros","Check zeros during factorizations at subinterval boundaries","EPSKrylovSchurSetDetectZeros",ctx->detect,&b,&flg);
1282:   if (flg) {
1283:     EPSKrylovSchurSetDetectZeros(eps,b);
1284:   }
1285:   i = 1;
1286:   j = k = PETSC_DECIDE;
1287:   PetscOptionsInt("-eps_krylovschur_nev","Number of eigenvalues to compute in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",40,&i,&f1);
1288:   PetscOptionsInt("-eps_krylovschur_ncv","Number of basis vectors in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&j,&f2);
1289:   PetscOptionsInt("-eps_krylovschur_mpd","Maximum dimension of projected problem in each subsolve (only for spectrum slicing)","EPSKrylovSchurSetDimensions",80,&k,&f3);
1290:   if (f1 || f2 || f3) {
1291:     EPSKrylovSchurSetDimensions(eps,i,j,k);
1292:   }
1293:   PetscOptionsTail();
1294:   return(0);
1295: }

1299: PetscErrorCode EPSView_KrylovSchur(EPS eps,PetscViewer viewer)
1300: {
1301:   PetscErrorCode  ierr;
1302:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
1303:   PetscBool       isascii;

1306:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1307:   if (isascii) {
1308:     PetscViewerASCIIPrintf(viewer,"  Krylov-Schur: %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
1309:     PetscViewerASCIIPrintf(viewer,"  Krylov-Schur: using the %slocking variant\n",ctx->lock?"":"non-");
1310:     if (eps->which==EPS_ALL) {
1311:       PetscViewerASCIIPrintf(viewer,"  Krylov-Schur: doing spectrum slicing with nev=%D, ncv=%D, mpd=%D\n",ctx->nev,ctx->ncv,ctx->mpd);
1312:       if (ctx->npart>1) {
1313:         PetscViewerASCIIPrintf(viewer,"  Krylov-Schur: multi-communicator spectrum slicing with %D partitions\n",ctx->npart);
1314:         if (ctx->detect) { PetscViewerASCIIPrintf(viewer,"  Krylov-Schur: detecting zeros when factorizing at subinterval boundaries\n"); }
1315:       }
1316:     }
1317:   }
1318:   return(0);
1319: }

1323: PetscErrorCode EPSDestroy_KrylovSchur(EPS eps)
1324: {

1328:   PetscFree(eps->data);
1329:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",NULL);
1330:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",NULL);
1331:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetLocking_C",NULL);
1332:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetLocking_C",NULL);
1333:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetPartitions_C",NULL);
1334:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetPartitions_C",NULL);
1335:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDetectZeros_C",NULL);
1336:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDetectZeros_C",NULL);
1337:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",NULL);
1338:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",NULL);
1339:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetSubintervals_C",NULL);
1340:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubintervals_C",NULL);
1341:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",NULL);
1342:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommInfo_C",NULL);
1343:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommPairs_C",NULL);
1344:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommMats_C",NULL);
1345:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurUpdateSubcommMats_C",NULL);
1346:   return(0);
1347: }

1351: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
1352: {
1355:   if (eps->which==EPS_ALL) {
1356:     EPSReset_KrylovSchur_Slice(eps);
1357:   }
1358:   return(0);
1359: }

1363: PETSC_EXTERN PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
1364: {
1365:   EPS_KRYLOVSCHUR *ctx;
1366:   PetscErrorCode  ierr;

1369:   PetscNewLog(eps,&ctx);
1370:   eps->data   = (void*)ctx;
1371:   ctx->lock   = PETSC_TRUE;
1372:   ctx->nev    = 1;
1373:   ctx->npart  = 1;
1374:   ctx->detect = PETSC_FALSE;
1375:   ctx->global = PETSC_TRUE;

1377:   eps->ops->setup          = EPSSetUp_KrylovSchur;
1378:   eps->ops->setfromoptions = EPSSetFromOptions_KrylovSchur;
1379:   eps->ops->destroy        = EPSDestroy_KrylovSchur;
1380:   eps->ops->reset          = EPSReset_KrylovSchur;
1381:   eps->ops->view           = EPSView_KrylovSchur;
1382:   eps->ops->backtransform  = EPSBackTransform_Default;
1383:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",EPSKrylovSchurSetRestart_KrylovSchur);
1384:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",EPSKrylovSchurGetRestart_KrylovSchur);
1385:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetLocking_C",EPSKrylovSchurSetLocking_KrylovSchur);
1386:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetLocking_C",EPSKrylovSchurGetLocking_KrylovSchur);
1387:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetPartitions_C",EPSKrylovSchurSetPartitions_KrylovSchur);
1388:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetPartitions_C",EPSKrylovSchurGetPartitions_KrylovSchur);
1389:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDetectZeros_C",EPSKrylovSchurSetDetectZeros_KrylovSchur);
1390:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDetectZeros_C",EPSKrylovSchurGetDetectZeros_KrylovSchur);
1391:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetDimensions_C",EPSKrylovSchurSetDimensions_KrylovSchur);
1392:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetDimensions_C",EPSKrylovSchurGetDimensions_KrylovSchur);
1393:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetSubintervals_C",EPSKrylovSchurSetSubintervals_KrylovSchur);
1394:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubintervals_C",EPSKrylovSchurGetSubintervals_KrylovSchur);
1395:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetInertias_C",EPSKrylovSchurGetInertias_KrylovSchur);
1396:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommInfo_C",EPSKrylovSchurGetSubcommInfo_KrylovSchur);
1397:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommPairs_C",EPSKrylovSchurGetSubcommPairs_KrylovSchur);
1398:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetSubcommMats_C",EPSKrylovSchurGetSubcommMats_KrylovSchur);
1399:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurUpdateSubcommMats_C",EPSKrylovSchurUpdateSubcommMats_KrylovSchur);
1400:   return(0);
1401: }