Actual source code: fnphi.c

slepc-3.18.0 2022-10-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: static const 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];

 35:   if (x==0.0) *y = rfactorial[ctx->k];
 36:   else {
 37:     phi[0] = PetscExpScalar(x);
 38:     for (i=1;i<=ctx->k;i++) phi[i] = (phi[i-1]-rfactorial[i-1])/x;
 39:     *y = phi[ctx->k];
 40:   }
 41:   return 0;
 42: }

 44: PetscErrorCode FNEvaluateDerivative_Phi(FN fn,PetscScalar x,PetscScalar *y)
 45: {
 46:   FN_PHI      *ctx = (FN_PHI*)fn->data;
 47:   PetscInt    i;
 48:   PetscScalar phi[MAX_INDEX+2];

 50:   if (x==0.0) *y = rfactorial[ctx->k+1];
 51:   else {
 52:     phi[0] = PetscExpScalar(x);
 53:     for (i=1;i<=ctx->k+1;i++) phi[i] = (phi[i-1]-rfactorial[i-1])/x;
 54:     *y = phi[ctx->k] - phi[ctx->k+1]*(PetscReal)ctx->k;
 55:   }
 56:   return 0;
 57: }

 59: PetscErrorCode FNEvaluateFunctionMatVec_Phi(FN fn,Mat A,Vec v)
 60: {
 61:   FN_PHI            *ctx = (FN_PHI*)fn->data;
 62:   PetscInt          i,j,m,n,nh;
 63:   PetscScalar       *Ha,*va,sfactor=1.0;
 64:   const PetscScalar *Aa,*Fa;

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

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

 93:   MatDenseGetArrayRead(ctx->F,&Fa);
 94:   VecGetArray(v,&va);
 95:   if (ctx->k) {
 96:     sfactor = PetscPowScalarInt(fn->alpha,-ctx->k);
 97:     for (i=0;i<m;i++) va[i] = sfactor*Fa[i+(n-1)*n];
 98:   } else {
 99:     for (i=0;i<m;i++) va[i] = sfactor*Fa[i+0*n];
100:   }
101:   VecRestoreArray(v,&va);
102:   MatDenseRestoreArrayRead(ctx->F,&Fa);
103:   return 0;
104: }

106: static PetscErrorCode FNPhiSetIndex_Phi(FN fn,PetscInt k)
107: {
108:   FN_PHI         *ctx = (FN_PHI*)fn->data;

112:   if (k!=ctx->k) {
113:     ctx->k = k;
114:     MatDestroy(&ctx->H);
115:     MatDestroy(&ctx->F);
116:   }
117:   return 0;
118: }

120: /*@
121:    FNPhiSetIndex - Sets the index of the phi-function.

123:    Logically Collective on fn

125:    Input Parameters:
126: +  fn - the math function context
127: -  k  - the index

129:    Notes:
130:    The phi-functions are defined as follows. The default is k=1.
131: .vb
132:       phi_0(x) = exp(x)
133:       phi_1(x) = (exp(x)-1)/x
134:       phi_k(x) = (phi_{k-1}(x)-1/(k-1)!)/x
135: .ve

137:    Level: intermediate

139: .seealso: FNPhiGetIndex()
140: @*/
141: PetscErrorCode FNPhiSetIndex(FN fn,PetscInt k)
142: {
145:   PetscTryMethod(fn,"FNPhiSetIndex_C",(FN,PetscInt),(fn,k));
146:   return 0;
147: }

149: static PetscErrorCode FNPhiGetIndex_Phi(FN fn,PetscInt *k)
150: {
151:   FN_PHI *ctx = (FN_PHI*)fn->data;

153:   *k = ctx->k;
154:   return 0;
155: }

157: /*@
158:    FNPhiGetIndex - Gets the index of the phi-function.

160:    Not Collective

162:    Input Parameter:
163: .  fn - the math function context

165:    Output Parameter:
166: .  k  - the index

168:    Level: intermediate

170: .seealso: FNPhiSetIndex()
171: @*/
172: PetscErrorCode FNPhiGetIndex(FN fn,PetscInt *k)
173: {
176:   PetscUseMethod(fn,"FNPhiGetIndex_C",(FN,PetscInt*),(fn,k));
177:   return 0;
178: }

