Actual source code: mfnkrylov.c

  1: /*

  3:    SLEPc matrix function solver: "krylov"

  5:    Method: Arnoldi

  7:    Algorithm:

  9:        Single-vector Arnoldi method to build a Krylov subspace, then
 10:        compute f(B) on the projected matrix B.

 12:    References:

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

 17:    Last update: Feb 2013

 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 21:    Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain

 23:    This file is part of SLEPc.

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

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

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

 39: #include <slepc-private/mfnimpl.h>                /*I "slepcmfn.h" I*/

 43: PetscErrorCode MFNSetUp_Krylov(MFN mfn)
 44: {
 45:   PetscErrorCode  ierr;

 48:   if (!mfn->ncv) mfn->ncv = PetscMin(30,mfn->n);
 49:   if (!mfn->max_it) mfn->max_it = PetscMax(100,2*mfn->n/mfn->ncv);
 50:   VecDuplicateVecs(mfn->t,mfn->ncv+1,&mfn->V);
 51:   PetscLogObjectParents(mfn,mfn->ncv+1,mfn->V);
 52:   mfn->allocated_ncv = mfn->ncv+1;
 53:   DSAllocate(mfn->ds,mfn->ncv+2);
 54:   DSSetType(mfn->ds,DSNHEP);
 55:   return(0);
 56: }

 60: static PetscErrorCode MFNBasicArnoldi(MFN mfn,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscReal *beta,PetscBool *breakdown)
 61: {
 63:   PetscInt       j,m = *M;
 64:   PetscReal      norm;

 67:   for (j=k;j<m-1;j++) {
 68:     MatMult(mfn->A,V[j],V[j+1]);
 69:     IPOrthogonalize(mfn->ip,0,NULL,j+1,NULL,V,V[j+1],H+ldh*j,&norm,breakdown);
 70:     H[j+1+ldh*j] = norm;
 71:     if (*breakdown) {
 72:       *M = j+1;
 73:       *beta = norm;
 74:       return(0);
 75:     } else {
 76:       VecScale(V[j+1],1/norm);
 77:     }
 78:   }
 79:   MatMult(mfn->A,V[m-1],f);
 80:   IPOrthogonalize(mfn->ip,0,NULL,m,NULL,V,f,H+ldh*(m-1),beta,NULL);
 81:   return(0);
 82: }

 86: PetscErrorCode MFNSolve_Krylov(MFN mfn,Vec b,Vec x)
 87: {
 89:   PetscInt       mxstep,mxrej,m,mb,ld,i,j,ireject,mx,k1;
 90:   Vec            r;
 91:   PetscScalar    *H,*B,*F,*betaF;
 92:   PetscReal      anorm,normb,avnorm,tol,err_loc,rndoff;
 93:   PetscReal      t,t_out,t_new,t_now,t_step;
 94:   PetscReal      xm,fact,s,sgn,p1,p2;
 95:   PetscReal      beta,beta2,gamma,delta;
 96:   PetscBool      breakdown;

 99:   m   = mfn->ncv;
100:   tol = mfn->tol;
101:   mxstep = mfn->max_it;
102:   mxrej = 10;
103:   gamma = 0.9;
104:   delta = 1.2;
105:   mb    = m;
106:   t     = PetscRealPart(mfn->sfactor);
107:   t_out = PetscAbsReal(t);
108:   t_new = 0.0;
109:   t_now = 0.0;
110:   MatNorm(mfn->A,NORM_INFINITY,&anorm);
111:   rndoff = anorm*PETSC_MACHINE_EPSILON;

113:   k1 = 2;
114:   xm = 1.0/(PetscReal)m;
115:   VecNorm(b,NORM_2,&normb);
116:   beta = normb;
117:   fact = PetscPowRealInt((m+1)/2.72,m+1)*PetscSqrtReal(2*PETSC_PI*(m+1));
118:   t_new = (1.0/anorm)*PetscPowReal((fact*tol)/(4.0*beta*anorm),xm);
119:   s = PetscPowReal(10,floor(log10(t_new))-1);
120:   t_new = ceil(t_new/s)*s;
121:   sgn = PetscSign(t);

123:   PetscMalloc((m+1)*sizeof(PetscScalar),&betaF);
124:   VecDuplicate(mfn->V[0],&r);
125:   VecCopy(b,x);
126:   DSGetLeadingDimension(mfn->ds,&ld);
127:   PetscMalloc(ld*ld*sizeof(PetscScalar),&B);

129:   while (mfn->reason == MFN_CONVERGED_ITERATING) {
130:     mfn->its++;
131:     if (PetscIsInfOrNanReal(t_new)) t_new = PETSC_MAX_REAL;
132:     t_step = PetscMin(t_out-t_now,t_new);

134:     VecCopy(x,mfn->V[0]);
135:     VecScale(mfn->V[0],1.0/beta);
136:     DSGetArray(mfn->ds,DS_MAT_A,&H);
137:     MFNBasicArnoldi(mfn,H,ld,mfn->V,0,&mb,r,&beta2,&breakdown);
138:     H[mb+(mb-1)*ld] = beta2;
139:     VecScale(r,1.0/beta2);
140:     VecCopy(r,mfn->V[m]);
141:     if (breakdown) {
142:       k1 = 0;
143:       t_step = t_out-t_now;
144:     }
145:     if (k1!=0) {
146:       H[m+1+ld*m] = 1.0;
147:       MatMult(mfn->A,mfn->V[m],r);
148:       VecNorm(r,NORM_2,&avnorm);
149:     }
150:     PetscMemcpy(B,H,ld*ld*sizeof(PetscScalar));
151:     DSRestoreArray(mfn->ds,DS_MAT_A,&H);

153:     ireject = 0;
154:     while (ireject <= mxrej) {
155:       mx = mb + k1;
156:       DSSetDimensions(mfn->ds,mx,0,0,0);
157:       DSGetArray(mfn->ds,DS_MAT_A,&H);
158:       for (i=0;i<mx;i++) {
159:         for (j=0;j<mx;j++) {
160:           H[i+j*ld] = sgn*B[i+j*ld]*t_step;
161:         }
162:       }
163:       DSRestoreArray(mfn->ds,DS_MAT_A,&H);
164:       DSComputeFunction(mfn->ds,mfn->function);

166:       if (k1==0) {
167:         err_loc = tol;
168:         break;
169:       } else {
170:         DSGetArray(mfn->ds,DS_MAT_F,&F);
171:         p1 = PetscAbsScalar(beta*F[m]);
172:         p2 = PetscAbsScalar(beta*F[m+1]*avnorm);
173:         DSRestoreArray(mfn->ds,DS_MAT_F,&F);
174:         if (p1 > 10*p2) {
175:           err_loc = p2;
176:           xm = 1.0/(PetscReal)m;
177:         } else if (p1 > p2) {
178:           err_loc = (p1*p2)/(p1-p2);
179:           xm = 1.0/(PetscReal)m;
180:         } else {
181:           err_loc = p1;
182:           xm = 1.0/(PetscReal)(m-1);
183:         }
184:       }
185:       if (err_loc <= delta*t_step*tol) break;
186:       else {
187:         t_step = gamma*t_step*PetscPowReal(t_step*tol/err_loc,xm);
188:         s = PetscPowReal(10,floor(log10(t_step))-1);
189:         t_step = ceil(t_step/s)*s;
190:         ireject = ireject+1;
191:       }
192:     }

194:     mx = mb + PetscMax(0,k1-1);
195:     DSGetArray(mfn->ds,DS_MAT_F,&F);
196:     for (j=0;j<mx;j++) betaF[j] = beta*F[j];
197:     DSRestoreArray(mfn->ds,DS_MAT_F,&F);
198:     VecSet(x,0.0);
199:     VecMAXPY(x,mx,betaF,mfn->V);
200:     VecNorm(x,NORM_2,&beta);

202:     t_now = t_now+t_step;
203:     if (t_now>=t_out) mfn->reason = MFN_CONVERGED_TOL;
204:     else {
205:       t_new = gamma*t_step*PetscPowReal((t_step*tol)/err_loc,xm);
206:       s = PetscPowReal(10,floor(log10(t_new))-1);
207:       t_new = ceil(t_new/s)*s;
208:     }
209:     err_loc = PetscMax(err_loc,rndoff);
210:     if (mfn->its==mxstep) mfn->reason = MFN_DIVERGED_ITS;
211:     MFNMonitor(mfn,mfn->its,t_now);
212:   }

214:   VecDestroy(&r);
215:   PetscFree(betaF);
216:   PetscFree(B);
217:   return(0);
218: }

222: PetscErrorCode MFNReset_Krylov(MFN mfn)
223: {
224:   PetscErrorCode  ierr;

227:   if (mfn->allocated_ncv > 0) {
228:     VecDestroyVecs(mfn->allocated_ncv,&mfn->V);
229:     mfn->allocated_ncv = 0;
230:   }
231:   return(0);
232: }

236: PETSC_EXTERN PetscErrorCode MFNCreate_Krylov(MFN mfn)
237: {
239:   mfn->ops->solve          = MFNSolve_Krylov;
240:   mfn->ops->setup          = MFNSetUp_Krylov;
241:   mfn->ops->reset          = MFNReset_Krylov;
242:   return(0);
243: }