Actual source code: mfnkrylov.c

slepc-3.15.1 2021-05-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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 matrix function solver: "krylov"

 13:    Method: Arnoldi with Eiermann-Ernst restart

 15:    Algorithm:

 17:        Build Arnoldi approximations using f(H) for the Hessenberg matrix H,
 18:        restart by discarding the Krylov basis but keeping H.

 20:    References:

 22:        [1] M. Eiermann and O. Ernst, "A restarted Krylov subspace method
 23:            for the evaluation of matrix functions", SIAM J. Numer. Anal.
 24:            44(6):2481-2504, 2006.
 25: */

 27: #include <slepc/private/mfnimpl.h>
 28: #include <slepcblaslapack.h>

 30: PetscErrorCode MFNSetUp_Krylov(MFN mfn)
 31: {
 33:   PetscInt       N;

 36:   MatGetSize(mfn->A,&N,NULL);
 37:   if (mfn->ncv==PETSC_DEFAULT) mfn->ncv = PetscMin(30,N);
 38:   if (mfn->max_it==PETSC_DEFAULT) mfn->max_it = 100;
 39:   MFNAllocateSolution(mfn,1);
 40:   return(0);
 41: }

 43: PetscErrorCode MFNSolve_Krylov(MFN mfn,Vec b,Vec x)
 44: {
 45:   PetscErrorCode    ierr;
 46:   PetscInt          n=0,m,ld,ldh,j;
 47:   PetscBLASInt      m_,inc=1;
 48:   Mat               G=NULL,H=NULL;
 49:   Vec               F=NULL;
 50:   PetscScalar       *array,*farray,*harray;
 51:   const PetscScalar *garray;
 52:   PetscReal         beta,betaold=0.0,nrm=1.0;
 53:   PetscBool         breakdown;

 56:   m  = mfn->ncv;
 57:   ld = m+1;
 58:   PetscCalloc1(ld*ld,&array);

 60:   /* set initial vector to b/||b|| */
 61:   BVInsertVec(mfn->V,0,b);
 62:   BVScaleColumn(mfn->V,0,1.0/mfn->bnorm);
 63:   VecSet(x,0.0);

 65:   /* Restart loop */
 66:   while (mfn->reason == MFN_CONVERGED_ITERATING) {
 67:     mfn->its++;

 69:     /* compute Arnoldi factorization */
 70:     BVMatArnoldi(mfn->V,mfn->transpose_solve?mfn->AT:mfn->A,array,ld,0,&m,&beta,&breakdown);

 72:     /* save previous Hessenberg matrix in G; allocate new storage for H and f(H) */
 73:     if (mfn->its>1) { G = H; H = NULL; }
 74:     ldh = n+m;
 75:     MFN_CreateVec(ldh,&F);
 76:     MFN_CreateDenseMat(ldh,&H);

 78:     /* glue together the previous H and the new H obtained with Arnoldi */
 79:     MatDenseGetArray(H,&harray);
 80:     for (j=0;j<m;j++) {
 81:       PetscArraycpy(harray+n+(j+n)*ldh,array+j*ld,m);
 82:     }
 83:     if (mfn->its>1) {
 84:       MatDenseGetArrayRead(G,&garray);
 85:       for (j=0;j<n;j++) {
 86:         PetscArraycpy(harray+j*ldh,garray+j*n,n);
 87:       }
 88:       MatDenseRestoreArrayRead(G,&garray);
 89:       MatDestroy(&G);
 90:       harray[n+(n-1)*ldh] = betaold;
 91:     }
 92:     MatDenseRestoreArray(H,&harray);

 94:     if (mfn->its==1) {
 95:       /* set symmetry flag of H from A */
 96:       MatPropagateSymmetryOptions(mfn->A,H);
 97:     }

 99:     /* evaluate f(H) */
100:     FNEvaluateFunctionMatVec(mfn->fn,H,F);

102:     /* x += ||b||*V*f(H)*e_1 */
103:     VecGetArray(F,&farray);
104:     PetscBLASIntCast(m,&m_);
105:     nrm = BLASnrm2_(&m_,farray+n,&inc);   /* relative norm of the update ||u||/||b|| */
106:     MFNMonitor(mfn,mfn->its,nrm);
107:     for (j=0;j<m;j++) farray[j+n] *= mfn->bnorm;
108:     BVSetActiveColumns(mfn->V,0,m);
109:     BVMultVec(mfn->V,1.0,1.0,x,farray+n);
110:     VecRestoreArray(F,&farray);

112:     /* check convergence */
113:     if (mfn->its >= mfn->max_it) mfn->reason = MFN_DIVERGED_ITS;
114:     if (mfn->its>1) {
115:       if (m<mfn->ncv || breakdown || beta==0.0 || nrm<mfn->tol) mfn->reason = MFN_CONVERGED_TOL;
116:     }

118:     /* restart with vector v_{m+1} */
119:     if (mfn->reason == MFN_CONVERGED_ITERATING) {
120:       BVCopyColumn(mfn->V,m,0);
121:       n += m;
122:       betaold = beta;
123:     }
124:   }

126:   MatDestroy(&H);
127:   MatDestroy(&G);
128:   VecDestroy(&F);
129:   PetscFree(array);
130:   return(0);
131: }

133: SLEPC_EXTERN PetscErrorCode MFNCreate_Krylov(MFN mfn)
134: {
136:   mfn->ops->solve          = MFNSolve_Krylov;
137:   mfn->ops->setup          = MFNSetUp_Krylov;
138:   return(0);
139: }