Actual source code: fnphi.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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: } FN_PHI;

 23: 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 };

 25: static void PhiFunction(PetscScalar x,PetscScalar *y,PetscInt k)
 26: {
 27:   PetscScalar phi;

 29:   if (!k) *y = PetscExpScalar(x);
 30:   else if (k==1) *y = (PetscExpScalar(x)-1.0)/x;
 31:   else {
 32:     /* phi_k(x) = (phi_{k-1}(x)-1/(k-1)!)/x */
 33:     PhiFunction(x,&phi,k-1);
 34:     *y = (phi-rfactorial[k-1])/x;
 35:   }
 36: }

 38: PetscErrorCode FNEvaluateFunction_Phi(FN fn,PetscScalar x,PetscScalar *y)
 39: {
 40:   FN_PHI *ctx = (FN_PHI*)fn->data;

 43:   PhiFunction(x,y,ctx->k);
 44:   return(0);
 45: }

 47: static void PhiDerivative(PetscScalar x,PetscScalar *y,PetscInt k)
 48: {
 49:   PetscScalar der,phi;

 51:   if (!k) *y = PetscExpScalar(x);
 52:   else if (k==1) {
 53:     der = PetscExpScalar(x);
 54:     phi = (der-1.0)/x;
 55:     *y = (der-phi)/x;
 56:   } else {
 57:     PhiDerivative(x,&der,k-1);
 58:     PhiFunction(x,&phi,k);
 59:     *y = (der-phi)/x;
 60:   }
 61: }

 63: PetscErrorCode FNEvaluateDerivative_Phi(FN fn,PetscScalar x,PetscScalar *y)
 64: {
 65:   FN_PHI *ctx = (FN_PHI*)fn->data;

 68:   PhiDerivative(x,y,ctx->k);
 69:   return(0);
 70: }

 72: static PetscErrorCode FNPhiSetIndex_Phi(FN fn,PetscInt k)
 73: {
 74:   FN_PHI *ctx = (FN_PHI*)fn->data;

 77:   if (k<0) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Index cannot be negative");
 78:   if (k>10) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Only implemented for k<=10");
 79:   ctx->k = k;
 80:   return(0);
 81: }

 83: /*@
 84:    FNPhiSetIndex - Sets the index of the phi-function.

 86:    Logically Collective on FN

 88:    Input Parameters:
 89: +  fn - the math function context
 90: -  k  - the index

 92:    Notes:
 93:    The phi-functions are defined as follows. The default is k=1.
 94: .vb
 95:       phi_0(x) = exp(x)
 96:       phi_1(x) = (exp(x)-1)/x
 97:       phi_k(x) = (phi_{k-1}(x)-1/(k-1)!)/x
 98: .ve

100:    Level: intermediate

102: .seealso: FNPhiGetIndex()
103: @*/
104: PetscErrorCode FNPhiSetIndex(FN fn,PetscInt k)
105: {

111:   PetscTryMethod(fn,"FNPhiSetIndex_C",(FN,PetscInt),(fn,k));
112:   return(0);
113: }

115: static PetscErrorCode FNPhiGetIndex_Phi(FN fn,PetscInt *k)
116: {
117:   FN_PHI *ctx = (FN_PHI*)fn->data;

120:   *k = ctx->k;
121:   return(0);
122: }

124: /*@
125:    FNPhiGetIndex - Gets the index of the phi-function.

127:    Not Collective

129:    Input Parameter:
130: .  fn - the math function context

132:    Output Parameter:
133: .  k  - the index

135:    Level: intermediate

137: .seealso: FNPhiSetIndex()
138: @*/
139: PetscErrorCode FNPhiGetIndex(FN fn,PetscInt *k)
140: {

146:   PetscUseMethod(fn,"FNPhiGetIndex_C",(FN,PetscInt*),(fn,k));
147:   return(0);
148: }

150: PetscErrorCode FNView_Phi(FN fn,PetscViewer viewer)
151: {
153:   FN_PHI         *ctx = (FN_PHI*)fn->data;
154:   PetscBool      isascii;
155:   char           str[50],strx[50];

158:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
159:   if (isascii) {
160:     PetscViewerASCIIPrintf(viewer,"  Phi_%D: ",ctx->k);
161:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
162:     if (fn->beta!=(PetscScalar)1.0) {
163:       SlepcSNPrintfScalar(str,50,fn->beta,PETSC_TRUE);
164:       PetscViewerASCIIPrintf(viewer,"%s*",str);
165:     }
166:     if (fn->alpha==(PetscScalar)1.0) {
167:       PetscSNPrintf(strx,50,"x");
168:     } else {
169:       SlepcSNPrintfScalar(str,50,fn->alpha,PETSC_TRUE);
170:       PetscSNPrintf(strx,50,"(%s*x)",str);
171:     }
172:     if (!ctx->k) {
173:       PetscViewerASCIIPrintf(viewer,"exp(%s)\n",strx);
174:     } else if (ctx->k==1) {
175:       PetscViewerASCIIPrintf(viewer,"(exp(%s)-1)/%s\n",strx,strx);
176:     } else {
177:       PetscViewerASCIIPrintf(viewer,"(phi_%D(%s)-1/%D!)/%s\n",ctx->k-1,strx,ctx->k-1,strx);
178:     }
179:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
180:   }
181:   return(0);
182: }

184: PetscErrorCode FNSetFromOptions_Phi(PetscOptionItems *PetscOptionsObject,FN fn)
185: {
187:   FN_PHI         *ctx = (FN_PHI*)fn->data;
188:   PetscInt       k;
189:   PetscBool      flag;

192:   PetscOptionsHead(PetscOptionsObject,"FN Phi Options");

194:     PetscOptionsInt("-fn_phi_index","Index of the phi-function","FNPhiSetIndex",ctx->k,&k,&flag);
195:     if (flag) { FNPhiSetIndex(fn,k); }

197:   PetscOptionsTail();
198:   return(0);
199: }

201: PetscErrorCode FNDuplicate_Phi(FN fn,MPI_Comm comm,FN *newfn)
202: {
203:   FN_PHI *ctx = (FN_PHI*)fn->data,*ctx2 = (FN_PHI*)(*newfn)->data;

206:   ctx2->k = ctx->k;
207:   return(0);
208: }

210: PetscErrorCode FNDestroy_Phi(FN fn)
211: {

215:   PetscFree(fn->data);
216:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",NULL);
217:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",NULL);
218:   return(0);
219: }

221: PETSC_EXTERN PetscErrorCode FNCreate_Phi(FN fn)
222: {
224:   FN_PHI         *ctx;

227:   PetscNewLog(fn,&ctx);
228:   fn->data = (void*)ctx;
229:   ctx->k   = 1;

231:   fn->ops->evaluatefunction    = FNEvaluateFunction_Phi;
232:   fn->ops->evaluatederivative  = FNEvaluateDerivative_Phi;
233:   fn->ops->setfromoptions      = FNSetFromOptions_Phi;
234:   fn->ops->view                = FNView_Phi;
235:   fn->ops->duplicate           = FNDuplicate_Phi;
236:   fn->ops->destroy             = FNDestroy_Phi;
237:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",FNPhiSetIndex_Phi);
238:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",FNPhiGetIndex_Phi);
239:   return(0);
240: }