Actual source code: mfnexpokit.c

slepc-3.16.0 2021-09-30
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: "expokit"

 13:    Method: Arnoldi method tailored for the matrix exponential

 15:    Algorithm:

 17:        Uses Arnoldi relations to compute exp(t_step*A)*v_last for
 18:        several time steps.

 20:    References:

 22:        [1] R. Sidje, "Expokit: a software package for computing matrix
 23:            exponentials", ACM Trans. Math. Softw. 24(1):130-156, 1998.
 24: */

 26: #include <slepc/private/mfnimpl.h>

 28: PetscErrorCode MFNSetUp_Expokit(MFN mfn)
 29: {
 31:   PetscInt       N;
 32:   PetscBool      isexp;

 35:   MatGetSize(mfn->A,&N,NULL);
 36:   if (mfn->ncv==PETSC_DEFAULT) mfn->ncv = PetscMin(30,N);
 37:   if (mfn->max_it==PETSC_DEFAULT) mfn->max_it = 100;
 38:   MFNAllocateSolution(mfn,2);

 40:   PetscObjectTypeCompare((PetscObject)mfn->fn,FNEXP,&isexp);
 41:   if (!isexp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This solver only supports the exponential function");
 42:   return(0);
 43: }

 45: PetscErrorCode MFNSolve_Expokit(MFN mfn,Vec b,Vec x)
 46: {
 47:   PetscErrorCode    ierr;
 48:   PetscInt          mxstep,mxrej,m,mb,ld,i,j,ireject,mx,k1;
 49:   Vec               v,r;
 50:   Mat               H,M=NULL,K=NULL;
 51:   FN                fn;
 52:   PetscScalar       *Harray,*B,*F,*betaF,t,sgn,sfactor;
 53:   const PetscScalar *pK;
 54:   PetscReal         anorm,avnorm,tol,err_loc,rndoff,t_out,t_new,t_now,t_step;
 55:   PetscReal         xm,fact,s,p1,p2,beta,beta2,gamma,delta;
 56:   PetscBool         breakdown;

 59:   m   = mfn->ncv;
 60:   tol = mfn->tol;
 61:   mxstep = mfn->max_it;
 62:   mxrej = 10;
 63:   gamma = 0.9;
 64:   delta = 1.2;
 65:   mb    = m;
 66:   FNGetScale(mfn->fn,&t,&sfactor);
 67:   FNDuplicate(mfn->fn,PetscObjectComm((PetscObject)mfn->fn),&fn);
 68:   FNSetScale(fn,1.0,1.0);
 69:   t_out = PetscAbsScalar(t);
 70:   t_now = 0.0;
 71:   MatNorm(mfn->A,NORM_INFINITY,&anorm);
 72:   rndoff = anorm*PETSC_MACHINE_EPSILON;

 74:   k1 = 2;
 75:   xm = 1.0/(PetscReal)m;
 76:   beta = mfn->bnorm;
 77:   fact = PetscPowRealInt((m+1)/2.72,m+1)*PetscSqrtReal(2.0*PETSC_PI*(m+1));
 78:   t_new = (1.0/anorm)*PetscPowReal((fact*tol)/(4.0*beta*anorm),xm);
 79:   s = PetscPowReal(10.0,PetscFloorReal(PetscLog10Real(t_new))-1);
 80:   t_new = PetscCeilReal(t_new/s)*s;
 81:   sgn = t/PetscAbsScalar(t);

 83:   VecCopy(b,x);
 84:   ld = m+2;
 85:   PetscCalloc2(m+1,&betaF,ld*ld,&B);
 86:   MatCreateSeqDense(PETSC_COMM_SELF,ld,ld,NULL,&H);
 87:   MatDenseGetArray(H,&Harray);

 89:   while (mfn->reason == MFN_CONVERGED_ITERATING) {
 90:     mfn->its++;
 91:     if (PetscIsInfOrNanReal(t_new)) t_new = PETSC_MAX_REAL;
 92:     t_step = PetscMin(t_out-t_now,t_new);
 93:     BVInsertVec(mfn->V,0,x);
 94:     BVScaleColumn(mfn->V,0,1.0/beta);
 95:     BVMatArnoldi(mfn->V,mfn->transpose_solve?mfn->AT:mfn->A,H,0,&mb,&beta2,&breakdown);
 96:     if (breakdown) {
 97:       k1 = 0;
 98:       t_step = t_out-t_now;
 99:     }
100:     if (k1!=0) {
101:       Harray[m+1+ld*m] = 1.0;
102:       BVGetColumn(mfn->V,m,&v);
103:       BVGetColumn(mfn->V,m+1,&r);
104:       MatMult(mfn->transpose_solve?mfn->AT:mfn->A,v,r);
105:       BVRestoreColumn(mfn->V,m,&v);
106:       BVRestoreColumn(mfn->V,m+1,&r);
107:       BVNormColumn(mfn->V,m+1,NORM_2,&avnorm);
108:     }
109:     PetscArraycpy(B,Harray,ld*ld);

111:     ireject = 0;
112:     while (ireject <= mxrej) {
113:       mx = mb + k1;
114:       for (i=0;i<mx;i++) {
115:         for (j=0;j<mx;j++) {
116:           Harray[i+j*ld] = sgn*B[i+j*ld]*t_step;
117:         }
118:       }
119:       MFN_CreateDenseMat(mx,&M);
120:       MFN_CreateDenseMat(mx,&K);
121:       MatDenseGetArray(M,&F);
122:       for (j=0;j<mx;j++) {
123:         PetscArraycpy(F+j*mx,Harray+j*ld,mx);
124:       }
125:       MatDenseRestoreArray(M,&F);
126:       FNEvaluateFunctionMat(fn,M,K);

128:       if (k1==0) {
129:         err_loc = tol;
130:         break;
131:       } else {
132:         MatDenseGetArrayRead(K,&pK);
133:         p1 = PetscAbsScalar(beta*pK[m]);
134:         p2 = PetscAbsScalar(beta*pK[m+1]*avnorm);
135:         MatDenseRestoreArrayRead(K,&pK);
136:         if (p1 > 10*p2) {
137:           err_loc = p2;
138:           xm = 1.0/(PetscReal)m;
139:         } else if (p1 > p2) {
140:           err_loc = (p1*p2)/(p1-p2);
141:           xm = 1.0/(PetscReal)m;
142:         } else {
143:           err_loc = p1;
144:           xm = 1.0/(PetscReal)(m-1);
145:         }
146:       }
147:       if (err_loc <= delta*t_step*tol) break;
148:       else {
149:         t_step = gamma*t_step*PetscPowReal(t_step*tol/err_loc,xm);
150:         s = PetscPowReal(10.0,PetscFloorReal(PetscLog10Real(t_step))-1);
151:         t_step = PetscCeilReal(t_step/s)*s;
152:         ireject = ireject+1;
153:       }
154:     }

156:     mx = mb + PetscMax(0,k1-1);
157:     MatDenseGetArrayRead(K,&pK);
158:     for (j=0;j<mx;j++) betaF[j] = beta*pK[j];
159:     MatDenseRestoreArrayRead(K,&pK);
160:     BVSetActiveColumns(mfn->V,0,mx);
161:     BVMultVec(mfn->V,1.0,0.0,x,betaF);
162:     VecNorm(x,NORM_2,&beta);

164:     t_now = t_now+t_step;
165:     if (t_now>=t_out) mfn->reason = MFN_CONVERGED_TOL;
166:     else {
167:       t_new = gamma*t_step*PetscPowReal((t_step*tol)/err_loc,xm);
168:       s = PetscPowReal(10.0,PetscFloorReal(PetscLog10Real(t_new))-1);
169:       t_new = PetscCeilReal(t_new/s)*s;
170:     }
171:     err_loc = PetscMax(err_loc,rndoff);
172:     if (mfn->its==mxstep) mfn->reason = MFN_DIVERGED_ITS;
173:     MFNMonitor(mfn,mfn->its,err_loc);
174:   }
175:   VecScale(x,sfactor);

177:   MatDestroy(&M);
178:   MatDestroy(&K);
179:   FNDestroy(&fn);
180:   MatDenseRestoreArray(H,&Harray);
181:   MatDestroy(&H);
182:   PetscFree2(betaF,B);
183:   return(0);
184: }

186: SLEPC_EXTERN PetscErrorCode MFNCreate_Expokit(MFN mfn)
187: {
189:   mfn->ops->solve          = MFNSolve_Expokit;
190:   mfn->ops->setup          = MFNSetUp_Expokit;
191:   return(0);
192: }