Actual source code: fnrational.c

slepc-main 2024-11-09
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:    Rational function  r(x) = p(x)/q(x), where p(x) and q(x) are polynomials
 12: */

 14: #include <slepc/private/fnimpl.h>

 16: typedef struct {
 17:   PetscScalar *pcoeff;    /* numerator coefficients */
 18:   PetscInt    np;         /* length of array pcoeff, p(x) has degree np-1 */
 19:   PetscScalar *qcoeff;    /* denominator coefficients */
 20:   PetscInt    nq;         /* length of array qcoeff, q(x) has degree nq-1 */
 21: } FN_RATIONAL;

 23: static PetscErrorCode FNEvaluateFunction_Rational(FN fn,PetscScalar x,PetscScalar *y)
 24: {
 25:   FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
 26:   PetscInt    i;
 27:   PetscScalar p,q;

 29:   PetscFunctionBegin;
 30:   if (!ctx->np) p = 1.0;
 31:   else {
 32:     p = ctx->pcoeff[0];
 33:     for (i=1;i<ctx->np;i++)
 34:       p = ctx->pcoeff[i]+x*p;
 35:   }
 36:   if (!ctx->nq) *y = p;
 37:   else {
 38:     q = ctx->qcoeff[0];
 39:     for (i=1;i<ctx->nq;i++)
 40:       q = ctx->qcoeff[i]+x*q;
 41:     PetscCheck(q!=0.0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Function not defined in the requested value");
 42:     *y = p/q;
 43:   }
 44:   PetscFunctionReturn(PETSC_SUCCESS);
 45: }

 47: /*
 48:    Horner evaluation of P=p(A)
 49:    d = degree of polynomial;   coeff = coefficients of polynomial;    W = workspace
 50: */
 51: static PetscErrorCode EvaluatePoly(Mat A,Mat P,Mat W,PetscInt d,PetscScalar *coeff)
 52: {
 53:   PetscInt j;

 55:   PetscFunctionBegin;
 56:   PetscCall(MatZeroEntries(P));
 57:   if (!d) PetscCall(MatShift(P,1.0));
 58:   else {
 59:     PetscCall(MatShift(P,coeff[0]));
 60:     for (j=1;j<d;j++) {
 61:       PetscCall(MatMatMult(P,A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&W));
 62:       PetscCall(MatCopy(W,P,SAME_NONZERO_PATTERN));
 63:       PetscCall(MatShift(P,coeff[j]));
 64:     }
 65:   }
 66:   PetscFunctionReturn(PETSC_SUCCESS);
 67: }

 69: static PetscErrorCode FNEvaluateFunctionMat_Rational(FN fn,Mat A,Mat B)
 70: {
 71:   FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
 72:   Mat         P,Q,W,F;
 73:   PetscBool   iscuda;

 75:   PetscFunctionBegin;
 76:   if (A==B) PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&P));
 77:   else P = B;
 78:   PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&W));

 80:   PetscCall(EvaluatePoly(A,P,W,ctx->np,ctx->pcoeff));
 81:   if (ctx->nq) {
 82:     PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q));
 83:     PetscCall(EvaluatePoly(A,Q,W,ctx->nq,ctx->qcoeff));
 84:     PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSECUDA,&iscuda));
 85:     PetscCall(MatGetFactor(Q,iscuda?MATSOLVERCUDA:MATSOLVERPETSC,MAT_FACTOR_LU,&F));
 86:     PetscCall(MatLUFactorSymbolic(F,Q,NULL,NULL,NULL));
 87:     PetscCall(MatLUFactorNumeric(F,Q,NULL));
 88:     PetscCall(MatMatSolve(F,P,P));
 89:     PetscCall(MatDestroy(&F));
 90:     PetscCall(MatDestroy(&Q));
 91:   }

 93:   if (A==B) {
 94:     PetscCall(MatCopy(P,B,SAME_NONZERO_PATTERN));
 95:     PetscCall(MatDestroy(&P));
 96:   }
 97:   PetscCall(MatDestroy(&W));
 98:   PetscFunctionReturn(PETSC_SUCCESS);
 99: }

101: static PetscErrorCode FNEvaluateFunctionMatVec_Rational(FN fn,Mat A,Vec v)
102: {
103:   FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
104:   Mat         P,Q,W,F;
105:   Vec         b;
106:   PetscBool   iscuda;

108:   PetscFunctionBegin;
109:   PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&P));
110:   PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&W));

