Actual source code: fnphi.c
slepc-3.7.0 2016-05-16
1: /*
2: Phi functions
3: phi_0(x) = exp(x)
4: phi_1(x) = (exp(x)-1)/x
5: phi_k(x) = (phi_{k-1}(x)-1/(k-1)!)/x
7: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
8: SLEPc - Scalable Library for Eigenvalue Problem Computations
9: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
11: This file is part of SLEPc.
13: SLEPc is free software: you can redistribute it and/or modify it under the
14: terms of version 3 of the GNU Lesser General Public License as published by
15: the Free Software Foundation.
17: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
18: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
20: more details.
22: You should have received a copy of the GNU Lesser General Public License
23: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
24: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
25: */
27: #include <slepc/private/fnimpl.h> /*I "slepcfn.h" I*/
29: typedef struct {
30: PetscInt k; /* index of the phi-function, defaults to k=1 */
31: } FN_PHI;
33: const static PetscReal rfactorial[] = { 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 };
35: static void PhiFunction(PetscScalar x,PetscScalar *y,PetscInt k)
36: {
37: PetscScalar phi;
39: if (!k) *y = PetscExpScalar(x);
40: else if (k==1) *y = (PetscExpScalar(x)-1.0)/x;
41: else {
42: /* phi_k(x) = (phi_{k-1}(x)-1/(k-1)!)/x */
43: PhiFunction(x,&phi,k-1);
44: *y = (phi-rfactorial[k-1])/x;
45: }
46: }
50: PetscErrorCode FNEvaluateFunction_Phi(FN fn,PetscScalar x,PetscScalar *y)
51: {
52: FN_PHI *ctx = (FN_PHI*)fn->data;
55: PhiFunction(x,y,ctx->k);
56: return(0);
57: }
59: static void PhiDerivative(PetscScalar x,PetscScalar *y,PetscInt k)
60: {
61: PetscScalar der,phi;
63: if (!k) *y = PetscExpScalar(x);
64: else if (k==1) {
65: der = PetscExpScalar(x);
66: phi = (der-1.0)/x;
67: *y = (der-phi)/x;
68: } else {
69: PhiDerivative(x,&der,k-1);
70: PhiFunction(x,&phi,k);
71: *y = (der-phi)/x;
72: }
73: }
77: PetscErrorCode FNEvaluateDerivative_Phi(FN fn,PetscScalar x,PetscScalar *y)
78: {
79: FN_PHI *ctx = (FN_PHI*)fn->data;
82: PhiDerivative(x,y,ctx->k);
83: return(0);
84: }
88: static PetscErrorCode FNPhiSetIndex_Phi(FN fn,PetscInt k)
89: {
90: FN_PHI *ctx = (FN_PHI*)fn->data;
93: ctx->k = k;
94: return(0);
95: }
99: /*@
100: FNPhiSetIndex - Sets the index of the phi-function.
102: Logically Collective on FN
104: Input Parameters:
105: + fn - the math function context
106: - k - the index
108: Notes:
109: The phi-functions are defined as follows. The default is k=1.
110: .vb
111: phi_0(x) = exp(x)
112: phi_1(x) = (exp(x)-1)/x
113: phi_k(x) = (phi_{k-1}(x)-1/(k-1)!)/x
114: .ve
116: Level: intermediate
118: .seealso: FNPhiGetIndex()
119: @*/
120: PetscErrorCode FNPhiSetIndex(FN fn,PetscInt k)
121: {
127: if (k<0) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Index cannot be negative");
128: if (k>10) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Only implemented for k<=10");
129: PetscTryMethod(fn,"FNPhiSetIndex_C",(FN,PetscInt),(fn,k));
130: return(0);
131: }
135: static PetscErrorCode FNPhiGetIndex_Phi(FN fn,PetscInt *k)
136: {
137: FN_PHI *ctx = (FN_PHI*)fn->data;
140: *k = ctx->k;
141: return(0);
142: }
146: /*@
147: FNPhiGetIndex - Gets the index of the phi-function.
149: Not Collective
151: Input Parameter:
152: . fn - the math function context
154: Output Parameter:
155: . k - the index
157: Level: intermediate
159: .seealso: FNPhiSetIndex()
160: @*/
161: PetscErrorCode FNPhiGetIndex(FN fn,PetscInt *k)
162: {
168: PetscUseMethod(fn,"FNPhiGetIndex_C",(FN,PetscInt*),(fn,k));
169: return(0);
170: }
174: PetscErrorCode FNView_Phi(FN fn,PetscViewer viewer)
175: {
177: FN_PHI *ctx = (FN_PHI*)fn->data;
178: PetscBool isascii;
179: char str[50],strx[50];
182: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
183: if (isascii) {
184: PetscViewerASCIIPrintf(viewer," Phi_%D: ",ctx->k);
185: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
186: if (fn->beta!=(PetscScalar)1.0) {
187: SlepcSNPrintfScalar(str,50,fn->beta,PETSC_TRUE);
188: PetscViewerASCIIPrintf(viewer,"%s*",str);
189: }
190: if (fn->alpha==(PetscScalar)1.0) {
191: PetscSNPrintf(strx,50,"x");
192: } else {
193: SlepcSNPrintfScalar(str,50,fn->alpha,PETSC_TRUE);
194: PetscSNPrintf(strx,50,"(%s*x)",str);
195: }
196: if (!ctx->k) {
197: PetscViewerASCIIPrintf(viewer,"exp(%s)\n",strx);
198: } else if (ctx->k==1) {
199: PetscViewerASCIIPrintf(viewer,"(exp(%s)-1)/%s\n",strx,strx);
200: } else {
201: PetscViewerASCIIPrintf(viewer,"(phi_%D(%s)-1/%D!)/%s\n",ctx->k-1,strx,ctx->k-1,strx);
202: }
203: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
204: }
205: return(0);
206: }
210: PetscErrorCode FNSetFromOptions_Phi(PetscOptionItems *PetscOptionsObject,FN fn)
211: {
213: FN_PHI *ctx = (FN_PHI*)fn->data;
214: PetscInt k;
215: PetscBool flag;
218: PetscOptionsHead(PetscOptionsObject,"FN Phi Options");
219: PetscOptionsInt("-fn_phi_index","Index of the phi-function","FNPhiSetIndex",ctx->k,&k,&flag);
220: if (flag) {
221: FNPhiSetIndex(fn,k);
222: }
223: PetscOptionsTail();
224: return(0);
225: }
229: PetscErrorCode FNDuplicate_Phi(FN fn,MPI_Comm comm,FN *newfn)
230: {
232: FN_PHI *ctx = (FN_PHI*)fn->data,*ctx2;
235: PetscNewLog(*newfn,&ctx2);
236: (*newfn)->data = (void*)ctx2;
237: ctx2->k = ctx->k;
238: return(0);
239: }
243: PetscErrorCode FNDestroy_Phi(FN fn)
244: {
248: PetscFree(fn->data);
249: PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",NULL);
250: PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",NULL);
251: return(0);
252: }
256: PETSC_EXTERN PetscErrorCode FNCreate_Phi(FN fn)
257: {
259: FN_PHI *ctx;
262: PetscNewLog(fn,&ctx);
263: fn->data = (void*)ctx;
264: ctx->k = 1;
266: fn->ops->evaluatefunction = FNEvaluateFunction_Phi;
267: fn->ops->evaluatederivative = FNEvaluateDerivative_Phi;
268: fn->ops->setfromoptions = FNSetFromOptions_Phi;
269: fn->ops->view = FNView_Phi;
270: fn->ops->duplicate = FNDuplicate_Phi;
271: fn->ops->destroy = FNDestroy_Phi;
272: PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",FNPhiSetIndex_Phi);
273: PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",FNPhiGetIndex_Phi);
274: return(0);
275: }