Actual source code: ks-indef.c

slepc-3.21.1 2024-04-26
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        k,l,ld,nv,t,nconv=0;
 22:   Mat             U,D;
 23:   Vec             vomega,w=eps->work[0];
 24:   PetscReal       *a,*b,beta,beta1=1.0,*omega;
 25:   PetscBool       breakdown=PETSC_FALSE,symmlost=PETSC_FALSE;

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

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

 33:   /* Extract sigma[0] from BV, computed during normalization */
 34:   PetscCall(DSSetDimensions(eps->ds,1,PETSC_DEFAULT,PETSC_DEFAULT));
 35:   PetscCall(BVSetActiveColumns(eps->V,0,1));
 36:   PetscCall(DSGetMatAndColumn(eps->ds,DS_MAT_D,0,&D,&vomega));
 37:   PetscCall(BVGetSignature(eps->V,vomega));
 38:   PetscCall(DSRestoreMatAndColumn(eps->ds,DS_MAT_D,0,&D,&vomega));
 39:   l = 0;

 41:   /* Restart loop */
 42:   while (eps->reason == EPS_CONVERGED_ITERATING) {
 43:     eps->its++;

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

 63:     /* Solve projected problem */
 64:     PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi));
 65:     PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL));
 66:     PetscCall(DSUpdateExtraRow(eps->ds));
 67:     PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));

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

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

 93:     if (eps->reason == EPS_CONVERGED_ITERATING) {
 94:       PetscCheck(!breakdown,PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Breakdown in Indefinite Krylov-Schur (beta=%g)",(double)beta);
 95:       /* Prepare the Rayleigh quotient for restart */
 96:       PetscCall(DSTruncate(eps->ds,k+l,PETSC_FALSE));
 97:     }
 98:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
 99:     PetscCall(DSGetMat(eps->ds,DS_MAT_Q,&U));
100:     PetscCall(BVMultInPlace(eps->V,U,eps->nconv,k+l));
101:     PetscCall(DSRestoreMat(eps->ds,DS_MAT_Q,&U));

103:     /* Move restart vector and update signature */
104:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
105:       PetscCall(BVCopyColumn(eps->V,nv,k+l));
106:       PetscCall(BVSetActiveColumns(eps->V,0,k+l));
107:       PetscCall(DSGetMatAndColumn(eps->ds,DS_MAT_D,0,&D,&vomega));
108:       PetscCall(BVSetSignature(eps->V,vomega));
109:       PetscCall(DSRestoreMatAndColumn(eps->ds,DS_MAT_D,0,&D,&vomega));
110:     }

112:     eps->nconv = k;
113:     PetscCall(EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv));
114:   }
115:   PetscCall(DSTruncate(eps->ds,eps->nconv,PETSC_TRUE));
116:   PetscFunctionReturn(PETSC_SUCCESS);
117: }