112:   PetscCall(EvaluatePoly(A,P,W,ctx->np,ctx->pcoeff));
113:   if (ctx->nq) {
114:     PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q));
115:     PetscCall(EvaluatePoly(A,Q,W,ctx->nq,ctx->qcoeff));
116:     PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSECUDA,&iscuda));
117:     PetscCall(MatGetFactor(Q,iscuda?MATSOLVERCUDA:MATSOLVERPETSC,MAT_FACTOR_LU,&F));
118:     PetscCall(MatLUFactorSymbolic(F,Q,NULL,NULL,NULL));
119:     PetscCall(MatLUFactorNumeric(F,Q,NULL));
120:     PetscCall(MatCreateVecs(P,&b,NULL));
121:     PetscCall(MatGetColumnVector(P,b,0));
122:     PetscCall(MatSolve(F,b,v));
123:     PetscCall(VecDestroy(&b));
124:     PetscCall(MatDestroy(&F));
125:     PetscCall(MatDestroy(&Q));
126:   } else PetscCall(MatGetColumnVector(P,v,0));

128:   PetscCall(MatDestroy(&P));
129:   PetscCall(MatDestroy(&W));
130:   PetscFunctionReturn(PETSC_SUCCESS);
131: }

133: static PetscErrorCode FNEvaluateDerivative_Rational(FN fn,PetscScalar x,PetscScalar *yp)
134: {
135:   FN_RATIONAL *ctx = (FN_RATIONAL*)fn->data;
136:   PetscInt    i;
137:   PetscScalar p,q,pp,qp;

139:   PetscFunctionBegin;
140:   if (!ctx->np) {
141:     p = 1.0;
142:     pp = 0.0;
143:   } else {
144:     p = ctx->pcoeff[0];
145:     pp = 0.0;
146:     for (i=1;i<ctx->np;i++) {
147:       pp = p+x*pp;
148:       p = ctx->pcoeff[i]+x*p;
149:     }
150:   }
151:   if (!ctx->nq) *yp = pp;
152:   else {
153:     q = ctx->qcoeff[0];
154:     qp = 0.0;
155:     for (i=1;i<ctx->nq;i++) {
156:       qp = q+x*qp;
157:       q = ctx->qcoeff[i]+x*q;
158:     }
159:     PetscCheck(q!=0.0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Derivative not defined in the requested value");
160:     *yp = (pp*q-p*qp)/(q*q);
161:   }
162:   PetscFunctionReturn(PETSC_SUCCESS);
163: }

165: static PetscErrorCode FNView_Rational(FN fn,PetscViewer viewer)
166: {
167:   FN_RATIONAL    *ctx = (FN_RATIONAL*)fn->data;
168:   PetscBool      isascii;
169:   PetscInt       i;
170:   char           str[50];

172:   PetscFunctionBegin;
173:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
174:   if (isascii) {
175:     if (fn->alpha!=(PetscScalar)1.0 || fn->beta!=(PetscScalar)1.0) {
176:       PetscCall(SlepcSNPrintfScalar(str,sizeof(str),fn->alpha,PETSC_FALSE));
177:       PetscCall(PetscViewerASCIIPrintf(viewer,"  scale factors: alpha=%s,",str));
178:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
179:       PetscCall(SlepcSNPrintfScalar(str,sizeof(str),fn->beta,PETSC_FALSE));
180:       PetscCall(PetscViewerASCIIPrintf(viewer," beta=%s\n",str));
181:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
182:     }
183:     if (!ctx->nq) {
184:       if (!ctx->np) PetscCall(PetscViewerASCIIPrintf(viewer,"  constant: 1.0\n"));
185:       else if (ctx->np==1) {
186:         PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[0],PETSC_FALSE));
187:         PetscCall(PetscViewerASCIIPrintf(viewer,"  constant: %s\n",str));
188:       } else {
189:         PetscCall(PetscViewerASCIIPrintf(viewer,"  polynomial: "));
190:         PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
191:         for (i=0;i<ctx->np-1;i++) {
192:           PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[i],PETSC_TRUE));
193:           PetscCall(PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->np-i-1));
194:         }
195:         PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[ctx->np-1],PETSC_TRUE));
196:         PetscCall(PetscViewerASCIIPrintf(viewer,"%s\n",str));
197:         PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
198:       }
199:     } else if (!ctx->np) {
200:       PetscCall(PetscViewerASCIIPrintf(viewer,"  inverse polynomial: 1 / ("));
201:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
202:       for (i=0;i<ctx->nq-1;i++) {
203:         PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[i],PETSC_TRUE));
204:         PetscCall(PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->nq-i-1));
205:       }
206:       PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[ctx->nq-1],PETSC_TRUE));
207:       PetscCall(PetscViewerASCIIPrintf(viewer,"%s)\n",str));
208:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
209:     } else {
210:       PetscCall(PetscViewerASCIIPrintf(viewer,"  rational function: ("));
211:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
212:       for (i=0;i<ctx->np-1;i++) {
213:         PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[i],PETSC_TRUE));
214:         PetscCall(PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->np-i-1));
215:       }
216:       PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->pcoeff[ctx->np-1],PETSC_TRUE));
217:       PetscCall(PetscViewerASCIIPrintf(viewer,"%s) / (",str));
218:       for (i=0;i<ctx->nq-1;i++) {
219:         PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[i],PETSC_TRUE));
220:         PetscCall(PetscViewerASCIIPrintf(viewer,"%s*x^%1" PetscInt_FMT,str,ctx->nq-i-1));
221:       }
222:       PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->qcoeff[ctx->nq-1],PETSC_TRUE));
223:       PetscCall(PetscViewerASCIIPrintf(viewer,"%s)\n",str));
224:       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
225:     }
226:   }
227:   PetscFunctionReturn(PETSC_SUCCESS);
228: }

