Actual source code: fnphi.c
slepc-3.11.2 2019-07-30
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: Phi functions
12: phi_0(x) = exp(x)
13: phi_1(x) = (exp(x)-1)/x
14: phi_k(x) = (phi_{k-1}(x)-1/(k-1)!)/x
15: */
17: #include <slepc/private/fnimpl.h> /*I "slepcfn.h" I*/
19: typedef struct {
20: PetscInt k; /* index of the phi-function, defaults to k=1 */
21: Mat H; /* auxiliary matrix of order m+k */
22: Mat F; /* auxiliary matrix to store exp(H) */
23: } FN_PHI;
25: #define MAX_INDEX 10
27: const static PetscReal rfactorial[MAX_INDEX+2] = { 1, 1, 0.5, 1.0/6, 1.0/24, 1.0/120, 1.0/720, 1.0/5040, 1.0/40320, 1.0/362880, 1.0/3628800, 1.0/39916800 };
29: PetscErrorCode FNEvaluateFunction_Phi(FN fn,PetscScalar x,PetscScalar *y)
30: {
31: FN_PHI *ctx = (FN_PHI*)fn->data;
32: PetscInt i;
33: PetscScalar phi[MAX_INDEX+1];
36: if (x==0.0) *y = rfactorial[ctx->k];
37: else {
38: phi[0] = PetscExpScalar(x);
39: for (i=1;i<=ctx->k;i++) phi[i] = (phi[i-1]-rfactorial[i-1])/x;
40: *y = phi[ctx->k];
41: }
42: return(0);
43: }
45: PetscErrorCode FNEvaluateDerivative_Phi(FN fn,PetscScalar x,PetscScalar *y)
46: {
47: FN_PHI *ctx = (FN_PHI*)fn->data;
48: PetscInt i;
49: PetscScalar phi[MAX_INDEX+2];
52: if (x==0.0) *y = rfactorial[ctx->k+1];
53: else {
54: phi[0] = PetscExpScalar(x);
55: for (i=1;i<=ctx->k+1;i++) phi[i] = (phi[i-1]-rfactorial[i-1])/x;
56: *y = phi[ctx->k] - phi[ctx->k+1]*(PetscReal)ctx->k;
57: }
58: return(0);
59: }
61: PetscErrorCode FNEvaluateFunctionMatVec_Phi(FN fn,Mat A,Vec v)
62: {
64: FN_PHI *ctx = (FN_PHI*)fn->data;
65: PetscInt i,j,m,n,nh;
66: PetscScalar *Aa,*Ha,*Fa,*va,sfactor=1.0;
69: MatGetSize(A,&m,NULL);
70: n = m+ctx->k;
71: if (ctx->H) {
72: MatGetSize(ctx->H,&nh,NULL);
73: if (n!=nh) {
74: MatDestroy(&ctx->H);
75: MatDestroy(&ctx->F);
76: }
77: }
78: if (!ctx->H) {
79: MatCreateDense(PETSC_COMM_SELF,n,n,n,n,NULL,&ctx->H);
80: MatCreateDense(PETSC_COMM_SELF,n,n,n,n,NULL,&ctx->F);
81: }
82: MatDenseGetArray(ctx->H,&Ha);
83: MatDenseGetArray(A,&Aa);
84: for (j=0;j<m;j++) {
85: PetscMemcpy(Ha+j*n,Aa+j*m,m*sizeof(PetscScalar));
86: }
87: MatDenseRestoreArray(A,&Aa);
88: if (ctx->k) {
89: for (j=0;j<m;j++) for (i=m;i<n;i++) Ha[i+j*n] = 0.0;
90: for (j=m;j<n;j++) for (i=0;i<n;i++) Ha[i+j*n] = 0.0;
91: Ha[0+m*n] = fn->alpha;
92: for (j=m+1;j<n;j++) Ha[j-1+j*n] = fn->alpha;
93: }
94: MatDenseRestoreArray(ctx->H,&Ha);
96: FNEvaluateFunctionMat_Exp_Higham(fn,ctx->H,ctx->F);
98: MatDenseGetArray(ctx->F,&Fa);
99: VecGetArray(v,&va);
100: if (ctx->k) {
101: sfactor = PetscPowScalarInt(fn->alpha,-ctx->k);
102: for (i=0;i<m;i++) va[i] = sfactor*Fa[i+(n-1)*n];
103: } else {
104: for (i=0;i<m;i++) va[i] = sfactor*Fa[i+0*n];
105: }
106: VecRestoreArray(v,&va);
107: MatDenseRestoreArray(ctx->F,&Fa);
108: return(0);
109: }
111: static PetscErrorCode FNPhiSetIndex_Phi(FN fn,PetscInt k)
112: {
114: FN_PHI *ctx = (FN_PHI*)fn->data;
117: if (k<0) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Index cannot be negative");
118: if (k>MAX_INDEX) SETERRQ1(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Phi functions only implemented for k<=%d",MAX_INDEX);
119: if (k!=ctx->k) {
120: ctx->k = k;
121: MatDestroy(&ctx->H);
122: MatDestroy(&ctx->F);
123: }
124: return(0);
125: }
127: /*@
128: FNPhiSetIndex - Sets the index of the phi-function.
130: Logically Collective on FN
132: Input Parameters:
133: + fn - the math function context
134: - k - the index
136: Notes:
137: The phi-functions are defined as follows. The default is k=1.
138: .vb
139: phi_0(x) = exp(x)
140: phi_1(x) = (exp(x)-1)/x
141: phi_k(x) = (phi_{k-1}(x)-1/(k-1)!)/x
142: .ve
144: Level: intermediate
146: .seealso: FNPhiGetIndex()
147: @*/
148: PetscErrorCode FNPhiSetIndex(FN fn,PetscInt k)
149: {
155: PetscTryMethod(fn,"FNPhiSetIndex_C",(FN,PetscInt),(fn,k));
156: return(0);
157: }
159: static PetscErrorCode FNPhiGetIndex_Phi(FN fn,PetscInt *k)
160: {
161: FN_PHI *ctx = (FN_PHI*)fn->data;
164: *k = ctx->k;
165: return(0);
166: }
168: /*@
169: FNPhiGetIndex - Gets the index of the phi-function.
171: Not Collective
173: Input Parameter:
174: . fn - the math function context
176: Output Parameter:
177: . k - the index
179: Level: intermediate
181: .seealso: FNPhiSetIndex()
182: @*/
183: PetscErrorCode FNPhiGetIndex(FN fn,PetscInt *k)
184: {
190: PetscUseMethod(fn,"FNPhiGetIndex_C",(FN,PetscInt*),(fn,k));
191: return(0);
192: }
194: PetscErrorCode FNView_Phi(FN fn,PetscViewer viewer)
195: {
197: FN_PHI *ctx = (FN_PHI*)fn->data;
198: PetscBool isascii;
199: char str[50],strx[50];
202: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
203: if (isascii) {
204: PetscViewerASCIIPrintf(viewer," Phi_%D: ",ctx->k);
205: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
206: if (fn->beta!=(PetscScalar)1.0) {
207: SlepcSNPrintfScalar(str,50,fn->beta,PETSC_TRUE);
208: PetscViewerASCIIPrintf(viewer,"%s*",str);
209: }
210: if (fn->alpha==(PetscScalar)1.0) {
211: PetscSNPrintf(strx,50,"x");
212: } else {
213: SlepcSNPrintfScalar(str,50,fn->alpha,PETSC_TRUE);
214: PetscSNPrintf(strx,50,"(%s*x)",str);
215: }
216: if (!ctx->k) {
217: PetscViewerASCIIPrintf(viewer,"exp(%s)\n",strx);
218: } else if (ctx->k==1) {
219: PetscViewerASCIIPrintf(viewer,"(exp(%s)-1)/%s\n",strx,strx);
220: } else {
221: PetscViewerASCIIPrintf(viewer,"(phi_%D(%s)-1/%D!)/%s\n",ctx->k-1,strx,ctx->k-1,strx);
222: }
223: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
224: }
225: return(0);
226: }
228: PetscErrorCode FNSetFromOptions_Phi(PetscOptionItems *PetscOptionsObject,FN fn)
229: {
231: FN_PHI *ctx = (FN_PHI*)fn->data;
232: PetscInt k;
233: PetscBool flag;
236: PetscOptionsHead(PetscOptionsObject,"FN Phi Options");
238: PetscOptionsInt("-fn_phi_index","Index of the phi-function","FNPhiSetIndex",ctx->k,&k,&flag);
239: if (flag) { FNPhiSetIndex(fn,k); }
241: PetscOptionsTail();
242: return(0);
243: }
245: PetscErrorCode FNDuplicate_Phi(FN fn,MPI_Comm comm,FN *newfn)
246: {
247: FN_PHI *ctx = (FN_PHI*)fn->data,*ctx2 = (FN_PHI*)(*newfn)->data;
250: ctx2->k = ctx->k;
251: return(0);
252: }
254: PetscErrorCode FNDestroy_Phi(FN fn)
255: {
257: FN_PHI *ctx = (FN_PHI*)fn->data;
260: MatDestroy(&ctx->H);
261: MatDestroy(&ctx->F);
262: PetscFree(fn->data);
263: PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",NULL);
264: PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",NULL);
265: return(0);
266: }
268: SLEPC_EXTERN PetscErrorCode FNCreate_Phi(FN fn)
269: {
271: FN_PHI *ctx;
274: PetscNewLog(fn,&ctx);
275: fn->data = (void*)ctx;
276: ctx->k = 1;
278: fn->ops->evaluatefunction = FNEvaluateFunction_Phi;
279: fn->ops->evaluatederivative = FNEvaluateDerivative_Phi;
280: fn->ops->evaluatefunctionmatvec[0] = FNEvaluateFunctionMatVec_Phi;
281: fn->ops->setfromoptions = FNSetFromOptions_Phi;
282: fn->ops->view = FNView_Phi;
283: fn->ops->duplicate = FNDuplicate_Phi;
284: fn->ops->destroy = FNDestroy_Phi;
285: PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",FNPhiSetIndex_Phi);
286: PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",FNPhiGetIndex_Phi);
287: return(0);
288: }