xref: /petsc/src/sys/tests/ex25.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests wrapping of math.h functions for real, complex, and scalar types \n";
2*c4762a1bSJed Brown #include <petscsys.h>
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown int main(int argc,char **argv)
5*c4762a1bSJed Brown {
6*c4762a1bSJed Brown   PetscErrorCode ierr;
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
9*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Real tests:\n");CHKERRQ(ierr);
10*c4762a1bSJed Brown   {
11*c4762a1bSJed Brown     PetscReal a,b,c;
12*c4762a1bSJed Brown     a = PetscRealConstant(0.5);
13*c4762a1bSJed Brown     c = PetscRealConstant(2.0);
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown     b = PetscSqrtReal(a);
16*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"sqrt(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
17*c4762a1bSJed Brown     b = PetscCbrtReal(a);
18*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"cbrt(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown     b = PetscHypotReal(a,c);
21*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"hypot(%f,%f) = %f\n",(double)a,(double)c,(double)b);CHKERRQ(ierr);
22*c4762a1bSJed Brown     b = PetscAtan2Real(a,c);
23*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"atan2(%f,%f) = %f\n",(double)a,(double)c,(double)b);CHKERRQ(ierr);
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown     b = PetscPowReal(a,c);
26*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"pow(%f,%f) = %f\n",(double)a,(double)c,(double)b);CHKERRQ(ierr);
27*c4762a1bSJed Brown     b = PetscExpReal(a);
28*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"exp(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
29*c4762a1bSJed Brown     b = PetscLogReal(a);
30*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"log(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
31*c4762a1bSJed Brown     b = PetscLog10Real(a);
32*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"log10(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
33*c4762a1bSJed Brown     b = PetscLog2Real(a);
34*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"log2(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown     b = PetscSinReal(a);
37*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"sin(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
38*c4762a1bSJed Brown     b = PetscCosReal(a);
39*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"cos(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
40*c4762a1bSJed Brown     b = PetscTanReal(a);
41*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"tan(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown     b = PetscAsinReal(a);
44*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"asin(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
45*c4762a1bSJed Brown     b = PetscAcosReal(a);
46*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"acos(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
47*c4762a1bSJed Brown     b = PetscAtanReal(a);
48*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"atan(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown     b = PetscSinhReal(a);
51*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"sinh(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
52*c4762a1bSJed Brown     b = PetscCoshReal(a);
53*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"cosh(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
54*c4762a1bSJed Brown     b = PetscTanhReal(a);
55*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"tanh(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
56*c4762a1bSJed Brown 
57*c4762a1bSJed Brown     b = PetscAsinhReal(a);
58*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"asinh(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
59*c4762a1bSJed Brown     b = PetscAcoshReal(c);
60*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"acosh(%f) = %f\n",(double)c,(double)b);CHKERRQ(ierr);
61*c4762a1bSJed Brown     b = PetscAtanhReal(a);
62*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"atanh(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown     b = PetscCeilReal(a);
65*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"ceil(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
66*c4762a1bSJed Brown     b = PetscFloorReal(a);
67*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"floor(%f) = %f\n",(double)a,(double)b);CHKERRQ(ierr);
68*c4762a1bSJed Brown     b = PetscFmodReal(a,c);
69*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"fmod(%f,%f) = %f\n",(double)a,(double)c,(double)b);CHKERRQ(ierr);
70*c4762a1bSJed Brown   }
71*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Scalar tests:\n");CHKERRQ(ierr);
72*c4762a1bSJed Brown   {
73*c4762a1bSJed Brown     PetscScalar a,b,c;
74*c4762a1bSJed Brown     a = PetscRealConstant(0.5);
75*c4762a1bSJed Brown     c = PetscRealConstant(2.0);
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown     b = PetscAbsScalar(a);
78*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"abs(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
79*c4762a1bSJed Brown     b = PetscArgScalar(a);
80*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"arg(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
81*c4762a1bSJed Brown     b = PetscConj(a);
82*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"conj(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
83*c4762a1bSJed Brown 
84*c4762a1bSJed Brown     b = PetscSqrtScalar(a);
85*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"sqrt(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown     b = PetscPowScalar(a,c);
88*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"pow(%f,%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(c),(double)PetscRealPart(b));CHKERRQ(ierr);
89*c4762a1bSJed Brown     b = PetscExpScalar(a);
90*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"exp(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
91*c4762a1bSJed Brown     b = PetscLogScalar(a);
92*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"log(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
93*c4762a1bSJed Brown 
94*c4762a1bSJed Brown     b = PetscSinScalar(a);
95*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"sin(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
96*c4762a1bSJed Brown     b = PetscCosScalar(a);
97*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"cos(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
98*c4762a1bSJed Brown     b = PetscTanScalar(a);
99*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"tan(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown     b = PetscAsinScalar(a);
102*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"asin(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
103*c4762a1bSJed Brown     b = PetscAcosScalar(a);
104*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"acos(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
105*c4762a1bSJed Brown     b = PetscAtanScalar(a);
106*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"atan(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
107*c4762a1bSJed Brown 
108*c4762a1bSJed Brown     b = PetscSinhScalar(a);
109*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"sinh(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
110*c4762a1bSJed Brown     b = PetscCoshScalar(a);
111*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"cosh(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
112*c4762a1bSJed Brown     b = PetscTanhScalar(a);
113*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"tanh(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
114*c4762a1bSJed Brown 
115*c4762a1bSJed Brown     b = PetscAsinhScalar(a);
116*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"asinh(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
117*c4762a1bSJed Brown     b = PetscAcoshScalar(c);
118*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"acosh(%f) = %f\n",(double)PetscRealPart(c),(double)PetscRealPart(b));CHKERRQ(ierr);
119*c4762a1bSJed Brown     b = PetscAtanhScalar(a);
120*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"atanh(%f) = %f\n",(double)PetscRealPart(a),(double)PetscRealPart(b));CHKERRQ(ierr);
121*c4762a1bSJed Brown   }
122*c4762a1bSJed Brown   ierr = PetscFinalize();
123*c4762a1bSJed Brown   return ierr;
124*c4762a1bSJed Brown }
125*c4762a1bSJed Brown 
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown /*TEST
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown    test:
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown TEST*/
132