xref: /petsc/include/petscmath.h (revision ff4e02e2c8501fafdc1e44c945f7f2b3d8ea9565)
1 /* $Id: petscmath.h,v 1.32 2001/08/30 20:37:06 bsmith Exp $ */
2 /*
3 
4       PETSc mathematics include file. Defines certain basic mathematical
5     constants and functions for working with single and double precision
6     floating point numbers as well as complex and integers.
7 
8     This file is included by petsc.h and should not be used directly.
9 
10 */
11 
12 #if !defined(__PETSCMATH_H)
13 #define __PETSCMATH_H
14 #include <math.h>
15 
16 /*
17 
18      Defines operations that are different for complex and real numbers;
19    note that one cannot really mix the use of complex and real in the same
20    PETSc program. All PETSc objects in one program are built around the object
21    PetscScalar which is either always a double or a complex.
22 
23 */
24 #if defined(PETSC_USE_COMPLEX)
25 
26 /*
27    PETSc now only supports std::complex
28 */
29 #include <complex>
30 
31 extern  MPI_Datatype        MPIU_COMPLEX;
32 #define MPIU_SCALAR         MPIU_COMPLEX
33 #if defined(PETSC_USE_MAT_SINGLE)
34 #define MPIU_MATSCALAR        ??Notdone
35 #else
36 #define MPIU_MATSCALAR      MPIU_COMPLEX
37 #endif
38 
39 #define PetscRealPart(a)        (a).real()
40 #define PetscImaginaryPart(a)   (a).imag()
41 #define PetscAbsScalar(a)   std::abs(a)
42 #define PetscConj(a)        std::conj(a)
43 #define PetscSqrtScalar(a)  std::sqrt(a)
44 #define PetscPowScalar(a,b) std::pow(a,b)
45 #define PetscExpScalar(a)   std::exp(a)
46 #define PetscSinScalar(a)   std::sin(a)
47 #define PetscCosScalar(a)   std::cos(a)
48 
49 typedef std::complex<double> PetscScalar;
50 
51 /* Compiling for real numbers only */
52 #else
53 #  if defined(PETSC_USE_SINGLE)
54 #    define MPIU_SCALAR           MPI_FLOAT
55 #  else
56 #    define MPIU_SCALAR           MPI_DOUBLE
57 #  endif
58 #  if defined(PETSC_USE_MAT_SINGLE) || defined(PETSC_USE_SINGLE)
59 #    define MPIU_MATSCALAR        MPI_FLOAT
60 #  else
61 #    define MPIU_MATSCALAR        MPI_DOUBLE
62 #  endif
63 #  define PetscRealPart(a)      (a)
64 #  define PetscImaginaryPart(a) (a)
65 #  define PetscAbsScalar(a)     (((a)<0.0)   ? -(a) : (a))
66 #  define PetscConj(a)          (a)
67 #  define PetscSqrtScalar(a)    sqrt(a)
68 #  define PetscPowScalar(a,b)   pow(a,b)
69 #  define PetscExpScalar(a)     exp(a)
70 #  define PetscSinScalar(a)     sin(a)
71 #  define PetscCosScalar(a)     cos(a)
72 
73 #  if defined(PETSC_USE_SINGLE)
74   typedef float PetscScalar;
75 #  else
76   typedef double PetscScalar;
77 #  endif
78 #endif
79 
80 #if defined(PETSC_USE_SINGLE)
81 #  define MPIU_REAL   MPI_FLOAT
82 #else
83 #  define MPIU_REAL   MPI_DOUBLE
84 #endif
85 
86 #define PetscSign(a) (((a) >= 0) ? ((a) == 0 ? 0 : 1) : -1)
87 #define PetscAbs(a)  (((a) >= 0) ? a : -a)
88 /*
89        Allows compiling PETSc so that matrix values are stored in
90    single precision but all other objects still use double
91    precision. This does not work for complex numbers in that case
92    it remains double
93 
94           EXPERIMENTAL! NOT YET COMPLETELY WORKING
95 */
96 
97 #if defined(PETSC_USE_MAT_SINGLE)
98 typedef float MatScalar;
99 #else
100 typedef PetscScalar MatScalar;
101 #endif
102 
103 #if defined(PETSC_USE_COMPLEX)
104 typedef double MatReal;
105 #elif defined(PETSC_USE_MAT_SINGLE) || defined(PETSC_USE_SINGLE)
106 typedef float MatReal;
107 #else
108 typedef double MatReal;
109 #endif
110 
111 #if defined(PETSC_USE_SINGLE)
112   typedef float PetscReal;
113 #else
114   typedef double PetscReal;
115 #endif
116 
117 /* --------------------------------------------------------------------------*/
118 
119 /*
120    Certain objects may be created using either single
121   or double precision.
122 */
123 typedef enum { PETSC_SCALAR_DOUBLE,PETSC_SCALAR_SINGLE } PetscScalarPrecision;
124 
125 /* PETSC_i is the imaginary number, i */
126 extern  PetscScalar       PETSC_i;
127 
128 #define PetscMin(a,b)   (((a)<(b)) ?  (a) : (b))
129 #define PetscMax(a,b)   (((a)<(b)) ?  (b) : (a))
130 #define PetscAbsInt(a)  (((a)<0)   ? -(a) : (a))
131 #define PetscAbsReal(a) (((a)<0)   ? -(a) : (a))
132 #define PetscSqr(a)     ((a)*(a))
133 
134 /* ----------------------------------------------------------------------------*/
135 /*
136      Basic constants
137 */
138 #define PETSC_PI                 3.14159265358979323846264
139 #define PETSC_DEGREES_TO_RADIANS 0.01745329251994
140 #define PETSC_MAX_INT            1000000000
141 #define PETSC_MIN_INT            -1000000000
142 
143 #if defined(PETSC_USE_SINGLE)
144 #  define PETSC_MAX                1.e30
145 #  define PETSC_MIN                -1.e30
146 #  define PETSC_MACHINE_EPSILON         1.e-7
147 #  define PETSC_SQRT_MACHINE_EPSILON    3.e-4
148 #  define PETSC_SMALL                   1.e-5
149 #else
150 #  define PETSC_MAX                1.e300
151 #  define PETSC_MIN                -1.e300
152 #  define PETSC_MACHINE_EPSILON         1.e-14
153 #  define PETSC_SQRT_MACHINE_EPSILON    1.e-7
154 #  define PETSC_SMALL                   1.e-10
155 #endif
156 
157 extern int PetscGlobalMax(PetscReal*,PetscReal*,MPI_Comm);
158 extern int PetscGlobalMin(PetscReal*,PetscReal*,MPI_Comm);
159 extern int PetscGlobalSum(PetscScalar*,PetscScalar*,MPI_Comm);
160 
161 
162 /* ----------------------------------------------------------------------------*/
163 /*
164     PetscLogDouble variables are used to contain double precision numbers
165   that are not used in the numerical computations, but rather in logging,
166   timing etc.
167 */
168 typedef double PetscLogDouble;
169 /*
170       Once PETSc is compiling with a ADIC enhanced version of MPI
171    we will create a new MPI_Datatype for the inactive double variables.
172 */
173 #if defined(AD_DERIV_H)
174 /* extern  MPI_Datatype  MPIU_PETSCLOGDOUBLE; */
175 #else
176 #if !defined(_petsc_mpi_uni)
177 #define MPIU_PETSCLOGDOUBLE MPI_DOUBLE
178 #endif
179 #endif
180 
181 #define PassiveReal PetscReal
182 #define PassiveScalar PetscScalar
183 
184 #define PETSCMAP1_a(a,b)  a ## _ ## b
185 #define PETSCMAP1_b(a,b)  PETSCMAP1_a(a,b)
186 #define PETSCMAP1(a)      PETSCMAP1_b(a,PetscScalar)
187 #endif
188