Actual source code: mfnexpokit.c
slepc-3.7.0 2016-05-16
1: /*
3: SLEPc matrix function solver: "expokit"
5: Method: Arnoldi method tailored for the matrix exponential
7: Algorithm:
9: Uses Arnoldi relations to compute exp(t_step*A)*v_last for
10: several time steps.
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: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
18: SLEPc - Scalable Library for Eigenvalue Problem Computations
19: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
21: This file is part of SLEPc.
23: SLEPc is free software: you can redistribute it and/or modify it under the
24: terms of version 3 of the GNU Lesser General Public License as published by
25: the Free Software Foundation.
27: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
28: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
29: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
30: more details.
32: You should have received a copy of the GNU Lesser General Public License
33: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
34: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35: */
37: #include <slepc/private/mfnimpl.h>
41: PetscErrorCode MFNSetUp_Expokit(MFN mfn)
42: {
44: PetscInt N;
45: PetscBool isexp;
48: MatGetSize(mfn->A,&N,NULL);
49: if (!mfn->ncv) mfn->ncv = PetscMin(30,N);
50: if (!mfn->max_it) mfn->max_it = 100;
51: MFNAllocateSolution(mfn,2);
53: PetscObjectTypeCompare((PetscObject)mfn->fn,FNEXP,&isexp);
54: if (!isexp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This solver only supports the exponential function");
55: return(0);
56: }
60: PetscErrorCode MFNSolve_Expokit(MFN mfn,Vec b,Vec x)
61: {
63: PetscInt mxstep,mxrej,m,mb,ld,i,j,ireject,mx,k1;
64: Vec v,r;
65: Mat M=NULL,K=NULL;
66: FN fn;
67: PetscScalar *H,*B,*F,*betaF,t,sgn,sfactor;
68: PetscReal anorm,avnorm,tol,err_loc,rndoff;
69: PetscReal t_out,t_new,t_now,t_step;
70: PetscReal xm,fact,s,p1,p2;
71: PetscReal beta,beta2,gamma,delta;
72: PetscBool breakdown;
75: m = mfn->ncv;
76: tol = mfn->tol;
77: mxstep = mfn->max_it;
78: mxrej = 10;
79: gamma = 0.9;
80: delta = 1.2;
81: mb = m;
82: FNGetScale(mfn->fn,&t,&sfactor);
83: FNDuplicate(mfn->fn,PetscObjectComm((PetscObject)mfn->fn),&fn);
84: FNSetScale(fn,1.0,1.0);
85: t_out = PetscAbsScalar(t);
86: t_now = 0.0;
87: MatNorm(mfn->A,NORM_INFINITY,&anorm);
88: rndoff = anorm*PETSC_MACHINE_EPSILON;
90: k1 = 2;
91: xm = 1.0/(PetscReal)m;
92: beta = mfn->bnorm;
93: fact = PetscPowRealInt((m+1)/2.72,m+1)*PetscSqrtReal(2*PETSC_PI*(m+1));
94: t_new = (1.0/anorm)*PetscPowReal((fact*tol)/(4.0*beta*anorm),xm);
95: s = PetscPowReal(10.0,PetscFloorReal(PetscLog10Real(t_new))-1);
96: t_new = PetscCeilReal(t_new/s)*s;
97: sgn = t/PetscAbsScalar(t);
99: VecCopy(b,x);
100: ld = m+2;
101: PetscCalloc3(m+1,&betaF,ld*ld,&H,ld*ld,&B);
103: while (mfn->reason == MFN_CONVERGED_ITERATING) {
104: mfn->its++;
105: if (PetscIsInfOrNanReal(t_new)) t_new = PETSC_MAX_REAL;
106: t_step = PetscMin(t_out-t_now,t_new);
107: BVInsertVec(mfn->V,0,x);
108: BVScaleColumn(mfn->V,0,1.0/beta);
109: MFNBasicArnoldi(mfn,H,ld,0,&mb,&beta2,&breakdown);
110: if (breakdown) {
111: k1 = 0;
112: t_step = t_out-t_now;
113: }
114: if (k1!=0) {
115: H[m+1+ld*m] = 1.0;
116: BVGetColumn(mfn->V,m,&v);
117: BVGetColumn(mfn->V,m+1,&r);
118: MatMult(mfn->A,v,r);
119: BVRestoreColumn(mfn->V,m,&v);
120: BVRestoreColumn(mfn->V,m+1,&r);
121: BVNormColumn(mfn->V,m+1,NORM_2,&avnorm);
122: }
123: PetscMemcpy(B,H,ld*ld*sizeof(PetscScalar));
125: ireject = 0;
126: while (ireject <= mxrej) {
127: mx = mb + k1;
128: for (i=0;i<mx;i++) {
129: for (j=0;j<mx;j++) {
130: H[i+j*ld] = sgn*B[i+j*ld]*t_step;
131: }
132: }
133: MFN_CreateDenseMat(mx,&M);
134: MFN_CreateDenseMat(mx,&K);
135: MatDenseGetArray(M,&F);
136: for (j=0;j<mx;j++) {
137: PetscMemcpy(F+j*mx,H+j*ld,mx*sizeof(PetscScalar));
138: }
139: MatDenseRestoreArray(M,&F);
140: FNEvaluateFunctionMat(fn,M,K);
142: if (k1==0) {
143: err_loc = tol;
144: break;
145: } else {
146: MatDenseGetArray(K,&F);
147: p1 = PetscAbsScalar(beta*F[m]);
148: p2 = PetscAbsScalar(beta*F[m+1]*avnorm);
149: MatDenseRestoreArray(K,&F);
150: if (p1 > 10*p2) {
151: err_loc = p2;
152: xm = 1.0/(PetscReal)m;
153: } else if (p1 > p2) {
154: err_loc = (p1*p2)/(p1-p2);
155: xm = 1.0/(PetscReal)m;
156: } else {
157: err_loc = p1;
158: xm = 1.0/(PetscReal)(m-1);
159: }
160: }
161: if (err_loc <= delta*t_step*tol) break;
162: else {
163: t_step = gamma*t_step*PetscPowReal(t_step*tol/err_loc,xm);
164: s = PetscPowReal(10.0,PetscFloorReal(PetscLog10Real(t_step))-1);
165: t_step = PetscCeilReal(t_step/s)*s;
166: ireject = ireject+1;
167: }
168: }
170: mx = mb + PetscMax(0,k1-1);
171: MatDenseGetArray(K,&F);
172: for (j=0;j<mx;j++) betaF[j] = beta*F[j];
173: MatDenseRestoreArray(K,&F);
174: BVSetActiveColumns(mfn->V,0,mx);
175: BVMultVec(mfn->V,1.0,0.0,x,betaF);
176: VecNorm(x,NORM_2,&beta);
178: t_now = t_now+t_step;
179: if (t_now>=t_out) mfn->reason = MFN_CONVERGED_TOL;
180: else {
181: t_new = gamma*t_step*PetscPowReal((t_step*tol)/err_loc,xm);
182: s = PetscPowReal(10.0,PetscFloorReal(PetscLog10Real(t_new))-1);
183: t_new = PetscCeilReal(t_new/s)*s;
184: }
185: err_loc = PetscMax(err_loc,rndoff);
186: if (mfn->its==mxstep) mfn->reason = MFN_DIVERGED_ITS;
187: MFNMonitor(mfn,mfn->its,err_loc);
188: }
189: VecScale(x,sfactor);
191: MatDestroy(&M);
192: MatDestroy(&K);
193: FNDestroy(&fn);
194: PetscFree3(betaF,H,B);
195: return(0);
196: }
200: PETSC_EXTERN PetscErrorCode MFNCreate_Expokit(MFN mfn)
201: {
203: mfn->ops->solve = MFNSolve_Expokit;
204: mfn->ops->setup = MFNSetUp_Expokit;
205: return(0);
206: }