230: static PetscErrorCode FNRationalSetNumerator_Rational(FN fn,PetscInt np,PetscScalar *pcoeff)
231: {
232:   FN_RATIONAL    *ctx = (FN_RATIONAL*)fn->data;
233:   PetscInt       i;

235:   PetscFunctionBegin;
236:   PetscCheck(np>=0,PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Argument np cannot be negative");
237:   ctx->np = np;
238:   PetscCall(PetscFree(ctx->pcoeff));
239:   if (np) {
240:     PetscCall(PetscMalloc1(np,&ctx->pcoeff));
241:     for (i=0;i<np;i++) ctx->pcoeff[i] = pcoeff[i];
242:   }
243:   PetscFunctionReturn(PETSC_SUCCESS);
244: }

246: /*@
247:    FNRationalSetNumerator - Sets the parameters defining the numerator of the
248:    rational function.

250:    Logically Collective

252:    Input Parameters:
253: +  fn     - the math function context
254: .  np     - number of coefficients
255: -  pcoeff - coefficients (array of scalar values)

257:    Notes:
258:    Let the rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials.
259:    This function provides the coefficients of the numerator p(x).
260:    Hence, p(x) is of degree np-1.
261:    If np is zero, then the numerator is assumed to be p(x)=1.

263:    In polynomials, high order coefficients are stored in the first positions
264:    of the array, e.g. to represent x^2-3 use {1,0,-3}.

266:    Level: intermediate

268: .seealso: FNRationalSetDenominator(), FNRationalGetNumerator()
269: @*/
270: PetscErrorCode FNRationalSetNumerator(FN fn,PetscInt np,PetscScalar pcoeff[])
271: {
272:   PetscFunctionBegin;
275:   if (np) PetscAssertPointer(pcoeff,3);
276:   PetscTryMethod(fn,"FNRationalSetNumerator_C",(FN,PetscInt,PetscScalar*),(fn,np,pcoeff));
277:   PetscFunctionReturn(PETSC_SUCCESS);
278: }

280: static PetscErrorCode FNRationalGetNumerator_Rational(FN fn,PetscInt *np,PetscScalar *pcoeff[])
281: {
282:   FN_RATIONAL    *ctx = (FN_RATIONAL*)fn->data;
283:   PetscInt       i;

285:   PetscFunctionBegin;
286:   if (np) *np = ctx->np;
287:   if (pcoeff) {
288:     if (!ctx->np) *pcoeff = NULL;
289:     else {
290:       PetscCall(PetscMalloc1(ctx->np,pcoeff));
291:       for (i=0;i<ctx->np;i++) (*pcoeff)[i] = ctx->pcoeff[i];
292:     }
293:   }
294:   PetscFunctionReturn(PETSC_SUCCESS);
295: }

297: /*@C
298:    FNRationalGetNumerator - Gets the parameters that define the numerator of the
299:    rational function.

301:    Not Collective

303:    Input Parameter:
304: .  fn     - the math function context

306:    Output Parameters:
307: +  np     - number of coefficients
308: -  pcoeff - coefficients (array of scalar values, length nq)

310:    Notes:
311:    The values passed by user with FNRationalSetNumerator() are returned (or null
312:    pointers otherwise).
313:    The pcoeff array should be freed by the user when no longer needed.

315:    Level: intermediate

317: .seealso: FNRationalSetNumerator()
318: @*/
319: PetscErrorCode FNRationalGetNumerator(FN fn,PetscInt *np,PetscScalar *pcoeff[])
320: {
321:   PetscFunctionBegin;
323:   PetscUseMethod(fn,"FNRationalGetNumerator_C",(FN,PetscInt*,PetscScalar**),(fn,np,pcoeff));
324:   PetscFunctionReturn(PETSC_SUCCESS);
325: }

327: static PetscErrorCode FNRationalSetDenominator_Rational(FN fn,PetscInt nq,PetscScalar *qcoeff)
328: {
329:   FN_RATIONAL    *ctx = (FN_RATIONAL*)fn->data;
330:   PetscInt       i;

332:   PetscFunctionBegin;
333:   PetscCheck(nq>=0,PetscObjectComm((PetscObject)fn),PETSC_ERR_ARG_OUTOFRANGE,"Argument nq cannot be negative");
334:   ctx->nq = nq;
335:   PetscCall(PetscFree(ctx->qcoeff));
336:   if (nq) {
337:     PetscCall(PetscMalloc1(nq,&ctx->qcoeff));
338:     for (i=0;i<nq;i++) ctx->qcoeff[i] = qcoeff[i];
339:   }
340:   PetscFunctionReturn(PETSC_SUCCESS);
341: }

343: /*@
344:    FNRationalSetDenominator - Sets the parameters defining the denominator of the
345:    rational function.

347:    Logically Collective

349:    Input Parameters:
350: +  fn     - the math function context
351: .  nq     - number of coefficients
352: -  qcoeff - coefficients (array of scalar values)

354:    Notes:
355:    Let the rational function r(x) = p(x)/q(x), where p(x) and q(x) are polynomials.
356:    This function provides the coefficients of the denominator q(x).
357:    Hence, q(x) is of degree nq-1.
358:    If nq is zero, then the function is assumed to be polynomial, r(x) = p(x).

360:    In polynomials, high order coefficients are stored in the first positions
361:    of the array, e.g. to represent x^2-3 use {1,0,-3}.

363:    Level: intermediate

365: .seealso: FNRationalSetNumerator(), FNRationalGetDenominator()
366: @*/
367: PetscErrorCode FNRationalSetDenominator(FN fn,PetscInt nq,PetscScalar qcoeff[])
368: {
369:   PetscFunctionBegin;
372:   if (nq) PetscAssertPointer(qcoeff,3);
373:   PetscTryMethod(fn,"FNRationalSetDenominator_C",(FN,PetscInt,PetscScalar*),(fn,nq,qcoeff));
374:   PetscFunctionReturn(PETSC_SUCCESS);
375: }

377: static PetscErrorCode FNRationalGetDenominator_Rational(FN fn,PetscInt *nq,PetscScalar *qcoeff[])
378: {
379:   FN_RATIONAL    *ctx = (FN_RATIONAL*)fn->data;
380:   PetscInt       i;

382:   PetscFunctionBegin;
383:   if (nq) *nq = ctx->nq;
384:   if (qcoeff) {
385:     if (!ctx->nq) *qcoeff = NULL;
386:     else {
387:       PetscCall(PetscMalloc1(ctx->nq,qcoeff));
388:       for (i=0;i<ctx->nq;i++) (*qcoeff)[i] = ctx->qcoeff[i];
389:     }
390:   }
391:   PetscFunctionReturn(PETSC_SUCCESS);
392: }

394: /*@C
395:    FNRationalGetDenominator - Gets the parameters that define the denominator of the
396:    rational function.

398:    Not Collective

400:    Input Parameter:
401: .  fn     - the math function context

403:    Output Parameters:
404: +  nq     - number of coefficients
405: -  qcoeff - coefficients (array of scalar values, length nq)

407:    Notes:
408:    The values passed by user with FNRationalSetDenominator() are returned (or a null
409:    pointer otherwise).
410:    The qcoeff array should be freed by the user when no longer needed.

412:    Level: intermediate

414: .seealso: FNRationalSetDenominator()
415: @*/
416: PetscErrorCode FNRationalGetDenominator(FN fn,PetscInt *nq,PetscScalar *qcoeff[])
417: {
418:   PetscFunctionBegin;
420:   PetscUseMethod(fn,"FNRationalGetDenominator_C",(FN,PetscInt*,PetscScalar**),(fn,nq,qcoeff));
421:   PetscFunctionReturn(PETSC_SUCCESS);
422: }

424: static PetscErrorCode FNSetFromOptions_Rational(FN fn,PetscOptionItems *PetscOptionsObject)
425: {
426: #define PARMAX 10
427:   PetscScalar    array[PARMAX];
428:   PetscInt       i,k;
429:   PetscBool      flg;

431:   PetscFunctionBegin;
432:   PetscOptionsHeadBegin(PetscOptionsObject,"FN Rational Options");

434:     k = PARMAX;
435:     for (i=0;i<k;i++) array[i] = 0;
436:     PetscCall(PetscOptionsScalarArray("-fn_rational_numerator","Numerator coefficients (one or more scalar values separated with a comma without spaces)","FNRationalSetNumerator",array,&k,&flg));
437:     if (flg) PetscCall(FNRationalSetNumerator(fn,k,array));

439:     k = PARMAX;
440:     for (i=0;i<k;i++) array[i] = 0;
441:     PetscCall(PetscOptionsScalarArray("-fn_rational_denominator","Denominator coefficients (one or more scalar values separated with a comma without spaces)","FNRationalSetDenominator",array,&k,&flg));
442:     if (flg) PetscCall(FNRationalSetDenominator(fn,k,array));

444:   PetscOptionsHeadEnd();
445:   PetscFunctionReturn(PETSC_SUCCESS);
446: }

448: static PetscErrorCode FNDuplicate_Rational(FN fn,MPI_Comm comm,FN *newfn)
449: {
450:   FN_RATIONAL    *ctx = (FN_RATIONAL*)fn->data,*ctx2 = (FN_RATIONAL*)(*newfn)->data;
451:   PetscInt       i;

453:   PetscFunctionBegin;
454:   ctx2->np = ctx->np;
455:   if (ctx->np) {
456:     PetscCall(PetscMalloc1(ctx->np,&ctx2->pcoeff));
457:     for (i=0;i<ctx->np;i++) ctx2->pcoeff[i] = ctx->pcoeff[i];
458:   }
459:   ctx2->nq = ctx->nq;
460:   if (ctx->nq) {
461:     PetscCall(PetscMalloc1(ctx->nq,&ctx2->qcoeff));
462:     for (i=0;i<ctx->nq;i++) ctx2->qcoeff[i] = ctx->qcoeff[i];
463:   }
464:   PetscFunctionReturn(PETSC_SUCCESS);
465: }

467: static PetscErrorCode FNDestroy_Rational(FN fn)
468: {
469:   FN_RATIONAL    *ctx = (FN_RATIONAL*)fn->data;

471:   PetscFunctionBegin;
472:   PetscCall(PetscFree(ctx->pcoeff));
473:   PetscCall(PetscFree(ctx->qcoeff));
474:   PetscCall(PetscFree(fn->data));
475:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",NULL));
476:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",NULL));
477:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",NULL));
478:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",NULL));
479:   PetscFunctionReturn(PETSC_SUCCESS);
480: }

482: SLEPC_EXTERN PetscErrorCode FNCreate_Rational(FN fn)
483: {
484:   FN_RATIONAL    *ctx;

486:   PetscFunctionBegin;
487:   PetscCall(PetscNew(&ctx));
488:   fn->data = (void*)ctx;

490:   fn->ops->evaluatefunction          = FNEvaluateFunction_Rational;
491:   fn->ops->evaluatederivative        = FNEvaluateDerivative_Rational;
492:   fn->ops->evaluatefunctionmat[0]    = FNEvaluateFunctionMat_Rational;
493:   fn->ops->evaluatefunctionmatvec[0] = FNEvaluateFunctionMatVec_Rational;
494: #if defined(PETSC_HAVE_CUDA)
495:   fn->ops->evaluatefunctionmatcuda[0]    = FNEvaluateFunctionMat_Rational;
496:   fn->ops->evaluatefunctionmatveccuda[0] = FNEvaluateFunctionMatVec_Rational;
497: #endif
498:   fn->ops->setfromoptions            = FNSetFromOptions_Rational;
499:   fn->ops->view                      = FNView_Rational;
500:   fn->ops->duplicate                 = FNDuplicate_Rational;
501:   fn->ops->destroy                   = FNDestroy_Rational;
502:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetNumerator_C",FNRationalSetNumerator_Rational));
503:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetNumerator_C",FNRationalGetNumerator_Rational));
504:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalSetDenominator_C",FNRationalSetDenominator_Rational));
505:   PetscCall(PetscObjectComposeFunction((PetscObject)fn,"FNRationalGetDenominator_C",FNRationalGetDenominator_Rational));
506:   PetscFunctionReturn(PETSC_SUCCESS);
507: }