Actual source code: fnphi.c

slepc-3.15.1 2021-05-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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>

 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: {
 63:   PetscErrorCode    ierr;
 64:   FN_PHI            *ctx = (FN_PHI*)fn->data;
 65:   PetscInt          i,j,m,n,nh;
 66:   PetscScalar       *Ha,*va,sfactor=1.0;
 67:   const PetscScalar *Aa,*Fa;

 70:   MatGetSize(A,&m,NULL);
 71:   n = m+ctx->k;
 72:   if (ctx->H) {
 73:     MatGetSize(ctx->H,&nh,NULL);
 74:     if (n!=nh) {
 75:       MatDestroy(&ctx->H);
 76:       MatDestroy(&ctx->F);
 77:     }
 78:   }
 79:   if (!ctx->H) {
 80:     MatCreateDense(PETSC_COMM_SELF,n,n,n,n,NULL,&ctx->H);
 81:     MatCreateDense(PETSC_COMM_SELF,n,n,n,n,NULL,&ctx->F);
 82:   }
 83:   MatDenseGetArray(ctx->H,&Ha);
 84:   MatDenseGetArrayRead(A,&Aa);
 85:   for (j=0;j<m;j++) {
 86:     PetscArraycpy(Ha+j*n,Aa+j*m,m);
 87:   }
 88:   MatDenseRestoreArrayRead(A,&Aa);
 89:   if (ctx->k) {
 90:     for (j=0;j<m;j++) for (i=m;i<n;i++) Ha[i+j*n] = 0.0;
 91:     for (j=m;j<n;j++) for (i=0;i<n;i++) Ha[i+j*n] = 0.0;
 92:     Ha[0+m*n] = fn->alpha;
 93:     for (j=m+1;j<n;j++) Ha[j-1+j*n] = fn->alpha;
 94:   }
 95:   MatDenseRestoreArray(ctx->H,&Ha);

 97:   FNEvaluateFunctionMat_Exp_Higham(fn,ctx->H,ctx->F);

 99:   MatDenseGetArrayRead(ctx->F,&Fa);
100:   VecGetArray(v,&va);
101:   if (ctx->k) {
102:     sfactor = PetscPowScalarInt(fn->alpha,-ctx->k);
103:     for (i=0;i<m;i++) va[i] = sfactor*Fa[i+(n-1)*n];
104:   } else {
105:     for (i=0;i<m;i++) va[i] = sfactor*Fa[i+0*n];
106:   }
107:   VecRestoreArray(v,&va);
108:   MatDenseRestoreArrayRead(ctx->F,&Fa);
109:   return(0);
110: }

