Actual source code: fnexp.c
slepc-3.12.1 2019-11-08
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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: Exponential function exp(x)
12: */
14: #include <slepc/private/fnimpl.h> /*I "slepcfn.h" I*/
15: #include <slepcblaslapack.h>
17: PetscErrorCode FNEvaluateFunction_Exp(FN fn,PetscScalar x,PetscScalar *y)
18: {
20: *y = PetscExpScalar(x);
21: return(0);
22: }
24: PetscErrorCode FNEvaluateDerivative_Exp(FN fn,PetscScalar x,PetscScalar *y)
25: {
27: *y = PetscExpScalar(x);
28: return(0);
29: }
31: #define MAX_PADE 6
32: #define SWAP(a,b,t) {t=a;a=b;b=t;}
34: PetscErrorCode FNEvaluateFunctionMat_Exp_Pade(FN fn,Mat A,Mat B)
35: {
36: #if defined(PETSC_MISSING_LAPACK_GESV) || defined(SLEPC_MISSING_LAPACK_LANGE)
38: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESV/LANGE - Lapack routines are unavailable");
39: #else
41: PetscBLASInt n,ld,ld2,*ipiv,info,inc=1;
42: PetscInt m,j,k,sexp;
43: PetscBool odd;
44: const PetscInt p=MAX_PADE;
45: PetscReal c[MAX_PADE+1],s,*rwork;
46: PetscScalar scale,mone=-1.0,one=1.0,two=2.0,zero=0.0;
47: PetscScalar *Aa,*Ba,*As,*A2,*Q,*P,*W,*aux;
50: MatDenseGetArray(A,&Aa);
51: MatDenseGetArray(B,&Ba);
52: MatGetSize(A,&m,NULL);
53: PetscBLASIntCast(m,&n);
54: ld = n;
55: ld2 = ld*ld;
56: P = Ba;
57: PetscMalloc6(m*m,&Q,m*m,&W,m*m,&As,m*m,&A2,ld,&rwork,ld,&ipiv);
58: PetscArraycpy(As,Aa,ld2);
60: /* Pade' coefficients */
61: c[0] = 1.0;
62: for (k=1;k<=p;k++) c[k] = c[k-1]*(p+1-k)/(k*(2*p+1-k));
64: /* Scaling */
65: s = LAPACKlange_("I",&n,&n,As,&ld,rwork);
66: PetscLogFlops(1.0*n*n);
67: if (s>0.5) {
68: sexp = PetscMax(0,(int)(PetscLogReal(s)/PetscLogReal(2.0))+2);
69: scale = PetscPowRealInt(2.0,-sexp);
70: PetscStackCallBLAS("BLASscal",BLASscal_(&ld2,&scale,As,&inc));
71: PetscLogFlops(1.0*n*n);
72: } else sexp = 0;
74: /* Horner evaluation */
75: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,As,&ld,As,&ld,&zero,A2,&ld));
76: PetscLogFlops(2.0*n*n*n);
77: PetscArrayzero(Q,ld2);
78: PetscArrayzero(P,ld2);
79: for (j=0;j<n;j++) {
80: Q[j+j*ld] = c[p];
81: P[j+j*ld] = c[p-1];
82: }
84: odd = PETSC_TRUE;
85: for (k=p-1;k>0;k--) {
86: if (odd) {
87: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,Q,&ld,A2,&ld,&zero,W,&ld));
88: SWAP(Q,W,aux);
89: for (j=0;j<n;j++) Q[j+j*ld] += c[k-1];
90: odd = PETSC_FALSE;
91: } else {
92: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,P,&ld,A2,&ld,&zero,W,&ld));
93: SWAP(P,W,aux);
94: for (j=0;j<n;j++) P[j+j*ld] += c[k-1];
95: odd = PETSC_TRUE;
96: }
97: PetscLogFlops(2.0*n*n*n);
98: }
99: /*if (odd) {
100: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,Q,&ld,As,&ld,&zero,W,&ld));
101: SWAP(Q,W,aux);
102: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&ld2,&mone,P,&inc,Q,&inc));
103: PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&n,&n,Q,&ld,ipiv,P,&ld,&info));
104: SlepcCheckLapackInfo("gesv",info);
105: PetscStackCallBLAS("BLASscal",BLASscal_(&ld2,&two,P,&inc));
106: for (j=0;j<n;j++) P[j+j*ld] += 1.0;
107: PetscStackCallBLAS("BLASscal",BLASscal_(&ld2,&mone,P,&inc));
108: } else {*/
109: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,P,&ld,As,&ld,&zero,W,&ld));
110: SWAP(P,W,aux);
111: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&ld2,&mone,P,&inc,Q,&inc));
112: PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&n,&n,Q,&ld,ipiv,P,&ld,&info));
113: SlepcCheckLapackInfo("gesv",info);
114: PetscStackCallBLAS("BLASscal",BLASscal_(&ld2,&two,P,&inc));
115: for (j=0;j<n;j++) P[j+j*ld] += 1.0;
116: /*}*/
117: PetscLogFlops(2.0*n*n*n+2.0*n*n*n/3.0+4.0*n*n);
119: for (k=1;k<=sexp;k++) {
120: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&one,P,&ld,P,&ld,&zero,W,&ld));
121: PetscArraycpy(P,W,ld2);
122: }
123: if (P!=Ba) { PetscArraycpy(Ba,P,ld2); }
124: PetscLogFlops(2.0*n*n*n*sexp);
126: PetscFree6(Q,W,As,A2,rwork,ipiv);
127: MatDenseRestoreArray(A,&Aa);
128: MatDenseRestoreArray(B,&Ba);
129: return(0);
130: #endif
131: }
133: /*
134: * Set scaling factor (s) and Pade degree (k,m)
135: */
136: static PetscErrorCode sexpm_params(PetscReal nrm,PetscInt *s,PetscInt *k,PetscInt *m)
137: {
139: if (nrm>1) {
140: if (nrm<200) {*s = 4; *k = 5; *m = *k-1;}
141: else if (nrm<1e4) {*s = 4; *k = 4; *m = *k+1;}
142: else if (nrm<1e6) {*s = 4; *k = 3; *m = *k+1;}
143: else if (nrm<1e9) {*s = 3; *k = 3; *m = *k+1;}
144: else if (nrm<1e11) {*s = 2; *k = 3; *m = *k+1;}
145: else if (nrm<1e12) {*s = 2; *k = 2; *m = *k+1;}
146: else if (nrm<1e14) {*s = 2; *k = 1; *m = *k+1;}
147: else {*s = 1; *k = 1; *m = *k+1;}
148: } else { /* nrm<1 */
149: if (nrm>0.5) {*s = 4; *k = 4; *m = *k-1;}
150: else if (nrm>0.3) {*s = 3; *k = 4; *m = *k-1;}
151: else if (nrm>0.15) {*s = 2; *k = 4; *m = *k-1;}
152: else if (nrm>0.07) {*s = 1; *k = 4; *m = *k-1;}
153: else if (nrm>0.01) {*s = 0; *k = 4; *m = *k-1;}
154: else if (nrm>3e-4) {*s = 0; *k = 3; *m = *k-1;}
155: else if (nrm>1e-5) {*s = 0; *k = 3; *m = 0;}
156: else if (nrm>1e-8) {*s = 0; *k = 2; *m = 0;}
157: else {*s = 0; *k = 1; *m = 0;}
158: }
159: return(0);
160: }
162: #if defined(PETSC_HAVE_COMPLEX)
163: /*
164: * Partial fraction form coefficients.
165: * If query, the function returns the size necessary to store the coefficients.
166: */
167: static PetscErrorCode getcoeffs(PetscInt k,PetscInt m,PetscComplex *r,PetscComplex *q,PetscComplex *remain,PetscBool query)
168: {
169: PetscInt i;
170: const PetscComplex /* m == k+1 */
171: p1r4[5] = {-1.582680186458572e+01 - 2.412564578224361e+01*PETSC_i,
172: -1.582680186458572e+01 + 2.412564578224361e+01*PETSC_i,
173: 1.499984465975511e+02 + 6.804227952202417e+01*PETSC_i,
174: 1.499984465975511e+02 - 6.804227952202417e+01*PETSC_i,
175: -2.733432894659307e+02 },
176: p1q4[5] = { 3.655694325463550e+00 + 6.543736899360086e+00*PETSC_i,
177: 3.655694325463550e+00 - 6.543736899360086e+00*PETSC_i,
178: 5.700953298671832e+00 + 3.210265600308496e+00*PETSC_i,
179: 5.700953298671832e+00 - 3.210265600308496e+00*PETSC_i,
180: 6.286704751729261e+00 },
181: p1r3[4] = {-1.130153999597152e+01 + 1.247167585025031e+01*PETSC_i,
182: -1.130153999597152e+01 - 1.247167585025031e+01*PETSC_i,
183: 1.330153999597152e+01 - 6.007173273704750e+01*PETSC_i,
184: 1.330153999597152e+01 + 6.007173273704750e+01*PETSC_i},
185: p1q3[4] = { 3.212806896871536e+00 + 4.773087433276636e+00*PETSC_i,
186: 3.212806896871536e+00 - 4.773087433276636e+00*PETSC_i,
187: 4.787193103128464e+00 + 1.567476416895212e+00*PETSC_i,
188: 4.787193103128464e+00 - 1.567476416895212e+00*PETSC_i},
189: p1r2[3] = { 7.648749087422928e+00 + 4.171640244747463e+00*PETSC_i,
190: 7.648749087422928e+00 - 4.171640244747463e+00*PETSC_i,
191: -1.829749817484586e+01 },
192: p1q2[3] = { 2.681082873627756e+00 + 3.050430199247411e+00*PETSC_i,
193: 2.681082873627756e+00 - 3.050430199247411e+00*PETSC_i,
194: 3.637834252744491e+00 },
195: p1r1[2] = { 1.000000000000000e+00 - 3.535533905932738e+00*PETSC_i,
196: 1.000000000000000e+00 + 3.535533905932738e+00*PETSC_i},
197: p1q1[2] = { 2.000000000000000e+00 + 1.414213562373095e+00*PETSC_i,
198: 2.000000000000000e+00 - 1.414213562373095e+00*PETSC_i};
199: const PetscComplex /* m == k-1 */
200: m1r5[4] = {-1.423367961376821e+02 - 1.385465094833037e+01*PETSC_i,
201: -1.423367961376821e+02 + 1.385465094833037e+01*PETSC_i,
202: 2.647367961376822e+02 - 4.814394493714596e+02*PETSC_i,
203: 2.647367961376822e+02 + 4.814394493714596e+02*PETSC_i},
204: m1q5[4] = { 5.203941240131764e+00 + 5.805856841805367e+00*PETSC_i,
205: 5.203941240131764e+00 - 5.805856841805367e+00*PETSC_i,
206: 6.796058759868242e+00 + 1.886649260140217e+00*PETSC_i,
207: 6.796058759868242e+00 - 1.886649260140217e+00*PETSC_i},
208: m1r4[3] = { 2.484269593165883e+01 + 7.460342395992306e+01*PETSC_i,
209: 2.484269593165883e+01 - 7.460342395992306e+01*PETSC_i,
210: -1.734353918633177e+02 },
211: m1q4[3] = { 4.675757014491557e+00 + 3.913489560603711e+00*PETSC_i,
212: 4.675757014491557e+00 - 3.913489560603711e+00*PETSC_i,
213: 5.648485971016893e+00 },
214: m1r3[2] = { 2.533333333333333e+01 - 2.733333333333333e+01*PETSC_i,
215: 2.533333333333333e+01 + 2.733333333333333e+01*PETSC_i},
216: m1q3[2] = { 4.000000000000000e+00 + 2.000000000000000e+00*PETSC_i,
217: 4.000000000000000e+00 - 2.000000000000000e+00*PETSC_i};
218: const PetscScalar /* m == k-1 */
219: m1remain5[2] = { 2.000000000000000e-01, 9.800000000000000e+00},
220: m1remain4[2] = {-2.500000000000000e-01, -7.750000000000000e+00},
221: m1remain3[2] = { 3.333333333333333e-01, 5.666666666666667e+00},
222: m1remain2[2] = {-0.5, -3.5},
223: remain3[4] = {1.0/6.0, 1.0/2.0, 1, 1},
224: remain2[3] = {1.0/2.0, 1, 1};
227: if (query) { /* query about buffer's size */
228: if (m==k+1) {
229: *remain = 0;
230: *r = *q = k+1;
231: return(0); /* quick return */
232: }
233: if (m==k-1) {
234: *remain = 2;
235: if (k==5) *r = *q = 4;
236: else if (k==4) *r = *q = 3;
237: else if (k==3) *r = *q = 2;
238: else if (k==2) *r = *q = 1;
239: }
240: if (m==0) {
241: *r = *q = 0;
242: *remain = k+1;
243: }
244: } else {
245: if (m==k+1) {
246: if (k==4) {
247: for (i=0;i<5;i++) { r[i] = p1r4[i]; q[i] = p1q4[i]; }
248: } else if (k==3) {
249: for (i=0;i<4;i++) { r[i] = p1r3[i]; q[i] = p1q3[i]; }
250: } else if (k==2) {
251: for (i=0;i<3;i++) { r[i] = p1r2[i]; q[i] = p1q2[i]; }
252: } else if (k==1) {
253: for (i=0;i<2;i++) { r[i] = p1r1[i]; q[i] = p1q1[i]; }
254: }
255: return(0); /* quick return */
256: }
257: if (m==k-1) {
258: if (k==5) {
259: for (i=0;i<4;i++) { r[i] = m1r5[i]; q[i] = m1q5[i]; }
260: for (i=0;i<2;i++) remain[i] = m1remain5[i];
261: } else if (k==4) {
262: for (i=0;i<3;i++) { r[i] = m1r4[i]; q[i] = m1q4[i]; }
263: for (i=0;i<2;i++) remain[i] = m1remain4[i];
264: } else if (k==3) {
265: for (i=0;i<2;i++) { r[i] = m1r3[i]; q[i] = m1q3[i]; remain[i] = m1remain3[i]; }
266: } else if (k==2) {
267: r[0] = -13.5; q[0] = 3;
268: for (i=0;i<2;i++) remain[i] = m1remain2[i];
269: }
270: }
271: if (m==0) {
272: r = q = 0;
273: if (k==3) {
274: for (i=0;i<4;i++) remain[i] = remain3[i];
275: } else if (k==2) {
276: for (i=0;i<3;i++) remain[i] = remain2[i];
277: }
278: }
279: }
280: return(0);
281: }
283: /*
284: * Product form coefficients.
285: * If query, the function returns the size necessary to store the coefficients.
286: */
287: static PetscErrorCode getcoeffsproduct(PetscInt k,PetscInt m,PetscComplex *p,PetscComplex *q,PetscComplex *mult,PetscBool query)
288: {
289: PetscInt i;
290: const PetscComplex /* m == k+1 */
291: p1p4[4] = {-5.203941240131764e+00 + 5.805856841805367e+00*PETSC_i,
292: -5.203941240131764e+00 - 5.805856841805367e+00*PETSC_i,
293: -6.796058759868242e+00 + 1.886649260140217e+00*PETSC_i,
294: -6.796058759868242e+00 - 1.886649260140217e+00*PETSC_i},
295: p1q4[5] = { 3.655694325463550e+00 + 6.543736899360086e+00*PETSC_i,
296: 3.655694325463550e+00 - 6.543736899360086e+00*PETSC_i,
297: 6.286704751729261e+00 ,
298: 5.700953298671832e+00 + 3.210265600308496e+00*PETSC_i,
299: 5.700953298671832e+00 - 3.210265600308496e+00*PETSC_i},
300: p1p3[3] = {-4.675757014491557e+00 + 3.913489560603711e+00*PETSC_i,
301: -4.675757014491557e+00 - 3.913489560603711e+00*PETSC_i,
302: -5.648485971016893e+00 },
303: p1q3[4] = { 3.212806896871536e+00 + 4.773087433276636e+00*PETSC_i,
304: 3.212806896871536e+00 - 4.773087433276636e+00*PETSC_i,
305: 4.787193103128464e+00 + 1.567476416895212e+00*PETSC_i,
306: 4.787193103128464e+00 - 1.567476416895212e+00*PETSC_i},
307: p1p2[2] = {-4.00000000000000e+00 + 2.000000000000000e+00*PETSC_i,
308: -4.00000000000000e+00 - 2.000000000000000e+00*PETSC_i},
309: p1q2[3] = { 2.681082873627756e+00 + 3.050430199247411e+00*PETSC_i,
310: 2.681082873627756e+00 - 3.050430199247411e+00*PETSC_i,
311: 3.637834252744491e+00 },
312: p1q1[2] = { 2.000000000000000e+00 + 1.414213562373095e+00*PETSC_i,
313: 2.000000000000000e+00 - 1.414213562373095e+00*PETSC_i};
314: const PetscComplex /* m == k-1 */
315: m1p5[5] = {-3.655694325463550e+00 + 6.543736899360086e+00*PETSC_i,
316: -3.655694325463550e+00 - 6.543736899360086e+00*PETSC_i,
317: -6.286704751729261e+00 ,
318: -5.700953298671832e+00 + 3.210265600308496e+00*PETSC_i,
319: -5.700953298671832e+00 - 3.210265600308496e+00*PETSC_i},
320: m1q5[4] = { 5.203941240131764e+00 + 5.805856841805367e+00*PETSC_i,
321: 5.203941240131764e+00 - 5.805856841805367e+00*PETSC_i,
322: 6.796058759868242e+00 + 1.886649260140217e+00*PETSC_i,
323: 6.796058759868242e+00 - 1.886649260140217e+00*PETSC_i},
324: m1p4[4] = {-3.212806896871536e+00 + 4.773087433276636e+00*PETSC_i,
325: -3.212806896871536e+00 - 4.773087433276636e+00*PETSC_i,
326: -4.787193103128464e+00 + 1.567476416895212e+00*PETSC_i,
327: -4.787193103128464e+00 - 1.567476416895212e+00*PETSC_i},
328: m1q4[3] = { 4.675757014491557e+00 + 3.913489560603711e+00*PETSC_i,
329: 4.675757014491557e+00 - 3.913489560603711e+00*PETSC_i,
330: 5.648485971016893e+00 },
331: m1p3[3] = {-2.681082873627756e+00 + 3.050430199247411e+00*PETSC_i,
332: -2.681082873627756e+00 - 3.050430199247411e+00*PETSC_i,
333: -3.637834252744491e+00 },
334: m1q3[2] = { 4.000000000000000e+00 + 2.000000000000000e+00*PETSC_i,
335: 4.000000000000000e+00 - 2.000000000000001e+00*PETSC_i},
336: m1p2[2] = {-2.000000000000000e+00 + 1.414213562373095e+00*PETSC_i,
337: -2.000000000000000e+00 - 1.414213562373095e+00*PETSC_i};
340: if (query) {
341: if (m == k+1) {
342: *mult = 1;
343: *p = k;
344: *q = k+1;
345: return(0);
346: }
347: if (m==k-1) {
348: *mult = 1;
349: *p = k;
350: *q = k-1;
351: }
352: } else {
353: if (m == k+1) {
354: *mult = PetscPowInt(-1,m);
355: *mult *= m;
356: if (k==4) {
357: for (i=0;i<4;i++) { p[i] = p1p4[i]; q[i] = p1q4[i]; }
358: q[4] = p1q4[4];
359: } else if (k==3) {
360: for (i=0;i<3;i++) { p[i] = p1p3[i]; q[i] = p1q3[i]; }
361: q[3] = p1q3[3];
362: } else if (k==2) {
363: for (i=0;i<2;i++) { p[i] = p1p2[i]; q[i] = p1q2[i]; }
364: q[2] = p1q2[2];
365: } else if (k==1) {
366: p[0] = -3;
367: for (i=0;i<2;i++) q[i] = p1q1[i];
368: }
369: return(0);
370: }
371: if (m==k-1) {
372: *mult = PetscPowInt(-1,m);
373: *mult /= k;
374: if (k==5) {
375: for (i=0;i<4;i++) { p[i] = m1p5[i]; q[i] = m1q5[i]; }
376: p[4] = m1p5[4];
377: } else if (k==4) {
378: for (i=0;i<3;i++) { p[i] = m1p4[i]; q[i] = m1q4[i]; }
379: p[3] = m1p4[3];
380: } else if (k==3) {
381: for (i=0;i<2;i++) { p[i] = m1p3[i]; q[i] = m1q3[i]; }
382: p[2] = m1p3[2];
383: } else if (k==2) {
384: for (i=0;i<2;i++) p[i] = m1p2[i];
385: q[0] = 3;
386: }
387: }
388: }
389: return(0);
390: }
391: #endif /* PETSC_HAVE_COMPLEX */
393: #if defined(PETSC_USE_COMPLEX)
394: static PetscErrorCode getisreal(PetscInt n,PetscComplex *a,PetscBool *result)
395: {
396: PetscInt i;
399: *result=PETSC_TRUE;
400: for (i=0;i<n&&*result;i++) {
401: if (PetscImaginaryPartComplex(a[i])) *result=PETSC_FALSE;
402: }
403: return(0);
404: }
405: #endif
407: /*
408: * Matrix exponential implementation based on algorithm and matlab code by Stefan Guettel
409: * and Yuji Nakatsukasa
410: *
411: * Stefan Guettel and Yuji Nakatsukasa, "Scaled and Squared Subdiagonal Pade
412: * Approximation for the Matrix Exponential",
413: * SIAM J. Matrix Anal. Appl. 37(1):145-170, 2016.
414: * https://doi.org/10.1137/15M1027553
415: */
416: PetscErrorCode FNEvaluateFunctionMat_Exp_GuettelNakatsukasa(FN fn,Mat A,Mat B)
417: {
418: #if !defined(PETSC_HAVE_COMPLEX)
420: SETERRQ(PETSC_COMM_SELF,1,"This function requires C99 or C++ complex support");
421: #elif defined(PETSC_MISSING_LAPACK_GEEV) || defined(PETSC_MISSING_LAPACK_GESV) || defined(SLEPC_MISSING_LAPACK_LANGE)
423: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEEV/GESV/LANGE - Lapack routines are unavailable");
424: #else
425: PetscInt i,j,n_,s,k,m,mod;
426: PetscBLASInt n,n2,irsize,rsizediv2,ipsize,iremainsize,query=-1,info,*piv,minlen,lwork,one=1;
427: PetscReal nrm,shift;
428: #if defined(PETSC_USE_COMPLEX)
429: PetscReal *rwork=NULL;
430: #endif
431: PetscComplex *As,*RR,*RR2,*expmA,*expmA2,*Maux,*Maux2,rsize,*r,psize,*p,remainsize,*remainterm,*rootp,*rootq,mult=0.0,scale,cone=1.0,czero=0.0,*aux;
432: PetscScalar *Aa,*Ba,*Ba2,*sMaux,*wr,*wi,expshift,sone=1.0,szero=0.0,*work,work1,*saux;
434: PetscBool isreal;
435: #if defined(PETSC_HAVE_ESSL)
436: PetscScalar sdummy;
437: PetscBLASInt idummy,io=0;
438: PetscScalar *wri;
439: #endif
442: MatGetSize(A,&n_,NULL);
443: PetscBLASIntCast(n_,&n);
444: MatDenseGetArray(A,&Aa);
445: MatDenseGetArray(B,&Ba);
446: Ba2 = Ba;
447: PetscBLASIntCast(n*n,&n2);
449: PetscMalloc2(n2,&sMaux,n2,&Maux);
450: Maux2 = Maux;
451: PetscMalloc2(n,&wr,n,&wi);
452: PetscArraycpy(sMaux,Aa,n2);
453: /* estimate rightmost eigenvalue and shift A with it */
454: #if !defined(PETSC_HAVE_ESSL)
455: #if !defined(PETSC_USE_COMPLEX)
456: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&n,sMaux,&n,wr,wi,NULL,&n,NULL,&n,&work1,&query,&info));
457: SlepcCheckLapackInfo("geev",info);
458: PetscBLASIntCast((PetscInt)PetscRealPart(work1),&lwork);
459: PetscMalloc1(lwork,&work);
460: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&n,sMaux,&n,wr,wi,NULL,&n,NULL,&n,work,&lwork,&info));
461: PetscFree(work);
462: #else
463: PetscArraycpy(Maux,Aa,n2);
464: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&n,Maux,&n,wr,NULL,&n,NULL,&n,&work1,&query,rwork,&info));
465: SlepcCheckLapackInfo("geev",info);
466: PetscBLASIntCast((PetscInt)PetscRealPart(work1),&lwork);
467: PetscMalloc2(2*n,&rwork,lwork,&work);
468: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&n,Maux,&n,wr,NULL,&n,NULL,&n,work,&lwork,rwork,&info));
469: PetscFree2(rwork,work);
470: #endif
471: SlepcCheckLapackInfo("geev",info);
472: #else /* defined(PETSC_HAVE_ESSL) */
473: PetscBLASIntCast(3*n,&lwork);
474: PetscMalloc2(lwork,&work,2*n,&wri);
475: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_(&io,Maux,&n,wri,&sdummy,&idummy,&idummy,&n,work,&lwork));
476: #if !defined(PETSC_USE_COMPLEX)
477: for (i=0;i<n;i++) {
478: wr[i] = wri[2*i];
479: wi[i] = wri[2*i+1];
480: }
481: #else
482: for (i=0;i<n;i++) wr[i] = wri[i];
483: #endif
484: PetscFree2(work,wri);
485: #endif
486: PetscLogFlops(25.0*n*n*n+(n*n*n)/3.0+1.0*n*n*n);
488: shift = PetscRealPart(wr[0]);
489: for (i=1;i<n;i++) {
490: if (PetscRealPart(wr[i]) > shift) shift = PetscRealPart(wr[i]);
491: }
492: PetscFree2(wr,wi);
493: /* shift so that largest real part is (about) 0 */
494: PetscArraycpy(sMaux,Aa,n2);
495: for (i=0;i<n;i++) {
496: sMaux[i+i*n] -= shift;
497: }
498: PetscLogFlops(1.0*n);
499: #if defined(PETSC_USE_COMPLEX)
500: PetscArraycpy(Maux,Aa,n2);
501: for (i=0;i<n;i++) {
502: Maux[i+i*n] -= shift;
503: }
504: PetscLogFlops(1.0*n);
505: #endif
507: /* estimate norm(A) and select the scaling factor */
508: nrm = LAPACKlange_("O",&n,&n,sMaux,&n,NULL);
509: PetscLogFlops(1.0*n*n);
510: sexpm_params(nrm,&s,&k,&m);
511: if (s==0 && k==1 && m==0) { /* exp(A) = I+A to eps! */
512: expshift = PetscExpReal(shift);
513: for (i=0;i<n;i++) sMaux[i+i*n] += 1.0;
514: PetscStackCallBLAS("BLASscal",BLASscal_(&n2,&expshift,sMaux,&one));
515: PetscLogFlops(1.0*(n+n2));
516: PetscArraycpy(Ba,sMaux,n2);
517: PetscFree2(sMaux,Maux);
518: MatDenseRestoreArray(A,&Aa);
519: MatDenseRestoreArray(B,&Ba);
520: return(0); /* quick return */
521: }
523: PetscMalloc4(n2,&expmA,n2,&As,n2,&RR,n,&piv);
524: expmA2 = expmA; RR2 = RR;
525: /* scale matrix */
526: #if !defined(PETSC_USE_COMPLEX)
527: for (i=0;i<n2;i++) {
528: As[i] = sMaux[i];
529: }
530: #else
531: PetscArraycpy(As,sMaux,n2);
532: #endif
533: scale = 1.0/PetscPowRealInt(2.0,s);
534: PetscStackCallBLAS("BLASCOMPLEXscal",BLASCOMPLEXscal_(&n2,&scale,As,&one));
535: SlepcLogFlopsComplex(1.0*n2);
537: /* evaluate Pade approximant (partial fraction or product form) */
538: if (fn->method==3 || !m) { /* partial fraction */
539: getcoeffs(k,m,&rsize,&psize,&remainsize,PETSC_TRUE);
540: PetscBLASIntCast((PetscInt)PetscRealPartComplex(rsize),&irsize);
541: PetscBLASIntCast((PetscInt)PetscRealPartComplex(psize),&ipsize);
542: PetscBLASIntCast((PetscInt)PetscRealPartComplex(remainsize),&iremainsize);
543: PetscMalloc3(irsize,&r,ipsize,&p,iremainsize,&remainterm);
544: getcoeffs(k,m,r,p,remainterm,PETSC_FALSE);
546: PetscArrayzero(expmA,n2);
547: #if !defined(PETSC_USE_COMPLEX)
548: isreal = PETSC_TRUE;
549: #else
550: getisreal(n2,Maux,&isreal);
551: #endif
552: if (isreal) {
553: rsizediv2 = irsize/2;
554: for (i=0;i<rsizediv2;i++) { /* use partial fraction to get R(As) */
555: PetscArraycpy(Maux,As,n2);
556: PetscArrayzero(RR,n2);
557: for (j=0;j<n;j++) {
558: Maux[j+j*n] -= p[2*i];
559: RR[j+j*n] = r[2*i];
560: }
561: PetscStackCallBLAS("LAPACKCOMPLEXgesv",LAPACKCOMPLEXgesv_(&n,&n,Maux,&n,piv,RR,&n,&info));
562: SlepcCheckLapackInfo("gesv",info);
563: for (j=0;j<n2;j++) {
564: expmA[j] += RR[j] + PetscConj(RR[j]);
565: }
566: /* loop(n) + gesv + loop(n2) */
567: SlepcLogFlopsComplex(1.0*n+(2.0*n*n*n/3.0+2.0*n*n*n)+2.0*n2);
568: }
570: mod = ipsize % 2;
571: if (mod) {
572: PetscArraycpy(Maux,As,n2);
573: PetscArrayzero(RR,n2);
574: for (j=0;j<n;j++) {
575: Maux[j+j*n] -= p[ipsize-1];
576: RR[j+j*n] = r[irsize-1];
577: }
578: PetscStackCallBLAS("LAPACKCOMPLEXgesv",LAPACKCOMPLEXgesv_(&n,&n,Maux,&n,piv,RR,&n,&info));
579: SlepcCheckLapackInfo("gesv",info);
580: for (j=0;j<n2;j++) {
581: expmA[j] += RR[j];
582: }
583: SlepcLogFlopsComplex(1.0*n+(2.0*n*n*n/3.0+2.0*n*n*n)+1.0*n2);
584: }
585: } else { /* complex */
586: for (i=0;i<irsize;i++) { /* use partial fraction to get R(As) */
587: PetscArraycpy(Maux,As,n2);
588: PetscArrayzero(RR,n2);
589: for (j=0;j<n;j++) {
590: Maux[j+j*n] -= p[i];
591: RR[j+j*n] = r[i];
592: }
593: PetscStackCallBLAS("LAPACKCOMPLEXgesv",LAPACKCOMPLEXgesv_(&n,&n,Maux,&n,piv,RR,&n,&info));
594: SlepcCheckLapackInfo("gesv",info);
595: for (j=0;j<n2;j++) {
596: expmA[j] += RR[j];
597: }
598: SlepcLogFlopsComplex(1.0*n+(2.0*n*n*n/3.0+2.0*n*n*n)+1.0*n2);
599: }
600: }
601: for (i=0;i<iremainsize;i++) {
602: if (!i) {
603: PetscArrayzero(RR,n2);
604: for (j=0;j<n;j++) {
605: RR[j+j*n] = remainterm[iremainsize-1];
606: }
607: } else {
608: PetscArraycpy(RR,As,n2);
609: for (j=1;j<i;j++) {
610: PetscStackCallBLAS("BLASCOMPLEXgemm",BLASCOMPLEXgemm_("N","N",&n,&n,&n,&cone,RR,&n,RR,&n,&czero,Maux,&n));
611: SWAP(RR,Maux,aux);
612: SlepcLogFlopsComplex(2.0*n*n*n);
613: }
614: PetscStackCallBLAS("BLASCOMPLEXscal",BLASCOMPLEXscal_(&n2,&remainterm[iremainsize-1-i],RR,&one));
615: SlepcLogFlopsComplex(1.0*n2);
616: }
617: for (j=0;j<n2;j++) {
618: expmA[j] += RR[j];
619: }
620: SlepcLogFlopsComplex(1.0*n2);
621: }
622: PetscFree3(r,p,remainterm);
623: } else { /* product form, default */
624: getcoeffsproduct(k,m,&rsize,&psize,&mult,PETSC_TRUE);
625: PetscBLASIntCast((PetscInt)PetscRealPartComplex(rsize),&irsize);
626: PetscBLASIntCast((PetscInt)PetscRealPartComplex(psize),&ipsize);
627: PetscMalloc2(irsize,&rootp,ipsize,&rootq);
628: getcoeffsproduct(k,m,rootp,rootq,&mult,PETSC_FALSE);
630: PetscArrayzero(expmA,n2);
631: for (i=0;i<n;i++) { /* initialize */
632: expmA[i+i*n] = 1.0;
633: }
634: minlen = PetscMin(irsize,ipsize);
635: for (i=0;i<minlen;i++) {
636: PetscArraycpy(RR,As,n2);
637: for (j=0;j<n;j++) {
638: RR[j+j*n] -= rootp[i];
639: }
640: PetscStackCallBLAS("BLASCOMPLEXgemm",BLASCOMPLEXgemm_("N","N",&n,&n,&n,&cone,RR,&n,expmA,&n,&czero,Maux,&n));
641: SWAP(expmA,Maux,aux);
642: PetscArraycpy(RR,As,n2);
643: for (j=0;j<n;j++) {
644: RR[j+j*n] -= rootq[i];
645: }
646: PetscStackCallBLAS("LAPACKCOMPLEXgesv",LAPACKCOMPLEXgesv_(&n,&n,RR,&n,piv,expmA,&n,&info));
647: SlepcCheckLapackInfo("gesv",info);
648: /* loop(n) + gemm + loop(n) + gesv */
649: SlepcLogFlopsComplex(1.0*n+(2.0*n*n*n)+1.0*n+(2.0*n*n*n/3.0+2.0*n*n*n));
650: }
651: /* extra numerator */
652: for (i=minlen;i<irsize;i++) {
653: PetscArraycpy(RR,As,n2);
654: for (j=0;j<n;j++) {
655: RR[j+j*n] -= rootp[i];
656: }
657: PetscStackCallBLAS("BLASCOMPLEXgemm",BLASCOMPLEXgemm_("N","N",&n,&n,&n,&cone,RR,&n,expmA,&n,&czero,Maux,&n));
658: SWAP(expmA,Maux,aux);
659: SlepcLogFlopsComplex(1.0*n+2.0*n*n*n);
660: }
661: /* extra denominator */
662: for (i=minlen;i<ipsize;i++) {
663: PetscArraycpy(RR,As,n2);
664: for (j=0;j<n;j++) RR[j+j*n] -= rootq[i];
665: PetscStackCallBLAS("LAPACKCOMPLEXgesv",LAPACKCOMPLEXgesv_(&n,&n,RR,&n,piv,expmA,&n,&info));
666: SlepcCheckLapackInfo("gesv",info);
667: SlepcLogFlopsComplex(1.0*n+(2.0*n*n*n/3.0+2.0*n*n*n));
668: }
669: PetscStackCallBLAS("BLASCOMPLEXscal",BLASCOMPLEXscal_(&n2,&mult,expmA,&one));
670: SlepcLogFlopsComplex(1.0*n2);
671: PetscFree2(rootp,rootq);
672: }
674: #if !defined(PETSC_USE_COMPLEX)
675: for (i=0;i<n2;i++) {
676: Ba2[i] = PetscRealPartComplex(expmA[i]);
677: }
678: #else
679: PetscArraycpy(Ba2,expmA,n2);
680: #endif
682: /* perform repeated squaring */
683: for (i=0;i<s;i++) { /* final squaring */
684: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&sone,Ba2,&n,Ba2,&n,&szero,sMaux,&n));
685: SWAP(Ba2,sMaux,saux);
686: PetscLogFlops(2.0*n*n*n);
687: }
688: if (Ba2!=Ba) {
689: PetscArraycpy(Ba,Ba2,n2);
690: sMaux = Ba2;
691: }
692: expshift = PetscExpReal(shift);
693: PetscStackCallBLAS("BLASscal",BLASscal_(&n2,&expshift,Ba,&one));
694: PetscLogFlops(1.0*n2);
696: /* restore pointers */
697: Maux = Maux2; expmA = expmA2; RR = RR2;
698: PetscFree2(sMaux,Maux);
699: PetscFree4(expmA,As,RR,piv);
700: MatDenseRestoreArray(A,&Aa);
701: MatDenseRestoreArray(B,&Ba);
702: return(0);
703: #endif
704: }
706: #define SMALLN 100
708: /*
709: * Function needed to compute optimal parameters (required workspace is 3*n*n)
710: */
711: static PetscInt ell(PetscBLASInt n,PetscScalar *A,PetscReal coeff,PetscInt m,PetscScalar *work,PetscRandom rand)
712: {
713: PetscScalar *Ascaled=work;
714: PetscReal nrm,alpha,beta,rwork[1];
715: PetscInt t;
716: PetscBLASInt i,j;
720: beta = PetscPowReal(coeff,1.0/(2*m+1));
721: for (i=0;i<n;i++)
722: for (j=0;j<n;j++)
723: Ascaled[i+j*n] = beta*PetscAbsScalar(A[i+j*n]);
724: nrm = LAPACKlange_("O",&n,&n,A,&n,rwork);
725: PetscLogFlops(2.0*n*n);
726: SlepcNormAm(n,Ascaled,2*m+1,work+n*n,rand,&alpha);
727: alpha /= nrm;
728: t = PetscMax((PetscInt)PetscCeilReal(PetscLogReal(2.0*alpha/PETSC_MACHINE_EPSILON)/PetscLogReal(2.0)/(2*m)),0);
729: PetscFunctionReturn(t);
730: }
732: /*
733: * Compute scaling parameter (s) and order of Pade approximant (m) (required workspace is 4*n*n)
734: */
735: static PetscErrorCode expm_params(PetscInt n,PetscScalar **Apowers,PetscInt *s,PetscInt *m,PetscScalar *work)
736: {
737: PetscErrorCode ierr;
738: PetscScalar sfactor,sone=1.0,szero=0.0,*A=Apowers[0],*Ascaled;
739: PetscReal d4,d6,d8,d10,eta1,eta3,eta4,eta5,rwork[1];
740: PetscBLASInt n_,n2,one=1;
741: PetscRandom rand;
742: const PetscReal coeff[5] = { 9.92063492063492e-06, 9.94131285136576e-11, /* backward error function */
743: 2.22819456055356e-16, 1.69079293431187e-22, 8.82996160201868e-36 };
744: const PetscReal theta[5] = { 1.495585217958292e-002, /* m = 3 */
745: 2.539398330063230e-001, /* m = 5 */
746: 9.504178996162932e-001, /* m = 7 */
747: 2.097847961257068e+000, /* m = 9 */
748: 5.371920351148152e+000 }; /* m = 13 */
751: *s = 0;
752: *m = 13;
753: PetscBLASIntCast(n,&n_);
754: PetscRandomCreate(PETSC_COMM_SELF,&rand);
755: d4 = PetscPowReal(LAPACKlange_("O",&n_,&n_,Apowers[2],&n_,rwork),1.0/4.0);
756: if (d4==0.0) { /* safeguard for the case A = 0 */
757: *m = 3;
758: goto done;
759: }
760: d6 = PetscPowReal(LAPACKlange_("O",&n_,&n_,Apowers[3],&n_,rwork),1.0/6.0);
761: PetscLogFlops(2.0*n*n);
762: eta1 = PetscMax(d4,d6);
763: if (eta1<=theta[0] && !ell(n_,A,coeff[0],3,work,rand)) {
764: *m = 3;
765: goto done;
766: }
767: if (eta1<=theta[1] && !ell(n_,A,coeff[1],5,work,rand)) {
768: *m = 5;
769: goto done;
770: }
771: if (n<SMALLN) {
772: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,Apowers[2],&n_,Apowers[2],&n_,&szero,work,&n_));
773: d8 = PetscPowReal(LAPACKlange_("O",&n_,&n_,work,&n_,rwork),1.0/8.0);
774: PetscLogFlops(2.0*n*n*n+1.0*n*n);
775: } else {
776: SlepcNormAm(n_,Apowers[2],2,work,rand,&d8);
777: d8 = PetscPowReal(d8,1.0/8.0);
778: }
779: eta3 = PetscMax(d6,d8);
780: if (eta3<=theta[2] && !ell(n_,A,coeff[2],7,work,rand)) {
781: *m = 7;
782: goto done;
783: }
784: if (eta3<=theta[3] && !ell(n_,A,coeff[3],9,work,rand)) {
785: *m = 9;
786: goto done;
787: }
788: if (n<SMALLN) {
789: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,Apowers[2],&n_,Apowers[3],&n_,&szero,work,&n_));
790: d10 = PetscPowReal(LAPACKlange_("O",&n_,&n_,work,&n_,rwork),1.0/10.0);
791: PetscLogFlops(2.0*n*n*n+1.0*n*n);
792: } else {
793: SlepcNormAm(n_,Apowers[1],5,work,rand,&d10);
794: d10 = PetscPowReal(d10,1.0/10.0);
795: }
796: eta4 = PetscMax(d8,d10);
797: eta5 = PetscMin(eta3,eta4);
798: *s = PetscMax((PetscInt)PetscCeilReal(PetscLogReal(eta5/theta[4])/PetscLogReal(2.0)),0);
799: if (*s) {
800: Ascaled = work+3*n*n;
801: n2 = n_*n_;
802: PetscStackCallBLAS("BLAScopy",BLAScopy_(&n2,A,&one,Ascaled,&one));
803: sfactor = PetscPowRealInt(2.0,-(*s));
804: PetscStackCallBLAS("BLASscal",BLASscal_(&n2,&sfactor,Ascaled,&one));
805: PetscLogFlops(1.0*n*n);
806: } else Ascaled = A;
807: *s += ell(n_,Ascaled,coeff[4],13,work,rand);
808: done:
809: PetscRandomDestroy(&rand);
810: return(0);
811: }
813: /*
814: * Matrix exponential implementation based on algorithm and matlab code by N. Higham and co-authors
815: *
816: * N. J. Higham, "The scaling and squaring method for the matrix exponential
817: * revisited", SIAM J. Matrix Anal. Appl. 26(4):1179-1193, 2005.
818: */
819: PetscErrorCode FNEvaluateFunctionMat_Exp_Higham(FN fn,Mat A,Mat B)
820: {
821: #if defined(PETSC_MISSING_LAPACK_GESV) || defined(SLEPC_MISSING_LAPACK_LANGE)
823: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESV/LANGE - Lapack routines are unavailable");
824: #else
825: PetscErrorCode ierr;
826: PetscBLASInt n_,n2,*ipiv,info,one=1;
827: PetscInt n,m,j,s;
828: PetscScalar scale,smone=-1.0,sone=1.0,stwo=2.0,szero=0.0;
829: PetscScalar *Aa,*Ba,*Apowers[5],*Q,*P,*W,*work,*aux;
830: const PetscScalar *c;
831: const PetscScalar c3[4] = { 120, 60, 12, 1 };
832: const PetscScalar c5[6] = { 30240, 15120, 3360, 420, 30, 1 };
833: const PetscScalar c7[8] = { 17297280, 8648640, 1995840, 277200, 25200, 1512, 56, 1 };
834: const PetscScalar c9[10] = { 17643225600, 8821612800, 2075673600, 302702400, 30270240,
835: 2162160, 110880, 3960, 90, 1 };
836: const PetscScalar c13[14] = { 64764752532480000, 32382376266240000, 7771770303897600,
837: 1187353796428800, 129060195264000, 10559470521600,
838: 670442572800, 33522128640, 1323241920,
839: 40840800, 960960, 16380, 182, 1 };
842: MatDenseGetArray(A,&Aa);
843: MatDenseGetArray(B,&Ba);
844: MatGetSize(A,&n,NULL);
845: PetscBLASIntCast(n,&n_);
846: n2 = n_*n_;
847: PetscMalloc2(8*n*n,&work,n,&ipiv);
849: /* Matrix powers */
850: Apowers[0] = work; /* Apowers[0] = A */
851: Apowers[1] = Apowers[0] + n*n; /* Apowers[1] = A^2 */
852: Apowers[2] = Apowers[1] + n*n; /* Apowers[2] = A^4 */
853: Apowers[3] = Apowers[2] + n*n; /* Apowers[3] = A^6 */
854: Apowers[4] = Apowers[3] + n*n; /* Apowers[4] = A^8 */
856: PetscArraycpy(Apowers[0],Aa,n2);
857: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,Apowers[0],&n_,Apowers[0],&n_,&szero,Apowers[1],&n_));
858: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,Apowers[1],&n_,Apowers[1],&n_,&szero,Apowers[2],&n_));
859: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,Apowers[1],&n_,Apowers[2],&n_,&szero,Apowers[3],&n_));
860: PetscLogFlops(6.0*n*n*n);
862: /* Compute scaling parameter and order of Pade approximant */
863: expm_params(n,Apowers,&s,&m,Apowers[4]);
865: if (s) { /* rescale */
866: for (j=0;j<4;j++) {
867: scale = PetscPowRealInt(2.0,-PetscMax(2*j,1)*s);
868: PetscStackCallBLAS("BLASscal",BLASscal_(&n2,&scale,Apowers[j],&one));
869: }
870: PetscLogFlops(4.0*n*n);
871: }
873: /* Evaluate the Pade approximant */
874: switch (m) {
875: case 3: c = c3; break;
876: case 5: c = c5; break;
877: case 7: c = c7; break;
878: case 9: c = c9; break;
879: case 13: c = c13; break;
880: default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong value of m %d",m);
881: }
882: P = Ba;
883: Q = Apowers[4] + n*n;
884: W = Q + n*n;
885: switch (m) {
886: case 3:
887: case 5:
888: case 7:
889: case 9:
890: if (m==9) PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,Apowers[1],&n_,Apowers[3],&n_,&szero,Apowers[4],&n_));
891: PetscArrayzero(P,n2);
892: PetscArrayzero(Q,n2);
893: for (j=0;j<n;j++) {
894: P[j+j*n] = c[1];
895: Q[j+j*n] = c[0];
896: }
897: for (j=m;j>=3;j-=2) {
898: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[j],Apowers[(j+1)/2-1],&one,P,&one));
899: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[j-1],Apowers[(j+1)/2-1],&one,Q,&one));
900: PetscLogFlops(4.0*n*n);
901: }
902: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,Apowers[0],&n_,P,&n_,&szero,W,&n_));
903: PetscLogFlops(2.0*n*n*n);
904: SWAP(P,W,aux);
905: break;
906: case 13:
907: /* P = A*(Apowers[3]*(c[13]*Apowers[3] + c[11]*Apowers[2] + c[9]*Apowers[1])
908: + c[7]*Apowers[3] + c[5]*Apowers[2] + c[3]*Apowers[1] + c[1]*I) */
909: PetscStackCallBLAS("BLAScopy",BLAScopy_(&n2,Apowers[3],&one,P,&one));
910: PetscStackCallBLAS("BLASscal",BLASscal_(&n2,&c[13],P,&one));
911: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[11],Apowers[2],&one,P,&one));
912: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[9],Apowers[1],&one,P,&one));
913: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,Apowers[3],&n_,P,&n_,&szero,W,&n_));
914: PetscLogFlops(5.0*n*n+2.0*n*n*n);
915: PetscArrayzero(P,n2);
916: for (j=0;j<n;j++) P[j+j*n] = c[1];
917: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[7],Apowers[3],&one,P,&one));
918: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[5],Apowers[2],&one,P,&one));
919: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[3],Apowers[1],&one,P,&one));
920: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&sone,P,&one,W,&one));
921: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,Apowers[0],&n_,W,&n_,&szero,P,&n_));
922: PetscLogFlops(7.0*n*n+2.0*n*n*n);
923: /* Q = Apowers[3]*(c[12]*Apowers[3] + c[10]*Apowers[2] + c[8]*Apowers[1])
924: + c[6]*Apowers[3] + c[4]*Apowers[2] + c[2]*Apowers[1] + c[0]*I */
925: PetscStackCallBLAS("BLAScopy",BLAScopy_(&n2,Apowers[3],&one,Q,&one));
926: PetscStackCallBLAS("BLASscal",BLASscal_(&n2,&c[12],Q,&one));
927: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[10],Apowers[2],&one,Q,&one));
928: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[8],Apowers[1],&one,Q,&one));
929: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,Apowers[3],&n_,Q,&n_,&szero,W,&n_));
930: PetscLogFlops(5.0*n*n+2.0*n*n*n);
931: PetscArrayzero(Q,n2);
932: for (j=0;j<n;j++) Q[j+j*n] = c[0];
933: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[6],Apowers[3],&one,Q,&one));
934: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[4],Apowers[2],&one,Q,&one));
935: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&c[2],Apowers[1],&one,Q,&one));
936: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&sone,W,&one,Q,&one));
937: PetscLogFlops(7.0*n*n);
938: break;
939: default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong value of m %d",m);
940: }
941: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&n2,&smone,P,&one,Q,&one));
942: PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&n_,&n_,Q,&n_,ipiv,P,&n_,&info));
943: SlepcCheckLapackInfo("gesv",info);
944: PetscStackCallBLAS("BLASscal",BLASscal_(&n2,&stwo,P,&one));
945: for (j=0;j<n;j++) P[j+j*n] += 1.0;
946: PetscLogFlops(2.0*n*n*n/3.0+4.0*n*n);
948: /* Squaring */
949: for (j=1;j<=s;j++) {
950: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,P,&n_,P,&n_,&szero,W,&n_));
951: SWAP(P,W,aux);
952: }
953: if (P!=Ba) { PetscArraycpy(Ba,P,n2); }
954: PetscLogFlops(2.0*n*n*n*s);
956: PetscFree2(work,ipiv);
957: MatDenseRestoreArray(A,&Aa);
958: MatDenseRestoreArray(B,&Ba);
959: return(0);
960: #endif
961: }
963: PetscErrorCode FNView_Exp(FN fn,PetscViewer viewer)
964: {
966: PetscBool isascii;
967: char str[50];
968: const char *methodname[] = {
969: "scaling & squaring, [m/m] Pade approximant (Higham)",
970: "scaling & squaring, [6/6] Pade approximant",
971: "scaling & squaring, subdiagonal Pade approximant (product form)",
972: "scaling & squaring, subdiagonal Pade approximant (partial fraction)"
973: };
974: const int nmeth=sizeof(methodname)/sizeof(methodname[0]);
977: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
978: if (isascii) {
979: if (fn->beta==(PetscScalar)1.0) {
980: if (fn->alpha==(PetscScalar)1.0) {
981: PetscViewerASCIIPrintf(viewer," Exponential: exp(x)\n");
982: } else {
983: SlepcSNPrintfScalar(str,50,fn->alpha,PETSC_TRUE);
984: PetscViewerASCIIPrintf(viewer," Exponential: exp(%s*x)\n",str);
985: }
986: } else {
987: SlepcSNPrintfScalar(str,50,fn->beta,PETSC_TRUE);
988: if (fn->alpha==(PetscScalar)1.0) {
989: PetscViewerASCIIPrintf(viewer," Exponential: %s*exp(x)\n",str);
990: } else {
991: PetscViewerASCIIPrintf(viewer," Exponential: %s",str);
992: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
993: SlepcSNPrintfScalar(str,50,fn->alpha,PETSC_TRUE);
994: PetscViewerASCIIPrintf(viewer,"*exp(%s*x)\n",str);
995: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
996: }
997: }
998: if (fn->method<nmeth) {
999: PetscViewerASCIIPrintf(viewer," computing matrix functions with: %s\n",methodname[fn->method]);
1000: }
1001: }
1002: return(0);
1003: }
1005: SLEPC_EXTERN PetscErrorCode FNCreate_Exp(FN fn)
1006: {
1008: fn->ops->evaluatefunction = FNEvaluateFunction_Exp;
1009: fn->ops->evaluatederivative = FNEvaluateDerivative_Exp;
1010: fn->ops->evaluatefunctionmat[0] = FNEvaluateFunctionMat_Exp_Higham;
1011: fn->ops->evaluatefunctionmat[1] = FNEvaluateFunctionMat_Exp_Pade;
1012: fn->ops->evaluatefunctionmat[2] = FNEvaluateFunctionMat_Exp_GuettelNakatsukasa; /* product form */
1013: fn->ops->evaluatefunctionmat[3] = FNEvaluateFunctionMat_Exp_GuettelNakatsukasa; /* partial fraction */
1014: fn->ops->view = FNView_Exp;
1015: return(0);
1016: }