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