Actual source code: fnphi.c

slepc-3.22.1 2024-10-28
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: static 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:   PetscFunctionBegin;
 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:   PetscFunctionReturn(PETSC_SUCCESS);
 43: }

 45: static 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];

 51:   PetscFunctionBegin;
 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:   PetscFunctionReturn(PETSC_SUCCESS);
 59: }

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

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

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

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

109: static PetscErrorCode FNPhiSetIndex_Phi(FN fn,PetscInt k)
110: {
111:   FN_PHI         *ctx = (FN_PHI*)fn->data;

113:   PetscFunctionBegin;
114:   PetscCheck(k>=0,PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Index cannot be negative");
115:   PetscCheck(k<=MAX_INDEX,PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Phi functions only implemented for k<=%d",MAX_INDEX);
116:   if (k!=ctx->k) {
117:     ctx->k = k;
118:     PetscCall(MatDestroy(&ctx->H));
119:     PetscCall(MatDestroy(&ctx->F));
120:   }
121:   PetscFunctionReturn(PETSC_SUCCESS);
122: }

124: /*@
125:    FNPhiSetIndex - Sets the index of the phi-function.

127:    Logically Collective

129:    Input Parameters:
130: +  fn - the math function context
131: -  k  - the index

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

141:    Level: intermediate

143: .seealso: FNPhiGetIndex()
144: @*/
145: PetscErrorCode FNPhiSetIndex(FN fn,PetscInt k)
146: {
147:   PetscFunctionBegin;
150:   PetscTryMethod(fn,"FNPhiSetIndex_C",(FN,PetscInt),(fn,k));
151:   PetscFunctionReturn(PETSC_SUCCESS);
152: }

154: static PetscErrorCode FNPhiGetIndex_Phi(FN fn,PetscInt *k)
155: {
156:   FN_PHI *ctx = (FN_PHI*)fn->data;

158:   PetscFunctionBegin;
159:   *k = ctx->k;
160:   PetscFunctionReturn(PETSC_SUCCESS);
161: }

163: /*@
164:    FNPhiGetIndex - Gets the index of the phi-function.

166:    Not Collective

168:    Input Parameter:
169: .  fn - the math function context

171:    Output Parameter:
172: .  k  - the index

174:    Level: intermediate

176: .seealso: FNPhiSetIndex()
177: @*/
178: PetscErrorCode FNPhiGetIndex(FN fn,PetscInt *k)
179: {
180:   PetscFunctionBegin;
182:   PetscAssertPointer(k,2);
183:   PetscUseMethod(fn,"FNPhiGetIndex_C",(FN,PetscInt*),(fn,k));
184:   PetscFunctionReturn(PETSC_SUCCESS);
185: }

187: static PetscErrorCode FNView_Phi(FN fn,PetscViewer viewer)
188: {
189:   FN_PHI         *ctx = (FN_PHI*)fn->data;
190:   PetscBool      isascii;
191:   char           str[50],strx[50];

193:   PetscFunctionBegin;
194:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
195:   if (isascii) {
196:     PetscCall(PetscViewerASCIIPrintf(viewer,"  phi_%" PetscInt_FMT ": ",ctx->k));
197:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
198:     if (fn->beta!=(PetscScalar)1.0) {
199:       PetscCall(SlepcSNPrintfScalar(str,sizeof(str),fn->beta,PETSC_TRUE));
200:       PetscCall(PetscViewerASCIIPrintf(viewer,"%s*",str));
201:     }
202:     if (fn->alpha==(PetscScalar)1.0) PetscCall(PetscSNPrintf(strx,sizeof(strx),"x"));
203:     else {
204:       PetscCall(SlepcSNPrintfScalar(str,sizeof(str),fn->alpha,PETSC_TRUE));
205:       PetscCall(PetscSNPrintf(strx,sizeof(strx),"(%s*x)",str));
206:     }
207:     if (!ctx->k) PetscCall(PetscViewerASCIIPrintf(viewer,"exp(%s)\n",strx));
208:     else if (ctx->k==1) PetscCall(PetscViewerASCIIPrintf(viewer,"(exp(%s)-1)/%s\n",strx,strx));
209:     else PetscCall(PetscViewerASCIIPrintf(viewer,"(phi_%" PetscInt_FMT "(%s)-1/%" PetscInt_FMT "!)/%s\n",ctx->k-1,strx,ctx->k-1,strx));
210:     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
211:   }
212:   PetscFunctionReturn(PETSC_SUCCESS);
213: }

215: static PetscErrorCode FNSetFromOptions_Phi(FN fn,PetscOptionItems *PetscOptionsObject)
216: {
217:   FN_PHI         *ctx = (FN_PHI*)fn->data;
218:   PetscInt       k;
219:   PetscBool      flag;

221:   PetscFunctionBegin;
222:   PetscOptionsHeadBegin(PetscOptionsObject,"FN Phi Options");

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

227:   PetscOptionsHeadEnd();
228:   PetscFunctionReturn(PETSC_SUCCESS);
229: }

231: static PetscErrorCode FNDuplicate_Phi(FN fn,MPI_Comm comm,FN *newfn)
232: {
233:   FN_PHI *ctx = (FN_PHI*)fn->data,*ctx2 = (FN_PHI*)(*newfn)->data;

235:   PetscFunctionBegin;
236:   ctx2->k = ctx->k;
237:   PetscFunctionReturn(PETSC_SUCCESS);
238: }

240: static PetscErrorCode FNDestroy_Phi(FN fn)
241: {
242:   FN_PHI         *ctx = (FN_PHI*)fn->data;

244:   PetscFunctionBegin;
245:   PetscCall(MatDestroy(&ctx->H));
246:   PetscCall(MatDestroy(&ctx->F));
247:   PetscCall(PetscFree(fn->data));
248:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",NULL));
249:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",NULL));
250:   PetscFunctionReturn(PETSC_SUCCESS);
251: }

253: SLEPC_EXTERN PetscErrorCode FNCreate_Phi(FN fn)
254: {
255:   FN_PHI         *ctx;

257:   PetscFunctionBegin;
258:   PetscCall(PetscNew(&ctx));
259:   fn->data = (void*)ctx;
260:   ctx->k   = 1;

262:   fn->ops->evaluatefunction          = FNEvaluateFunction_Phi;
263:   fn->ops->evaluatederivative        = FNEvaluateDerivative_Phi;
264:   fn->ops->evaluatefunctionmatvec[0] = FNEvaluateFunctionMatVec_Phi;
265:   fn->ops->setfromoptions            = FNSetFromOptions_Phi;
266:   fn->ops->view                      = FNView_Phi;
267:   fn->ops->duplicate                 = FNDuplicate_Phi;
268:   fn->ops->destroy                   = FNDestroy_Phi;
269:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNPhiSetIndex_C",FNPhiSetIndex_Phi));
270:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNPhiGetIndex_C",FNPhiGetIndex_Phi));
271:   PetscFunctionReturn(PETSC_SUCCESS);
272: }