Actual source code: test1.c
slepc-3.21.2 2024-09-25
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: */
11: static char help[] = "Test rational function.\n\n";
13: #include <slepcfn.h>
15: int main(int argc,char **argv)
16: {
17: FN fn;
18: PetscInt i,na,nb;
19: PetscScalar x,y,yp,p[10],q[10],five=5.0,*pp,*qq;
20: char strx[50],str[50];
22: PetscFunctionBeginUser;
23: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
24: PetscCall(FNCreate(PETSC_COMM_WORLD,&fn));
26: /* polynomial p(x) */
27: na = 5;
28: p[0] = -3.1; p[1] = 1.1; p[2] = 1.0; p[3] = -2.0; p[4] = 3.5;
29: PetscCall(FNSetType(fn,FNRATIONAL));
30: PetscCall(FNRationalSetNumerator(fn,na,p));
31: PetscCall(FNView(fn,NULL));
32: x = 2.2;
33: PetscCall(SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE));
34: PetscCall(FNEvaluateFunction(fn,x,&y));
35: PetscCall(FNEvaluateDerivative(fn,x,&yp));
36: PetscCall(SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE));
37: PetscCall(PetscPrintf(PETSC_COMM_WORLD," f(%s)=%s\n",strx,str));
38: PetscCall(SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE));
39: PetscCall(PetscPrintf(PETSC_COMM_WORLD," f'(%s)=%s\n",strx,str));
41: /* inverse of polynomial 1/q(x) */
42: nb = 3;
43: q[0] = -3.1; q[1] = 1.1; q[2] = 1.0;
44: PetscCall(FNSetType(fn,FNRATIONAL));
45: PetscCall(FNRationalSetNumerator(fn,0,NULL)); /* reset previous values */
46: PetscCall(FNRationalSetDenominator(fn,nb,q));
47: PetscCall(FNView(fn,NULL));
48: x = 2.2;
49: PetscCall(SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE));
50: PetscCall(FNEvaluateFunction(fn,x,&y));
51: PetscCall(FNEvaluateDerivative(fn,x,&yp));
52: PetscCall(SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE));
53: PetscCall(PetscPrintf(PETSC_COMM_WORLD," f(%s)=%s\n",strx,str));
54: PetscCall(SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE));
55: PetscCall(PetscPrintf(PETSC_COMM_WORLD," f'(%s)=%s\n",strx,str));
57: /* rational p(x)/q(x) */
58: na = 2; nb = 3;
59: p[0] = 1.1; p[1] = 1.1;
60: q[0] = 1.0; q[1] = -2.0; q[2] = 3.5;
61: PetscCall(FNSetType(fn,FNRATIONAL));
62: PetscCall(FNRationalSetNumerator(fn,na,p));
63: PetscCall(FNRationalSetDenominator(fn,nb,q));
64: PetscCall(FNSetScale(fn,1.2,0.5));
65: PetscCall(FNView(fn,NULL));
66: x = 2.2;
67: PetscCall(SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE));
68: PetscCall(FNEvaluateFunction(fn,x,&y));
69: PetscCall(FNEvaluateDerivative(fn,x,&yp));
70: PetscCall(SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE));
71: PetscCall(PetscPrintf(PETSC_COMM_WORLD," f(%s)=%s\n",strx,str));
72: PetscCall(SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE));
73: PetscCall(PetscPrintf(PETSC_COMM_WORLD," f'(%s)=%s\n",strx,str));
75: PetscCall(FNRationalGetNumerator(fn,&na,&pp));
76: PetscCall(FNRationalGetDenominator(fn,&nb,&qq));
77: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Coefficients:\n Numerator: "));
78: for (i=0;i<na;i++) {
79: PetscCall(SlepcSNPrintfScalar(str,sizeof(str),pp[i],PETSC_FALSE));
80: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s ",str));
81: }
82: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n Denominator: "));
83: for (i=0;i<nb;i++) {
84: PetscCall(SlepcSNPrintfScalar(str,sizeof(str),qq[i],PETSC_FALSE));
85: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s ",str));
86: }
87: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
88: PetscCall(PetscFree(pp));
89: PetscCall(PetscFree(qq));
91: /* constant */
92: PetscCall(FNSetType(fn,FNRATIONAL));
93: PetscCall(FNRationalSetNumerator(fn,1,&five));
94: PetscCall(FNRationalSetDenominator(fn,0,NULL)); /* reset previous values */
95: PetscCall(FNView(fn,NULL));
96: x = 2.2;
97: PetscCall(SlepcSNPrintfScalar(strx,sizeof(strx),x,PETSC_FALSE));
98: PetscCall(FNEvaluateFunction(fn,x,&y));
99: PetscCall(FNEvaluateDerivative(fn,x,&yp));
100: PetscCall(SlepcSNPrintfScalar(str,sizeof(str),y,PETSC_FALSE));
101: PetscCall(PetscPrintf(PETSC_COMM_WORLD," f(%s)=%s\n",strx,str));
102: PetscCall(SlepcSNPrintfScalar(str,sizeof(str),yp,PETSC_FALSE));
103: PetscCall(PetscPrintf(PETSC_COMM_WORLD," f'(%s)=%s\n",strx,str));
105: PetscCall(FNDestroy(&fn));
106: PetscCall(SlepcFinalize());
107: return 0;
108: }
110: /*TEST
112: test:
113: suffix: 1
114: nsize: 1
116: TEST*/