112: static PetscErrorCode FNPhiSetIndex_Phi(FN fn,PetscInt k)
113: {
115:   FN_PHI         *ctx = (FN_PHI*)fn->data;

118:   if (k<0) SETERRQ(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Index cannot be negative");
119:   if (k>MAX_INDEX) SETERRQ1(PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Phi functions only implemented for k<=%d",MAX_INDEX);
120:   if (k!=ctx->k) {
121:     ctx->k = k;
122:     MatDestroy(&ctx->H);
123:     MatDestroy(&ctx->F);
124:   }
125:   return(0);
126: }

128: /*@
129:    FNPhiSetIndex - Sets the index of the phi-function.

131:    Logically Collective on fn

133:    Input Parameters:
134: +  fn - the math function context
135: -  k  - the index

137:    Notes:
138:    The phi-functions are defined as follows. The default is k=1.
139: .vb
140:       phi_0(x) = exp(x)
141:       phi_1(x) = (exp(x)-1)/x
142:       phi_k(x) = (phi_{k-1}(x)-1/(k-1)!)/x
143: .ve

145:    Level: intermediate

147: .seealso: FNPhiGetIndex()
148: @*/
149: PetscErrorCode FNPhiSetIndex(FN fn,PetscInt k)
150: {

156:   PetscTryMethod(fn,"FNPhiSetIndex_C",(FN,PetscInt),(fn,k));
157:   return(0);
158: }

160: static PetscErrorCode FNPhiGetIndex_Phi(FN fn,PetscInt *k)
161: {
162:   FN_PHI *ctx = (FN_PHI*)fn->data;

165:   *k = ctx->k;
166:   return(0);
167: }

169: /*@
170:    FNPhiGetIndex - Gets the index of the phi-function.

172:    Not Collective

174:    Input Parameter:
175: .  fn - the math function context

177:    Output Parameter:
178: .  k  - the index

180:    Level: intermediate

182: .seealso: FNPhiSetIndex()
183: @*/
184: PetscErrorCode FNPhiGetIndex(FN fn,PetscInt *k)
185: {

191:   PetscUseMethod(fn,"FNPhiGetIndex_C",(FN,PetscInt*),(fn,k));
192:   return(0);
193: }

195: PetscErrorCode FNView_Phi(FN fn,PetscViewer viewer)
196: {
198:   FN_PHI         *ctx = (FN_PHI*)fn->data;
199:   PetscBool      isascii;
200:   char           str[50],strx[50];

203:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
204:   if (isascii) {
205:     PetscViewerASCIIPrintf(viewer,"  Phi_%D: ",ctx->k);
206:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
207:     if (fn->beta!=(PetscScalar)1.0) {
208:       SlepcSNPrintfScalar(str,sizeof(str),fn->beta,PETSC_TRUE);
209:       PetscViewerASCIIPrintf(viewer,"%s*",str);
210:     }
211:     if (fn->alpha==(PetscScalar)1.0) {
212:       PetscSNPrintf(strx,sizeof(strx),"x");
213:     } else {
214:       SlepcSNPrintfScalar(str,sizeof(str),fn->alpha,PETSC_TRUE);
215:       PetscSNPrintf(strx,sizeof(strx),"(%s*x)",str);
216:     }
217:     if (!ctx->k) {
218:       PetscViewerASCIIPrintf(viewer,"exp(%s)\n",strx);
219:     } else if (ctx->k==1) {
220:       PetscViewerASCIIPrintf(viewer,"(exp(%s)-1)/%s\n",strx,strx);
221:     } else {
222:       PetscViewerASCIIPrintf(viewer,"(phi_%D(%s)-1/%D!)/%s\n",ctx->k-1,strx,ctx->k-1,strx);
223:     }
224:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
225:   }
226:   return(0);
227: }

229: PetscErrorCode FNSetFromOptions_Phi(PetscOptionItems *PetscOptionsObject,FN fn)
230: {
232:   FN_PHI         *ctx = (FN_PHI*)fn->data;
233:   PetscInt       k;
234:   PetscBool      flag;

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

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

242:   PetscOptionsTail();
243:   return(0);
244: }

246: PetscErrorCode FNDuplicate_Phi(FN fn,MPI_Comm comm,FN *newfn)
247: {
248:   FN_PHI *ctx = (FN_PHI*)fn->data,*ctx2 = (FN_PHI*)(*newfn)->data;

251:   ctx2->k = ctx->k;
252:   return(0);
253: }

255: PetscErrorCode FNDestroy_Phi(FN fn)
256: {
258:   FN_PHI         *ctx = (FN_PHI*)fn->data;

261:   MatDestroy(&ctx->H);
262:   MatDestroy(&ctx->F);
263:   PetscFree(fn->data);
264:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",NULL);
265:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",NULL);
266:   return(0);
267: }

269: SLEPC_EXTERN PetscErrorCode FNCreate_Phi(FN fn)
270: {
272:   FN_PHI         *ctx;

275:   PetscNewLog(fn,&ctx);
276:   fn->data = (void*)ctx;
277:   ctx->k   = 1;

279:   fn->ops->evaluatefunction          = FNEvaluateFunction_Phi;
280:   fn->ops->evaluatederivative        = FNEvaluateDerivative_Phi;
281:   fn->ops->evaluatefunctionmatvec[0] = FNEvaluateFunctionMatVec_Phi;
282:   fn->ops->setfromoptions            = FNSetFromOptions_Phi;
283:   fn->ops->view                      = FNView_Phi;
284:   fn->ops->duplicate                 = FNDuplicate_Phi;
285:   fn->ops->destroy                   = FNDestroy_Phi;
286:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",FNPhiSetIndex_Phi);
287:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",FNPhiGetIndex_Phi);
288:   return(0);
289: }