Actual source code: ks-indef.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:    SLEPc eigensolver: "krylovschur"

 13:    Method: Krylov-Schur for symmetric-indefinite eigenproblems
 14: */
 15: #include <slepc/private/epsimpl.h>
 16: #include "krylovschur.h"

 18: PetscErrorCode EPSSolve_KrylovSchur_Indefinite(EPS eps)
 19: {
 20:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
 21:   PetscInt        i,k,l,ld,nv,t,nconv=0;
 22:   Mat             U;
 23:   Vec             vomega,w=eps->work[0];
 24:   PetscScalar     *aux;
 25:   PetscReal       *a,*b,beta,beta1=1.0,*omega;
 26:   PetscBool       breakdown=PETSC_FALSE,symmlost=PETSC_FALSE;

 28:   DSGetLeadingDimension(eps->ds,&ld);

 30:   /* Get the starting Lanczos vector */
 31:   EPSGetStartVector(eps,0,NULL);

 33:   /* Extract sigma[0] from BV, computed during normalization */
 34:   BVSetActiveColumns(eps->V,0,1);
 35:   VecCreateSeq(PETSC_COMM_SELF,1,&vomega);
 36:   BVGetSignature(eps->V,vomega);
 37:   VecGetArray(vomega,&aux);
 38:   DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
 39:   omega[0] = PetscRealPart(aux[0]);
 40:   DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
 41:   VecRestoreArray(vomega,&aux);
 42:   VecDestroy(&vomega);
 43:   l = 0;

 45:   /* Restart loop */
 46:   while (eps->reason == EPS_CONVERGED_ITERATING) {
 47:     eps->its++;

 49:     /* Compute an nv-step Lanczos factorization */
 50:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
 51:     DSGetArrayReal(eps->ds,DS_MAT_T,&a);
 52:     b = a + ld;
 53:     DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
 54:     EPSPseudoLanczos(eps,a,b,omega,eps->nconv+l,&nv,&breakdown,&symmlost,NULL,w);
 55:     if (symmlost) {
 56:       eps->reason = EPS_DIVERGED_SYMMETRY_LOST;
 57:       if (nv==eps->nconv+l+1) { eps->nconv = nconv; break; }
 58:     }
 59:     beta = b[nv-1];
 60:     DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
 61:     DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
 62:     DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l);
 63:     if (l==0) DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
 64:     else DSSetState(eps->ds,DS_STATE_RAW);
 65:     BVSetActiveColumns(eps->V,eps->nconv,nv);

 67:     /* Solve projected problem */
 68:     DSSolve(eps->ds,eps->eigr,eps->eigi);
 69:     DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
 70:     DSUpdateExtraRow(eps->ds);
 71:     DSSynchronize(eps->ds,eps->eigr,eps->eigi);

 73:     /* Check convergence */
 74:     DSGetDimensions(eps->ds,NULL,NULL,NULL,&t);
 75: #if 0
 76:     /* take into account also left residual */
 77:     BVGetColumn(eps->V,nv,&u);
 78:     VecNorm(u,NORM_2,&beta1);
 79:     BVRestoreColumn(eps->V,nv,&u);
 80:     VecNorm(w,NORM_2,&beta2);  /* w contains B*V[nv] */
 81:     beta1 = PetscMax(beta1,beta2);
 82: #endif
 83:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,t-eps->nconv,beta*beta1,0.0,1.0,&k);
 84:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
 85:     nconv = k;

 87:     /* Update l */
 88:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
 89:     else {
 90:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
 91:       l = PetscMin(l,t);
 92:       DSGetTruncateSize(eps->ds,k,t,&l);
 93:     }
 94:     if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
 95:     if (l) PetscInfo(eps,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l);

 97:     if (eps->reason == EPS_CONVERGED_ITERATING) {
 99:       /* Prepare the Rayleigh quotient for restart */
100:       DSTruncate(eps->ds,k+l,PETSC_FALSE);
101:     }
102:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
103:     DSGetMat(eps->ds,DS_MAT_Q,&U);
104:     BVMultInPlace(eps->V,U,eps->nconv,k+l);
105:     MatDestroy(&U);

107:     /* Move restart vector and update signature */
108:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
109:       BVCopyColumn(eps->V,nv,k+l);
110:       DSGetArrayReal(eps->ds,DS_MAT_D,&omega);
111:       VecCreateSeq(PETSC_COMM_SELF,k+l,&vomega);
112:       VecGetArray(vomega,&aux);
113:       for (i=0;i<k+l;i++) aux[i] = omega[i];
114:       VecRestoreArray(vomega,&aux);
115:       BVSetActiveColumns(eps->V,0,k+l);
116:       BVSetSignature(eps->V,vomega);
117:       VecDestroy(&vomega);
118:       DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega);
119:     }

121:     eps->nconv = k;
122:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv);
123:   }
124:   DSTruncate(eps->ds,eps->nconv,PETSC_TRUE);
125:   PetscFunctionReturn(0);
126: }