180: PetscErrorCode FNView_Phi(FN fn,PetscViewer viewer)
181: {
182:   FN_PHI         *ctx = (FN_PHI*)fn->data;
183:   PetscBool      isascii;
184:   char           str[50],strx[50];

186:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
187:   if (isascii) {
188:     PetscViewerASCIIPrintf(viewer,"  phi_%" PetscInt_FMT ": ",ctx->k);
189:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
190:     if (fn->beta!=(PetscScalar)1.0) {
191:       SlepcSNPrintfScalar(str,sizeof(str),fn->beta,PETSC_TRUE);
192:       PetscViewerASCIIPrintf(viewer,"%s*",str);
193:     }
194:     if (fn->alpha==(PetscScalar)1.0) PetscSNPrintf(strx,sizeof(strx),"x");
195:     else {
196:       SlepcSNPrintfScalar(str,sizeof(str),fn->alpha,PETSC_TRUE);
197:       PetscSNPrintf(strx,sizeof(strx),"(%s*x)",str);
198:     }
199:     if (!ctx->k) PetscViewerASCIIPrintf(viewer,"exp(%s)\n",strx);
200:     else if (ctx->k==1) PetscViewerASCIIPrintf(viewer,"(exp(%s)-1)/%s\n",strx,strx);
201:     else PetscViewerASCIIPrintf(viewer,"(phi_%" PetscInt_FMT "(%s)-1/%" PetscInt_FMT "!)/%s\n",ctx->k-1,strx,ctx->k-1,strx);
202:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
203:   }
204:   return 0;
205: }

207: PetscErrorCode FNSetFromOptions_Phi(FN fn,PetscOptionItems *PetscOptionsObject)
208: {
209:   FN_PHI         *ctx = (FN_PHI*)fn->data;
210:   PetscInt       k;
211:   PetscBool      flag;

213:   PetscOptionsHeadBegin(PetscOptionsObject,"FN Phi Options");

215:     PetscOptionsInt("-fn_phi_index","Index of the phi-function","FNPhiSetIndex",ctx->k,&k,&flag);
216:     if (flag) FNPhiSetIndex(fn,k);

218:   PetscOptionsHeadEnd();
219:   return 0;
220: }

222: PetscErrorCode FNDuplicate_Phi(FN fn,MPI_Comm comm,FN *newfn)
223: {
224:   FN_PHI *ctx = (FN_PHI*)fn->data,*ctx2 = (FN_PHI*)(*newfn)->data;

226:   ctx2->k = ctx->k;
227:   return 0;
228: }

230: PetscErrorCode FNDestroy_Phi(FN fn)
231: {
232:   FN_PHI         *ctx = (FN_PHI*)fn->data;

234:   MatDestroy(&ctx->H);
235:   MatDestroy(&ctx->F);
236:   PetscFree(fn->data);
237:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",NULL);
238:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",NULL);
239:   return 0;
240: }

242: SLEPC_EXTERN PetscErrorCode FNCreate_Phi(FN fn)
243: {
244:   FN_PHI         *ctx;

246:   PetscNew(&ctx);
247:   fn->data = (void*)ctx;
248:   ctx->k   = 1;

250:   fn->ops->evaluatefunction          = FNEvaluateFunction_Phi;
251:   fn->ops->evaluatederivative        = FNEvaluateDerivative_Phi;
252:   fn->ops->evaluatefunctionmatvec[0] = FNEvaluateFunctionMatVec_Phi;
253:   fn->ops->setfromoptions            = FNSetFromOptions_Phi;
254:   fn->ops->view                      = FNView_Phi;
255:   fn->ops->duplicate                 = FNDuplicate_Phi;
256:   fn->ops->destroy                   = FNDestroy_Phi;
257:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",FNPhiSetIndex_Phi);
258:   PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",FNPhiGetIndex_Phi);
259:   return 0;
260: }