xref: /petsc/include/petscmath.h (revision feff33ee0b5b037fa8f9f294dede656a2f85cc47)
1 /*
2 
3       PETSc mathematics include file. Defines certain basic mathematical
4     constants and functions for working with single, double, and quad precision
5     floating point numbers as well as complex single and double.
6 
7     This file is included by petscsys.h and should not be used directly.
8 
9 */
10 
11 #if !defined(__PETSCMATH_H)
12 #define __PETSCMATH_H
13 #include <math.h>
14 
15 /*
16 
17      Defines operations that are different for complex and real numbers;
18    note that one cannot mix the use of complex and real in the same
19    PETSc program. All PETSc objects in one program are built around the object
20    PetscScalar which is either always a real or a complex.
21 
22 */
23 
24 #if defined(PETSC_USE_REAL_SINGLE)
25 #define MPIU_REAL   MPI_FLOAT
26 typedef float PetscReal;
27 #define PetscRoundReal(a)   roundf(a)
28 #define PetscSqrtReal(a)    sqrtf(a)
29 #define PetscExpReal(a)     expf(a)
30 #define PetscLogReal(a)     logf(a)
31 #define PetscLog10Real(a)   log10f(a)
32 #ifdef PETSC_HAVE_LOG2
33 #define PetscLog2Real(a)    log2f(a)
34 #endif
35 #define PetscSinReal(a)     sinf(a)
36 #define PetscCosReal(a)     cosf(a)
37 #define PetscTanReal(a)     tanf(a)
38 #define PetscAsinReal(a)    asinf(a)
39 #define PetscAcosReal(a)    acosf(a)
40 #define PetscAtanReal(a)    atanf(a)
41 #define PetscAtan2Real(a,b) atan2f(a,b)
42 #define PetscSinhReal(a)    sinhf(a)
43 #define PetscCoshReal(a)    coshf(a)
44 #define PetscTanhReal(a)    tanhf(a)
45 #define PetscPowReal(a,b)   powf(a,b)
46 #define PetscCeilReal(a)    ceilf(a)
47 #define PetscFloorReal(a)   floorf(a)
48 #define PetscFmodReal(a,b)  fmodf(a,b)
49 #define PetscTGamma(a)      tgammaf(a)
50 #elif defined(PETSC_USE_REAL_DOUBLE)
51 #define MPIU_REAL   MPI_DOUBLE
52 typedef double PetscReal;
53 #define PetscRoundReal(a)   round(a)
54 #define PetscSqrtReal(a)    sqrt(a)
55 #define PetscExpReal(a)     exp(a)
56 #define PetscLogReal(a)     log(a)
57 #define PetscLog10Real(a)   log10(a)
58 #ifdef PETSC_HAVE_LOG2
59 #define PetscLog2Real(a)    log2(a)
60 #endif
61 #define PetscSinReal(a)     sin(a)
62 #define PetscCosReal(a)     cos(a)
63 #define PetscTanReal(a)     tan(a)
64 #define PetscAsinReal(a)    asin(a)
65 #define PetscAcosReal(a)    acos(a)
66 #define PetscAtanReal(a)    atan(a)
67 #define PetscAtan2Real(a,b) atan2(a,b)
68 #define PetscSinhReal(a)    sinh(a)
69 #define PetscCoshReal(a)    cosh(a)
70 #define PetscTanhReal(a)    tanh(a)
71 #define PetscPowReal(a,b)   pow(a,b)
72 #define PetscCeilReal(a)    ceil(a)
73 #define PetscFloorReal(a)   floor(a)
74 #define PetscFmodReal(a,b)  fmod(a,b)
75 #define PetscTGamma(a)      tgamma(a)
76 #elif defined(PETSC_USE_REAL___FLOAT128)
77 #if defined(__cplusplus)
78 extern "C" {
79 #endif
80 #include <quadmath.h>
81 #if defined(__cplusplus)
82 }
83 #endif
84 PETSC_EXTERN MPI_Datatype MPIU___FLOAT128 PetscAttrMPITypeTag(__float128);
85 #define MPIU_REAL MPIU___FLOAT128
86 typedef __float128 PetscReal;
87 #define PetscRoundReal(a)   roundq(a)
88 #define PetscSqrtReal(a)    sqrtq(a)
89 #define PetscExpReal(a)     expq(a)
90 #define PetscLogReal(a)     logq(a)
91 #define PetscLog10Real(a)   log10q(a)
92 #ifdef PETSC_HAVE_LOG2
93 #define PetscLog2Real(a)    log2q(a)
94 #endif
95 #define PetscSinReal(a)     sinq(a)
96 #define PetscCosReal(a)     cosq(a)
97 #define PetscTanReal(a)     tanq(a)
98 #define PetscAsinReal(a)    asinq(a)
99 #define PetscAcosReal(a)    acosq(a)
100 #define PetscAtanReal(a)    atanq(a)
101 #define PetscAtan2Real(a,b) atan2q(a,b)
102 #define PetscSinhReal(a)    sinhq(a)
103 #define PetscCoshReal(a)    coshq(a)
104 #define PetscTanhReal(a)    tanhq(a)
105 #define PetscPowReal(a,b)   powq(a,b)
106 #define PetscCeilReal(a)    ceilq(a)
107 #define PetscFloorReal(a)   floorq(a)
108 #define PetscFmodReal(a,b)  fmodq(a,b)
109 #define PetscTGamma(a)      tgammaq(a)
110 #elif defined(PETSC_USE_REAL___FP16)
111 PETSC_EXTERN MPI_Datatype MPIU___FP16 PetscAttrMPITypeTag(__fp16);
112 #define MPIU_REAL MPIU___FP16
113 typedef __fp16 PetscReal;
114 #define PetscRound(a)       roundf(a)
115 #define PetscSqrtReal(a)    sqrtf(a)
116 #define PetscExpReal(a)     expf(a)
117 #define PetscLogReal(a)     logf(a)
118 #define PetscLog10Real(a)   log10f(a)
119 #ifdef PETSC_HAVE_LOG2
120 #define PetscLog2Real(a)    log2f(a)
121 #endif
122 #define PetscSinReal(a)     sinf(a)
123 #define PetscCosReal(a)     cosf(a)
124 #define PetscTanReal(a)     tanf(a)
125 #define PetscAsinReal(a)    asinf(a)
126 #define PetscAcosReal(a)    acosf(a)
127 #define PetscAtanReal(a)    atanf(a)
128 #define PetscAtan2Real(a,b) atan2f(a,b)
129 #define PetscSinhReal(a)    sinhf(a)
130 #define PetscCoshReal(a)    coshf(a)
131 #define PetscTanhReal(a)    tanhf(a)
132 #define PetscPowReal(a,b)   powf(a,b)
133 #define PetscCeilReal(a)    ceilf(a)
134 #define PetscFloorReal(a)   floorf(a)
135 #define PetscFmodReal(a,b)  fmodf(a,b)
136 #define PetscTGamma(a)      tgammaf(a)
137 #endif /* PETSC_USE_REAL_* */
138 
139 /*
140     Complex number definitions
141  */
142 #if defined(__cplusplus) && defined(PETSC_HAVE_CXX_COMPLEX) && !defined(PETSC_USE_REAL___FLOAT128)
143 #if !defined(PETSC_SKIP_COMPLEX)
144 #define PETSC_HAVE_COMPLEX 1
145 /* C++ support of complex number */
146 #if defined(PETSC_HAVE_VECCUDA) && __CUDACC_VER_MAJOR__ > 6
147 /* complex headers in thrust only available in CUDA 7.0 and above */
148 #define complexlib thrust
149 #include <thrust/complex.h>
150 #else
151 #define complexlib std
152 #include <complex>
153 #endif
154 
155 #define PetscRealPartComplex(a)      (a).real()
156 #define PetscImaginaryPartComplex(a) (a).imag()
157 #define PetscAbsComplex(a)           complexlib::abs(a)
158 #define PetscConjComplex(a)          complexlib::conj(a)
159 #define PetscSqrtComplex(a)          complexlib::sqrt(a)
160 #define PetscPowComplex(a,b)         complexlib::pow(a,b)
161 #define PetscExpComplex(a)           complexlib::exp(a)
162 #define PetscLogComplex(a)           complexlib::log(a)
163 #define PetscSinComplex(a)           complexlib::sin(a)
164 #define PetscCosComplex(a)           complexlib::cos(a)
165 #define PetscAsinComplex(a)          complexlib::asin(a)
166 #define PetscAcosComplex(a)          complexlib::acos(a)
167 #if defined(PETSC_HAVE_TANCOMPLEX)
168 #define PetscTanComplex(a)           complexlib::tan(a)
169 #else
170 #define PetscTanComplex(a)           PetscSinComplex(a)/PetscCosComplex(a)
171 #endif
172 #define PetscSinhComplex(a)          complexlib::sinh(a)
173 #define PetscCoshComplex(a)          complexlib::cosh(a)
174 #if defined(PETSC_HAVE_TANHCOMPLEX)
175 #define PetscTanhComplex(a)          complexlib::tanh(a)
176 #else
177 #define PetscTanhComplex(a)          PetscSinhComplex(a)/PetscCoshComplex(a)
178 #endif
179 
180 #if defined(PETSC_USE_REAL_SINGLE)
181 typedef complexlib::complex<float> PetscComplex;
182 #if defined(PETSC_USE_CXX_COMPLEX_FLOAT_WORKAROUND)
183 static inline PetscComplex operator+(const PetscComplex& lhs, const double& rhs) { return lhs + float(rhs); }
184 static inline PetscComplex operator+(const double& lhs, const PetscComplex& rhs) { return float(lhs) + rhs; }
185 static inline PetscComplex operator-(const PetscComplex& lhs, const double& rhs) { return lhs - float(rhs); }
186 static inline PetscComplex operator-(const double& lhs, const PetscComplex& rhs) { return float(lhs) - rhs; }
187 static inline PetscComplex operator*(const PetscComplex& lhs, const double& rhs) { return lhs * float(rhs); }
188 static inline PetscComplex operator*(const double& lhs, const PetscComplex& rhs) { return float(lhs) * rhs; }
189 static inline PetscComplex operator/(const PetscComplex& lhs, const double& rhs) { return lhs / float(rhs); }
190 static inline PetscComplex operator/(const double& lhs, const PetscComplex& rhs) { return float(lhs) / rhs; }
191 static inline bool operator==(const PetscComplex& lhs, const double& rhs) { return lhs.imag() == float(0) && lhs.real() == float(rhs); }
192 static inline bool operator==(const double& lhs, const PetscComplex& rhs) { return rhs.imag() == float(0) && rhs.real() == float(lhs); }
193 static inline bool operator!=(const PetscComplex& lhs, const double& rhs) { return lhs.imag() != float(0) || lhs.real() != float(rhs); }
194 static inline bool operator!=(const double& lhs, const PetscComplex& rhs) { return rhs.imag() != float(0) || rhs.real() != float(lhs); }
195 #endif  /* PETSC_USE_CXX_COMPLEX_FLOAT_WORKAROUND */
196 #elif defined(PETSC_USE_REAL_DOUBLE)
197 typedef complexlib::complex<double> PetscComplex;
198 #if defined(PETSC_USE_CXX_COMPLEX_FLOAT_WORKAROUND)
199 static inline PetscComplex operator+(const PetscComplex& lhs, const PetscInt& rhs) { return lhs + double(rhs); }
200 static inline PetscComplex operator+(const PetscInt& lhs, const PetscComplex& rhs) { return double(lhs) + rhs; }
201 static inline PetscComplex operator-(const PetscComplex& lhs, const PetscInt& rhs) { return lhs - double(rhs); }
202 static inline PetscComplex operator-(const PetscInt& lhs, const PetscComplex& rhs) { return double(lhs) - rhs; }
203 static inline PetscComplex operator*(const PetscComplex& lhs, const PetscInt& rhs) { return lhs * double(rhs); }
204 static inline PetscComplex operator*(const PetscInt& lhs, const PetscComplex& rhs) { return double(lhs) * rhs; }
205 static inline PetscComplex operator/(const PetscComplex& lhs, const PetscInt& rhs) { return lhs / double(rhs); }
206 static inline PetscComplex operator/(const PetscInt& lhs, const PetscComplex& rhs) { return double(lhs) / rhs; }
207 static inline bool operator==(const PetscComplex& lhs, const PetscInt& rhs) { return lhs.imag() == double(0) && lhs.real() == double(rhs); }
208 static inline bool operator==(const PetscInt& lhs, const PetscComplex& rhs) { return rhs.imag() == double(0) && rhs.real() == double(lhs); }
209 static inline bool operator!=(const PetscComplex& lhs, const PetscInt& rhs) { return lhs.imag() != double(0) || lhs.real() != double(rhs); }
210 static inline bool operator!=(const PetscInt& lhs, const PetscComplex& rhs) { return rhs.imag() != double(0) || rhs.real() != double(lhs); }
211 #endif  /* PETSC_USE_CXX_COMPLEX_FLOAT_WORKAROUND */
212 #elif defined(PETSC_USE_REAL___FLOAT128)
213 typedef complexlib::complex<__float128> PetscComplex; /* Notstandard and not expected to work, use __complex128 */
214 PETSC_EXTERN MPI_Datatype MPIU___COMPLEX128;
215 #endif  /* PETSC_USE_REAL_ */
216 #endif  /* ! PETSC_SKIP_COMPLEX */
217 
218 #elif defined(PETSC_HAVE_C99_COMPLEX) && !defined(PETSC_USE_REAL___FP16)
219 #if !defined(PETSC_SKIP_COMPLEX)
220 #define PETSC_HAVE_COMPLEX 1
221 #include <complex.h>
222 
223 #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
224 typedef float _Complex PetscComplex;
225 
226 #define PetscRealPartComplex(a)      crealf(a)
227 #define PetscImaginaryPartComplex(a) cimagf(a)
228 #define PetscAbsComplex(a)           cabsf(a)
229 #define PetscConjComplex(a)          conjf(a)
230 #define PetscSqrtComplex(a)          csqrtf(a)
231 #define PetscPowComplex(a,b)         cpowf(a,b)
232 #define PetscExpComplex(a)           cexpf(a)
233 #define PetscLogComplex(a)           clogf(a)
234 #define PetscSinComplex(a)           csinf(a)
235 #define PetscCosComplex(a)           ccosf(a)
236 #define PetscAsinComplex(a)          casinf(a)
237 #define PetscAcosComplex(a)          cacosf(a)
238 #define PetscTanComplex(a)           ctanf(a)
239 #define PetscSinhComplex(a)          csinhf(a)
240 #define PetscCoshComplex(a)          ccoshf(a)
241 #define PetscTanhComplex(a)          ctanhf(a)
242 
243 #elif defined(PETSC_USE_REAL_DOUBLE)
244 typedef double _Complex PetscComplex;
245 
246 #define PetscRealPartComplex(a)      creal(a)
247 #define PetscImaginaryPartComplex(a) cimag(a)
248 #define PetscAbsComplex(a)           cabs(a)
249 #define PetscConjComplex(a)          conj(a)
250 #define PetscSqrtComplex(a)          csqrt(a)
251 #define PetscPowComplex(a,b)         cpow(a,b)
252 #define PetscExpComplex(a)           cexp(a)
253 #define PetscLogComplex(a)           clog(a)
254 #define PetscSinComplex(a)           csin(a)
255 #define PetscCosComplex(a)           ccos(a)
256 #define PetscAsinComplex(a)          casin(a)
257 #define PetscAcosComplex(a)          cacos(a)
258 #define PetscTanComplex(a)           ctan(a)
259 #define PetscSinhComplex(a)          csinh(a)
260 #define PetscCoshComplex(a)          ccosh(a)
261 #define PetscTanhComplex(a)          ctanh(a)
262 
263 #elif defined(PETSC_USE_REAL___FLOAT128)
264 typedef __complex128 PetscComplex;
265 PETSC_EXTERN MPI_Datatype MPIU___COMPLEX128 PetscAttrMPITypeTag(__complex128);
266 
267 #define PetscRealPartComplex(a)      crealq(a)
268 #define PetscImaginaryPartComplex(a) cimagq(a)
269 #define PetscAbsComplex(a)           cabsq(a)
270 #define PetscConjComplex(a)          conjq(a)
271 #define PetscSqrtComplex(a)          csqrtq(a)
272 #define PetscPowComplex(a,b)         cpowq(a,b)
273 #define PetscExpComplex(a)           cexpq(a)
274 #define PetscLogComplex(a)           clogq(a)
275 #define PetscSinComplex(a)           csinq(a)
276 #define PetscCosComplex(a)           ccosq(a)
277 #define PetscAsinComplex(a)          casinq(a)
278 #define PetscAcosComplex(a)          cacosq(a)
279 #define PetscTanComplex(a)           ctanq(a)
280 #define PetscSinhComplex(a)          csinhq(a)
281 #define PetscCoshComplex(a)          ccoshq(a)
282 #define PetscTanhComplex(a)          ctanhq(a)
283 
284 #endif /* PETSC_USE_REAL_* */
285 #elif (defined(PETSC_USE_COMPLEX) && !defined(PETSC_SKIP_COMPLEX))
286 #error "PETSc was configured --with-scalar-type=complex, but a language-appropriate complex library is not available"
287 #endif /* !PETSC_SKIP_COMPLEX */
288 #endif /* (__cplusplus && PETSC_HAVE_CXX_COMPLEX) else-if (!__cplusplus && PETSC_HAVE_C99_COMPLEX) */
289 
290 #if defined(PETSC_HAVE_COMPLEX)
291 #if defined(PETSC_HAVE_MPI_C_DOUBLE_COMPLEX)
292 #define MPIU_C_DOUBLE_COMPLEX MPI_C_DOUBLE_COMPLEX
293 #define MPIU_C_COMPLEX MPI_C_COMPLEX
294 #else
295 # if defined(__cplusplus) && defined(PETSC_HAVE_CXX_COMPLEX) && !defined(PETSC_USE_REAL___FLOAT128)
296   typedef complexlib::complex<double> petsc_mpiu_c_double_complex;
297   typedef complexlib::complex<float> petsc_mpiu_c_complex;
298 # elif !defined(__cplusplus) && defined(PETSC_HAVE_C99_COMPLEX)
299   typedef double _Complex petsc_mpiu_c_double_complex;
300   typedef float _Complex petsc_mpiu_c_complex;
301 # else
302   typedef struct {double real,imag;} petsc_mpiu_c_double_complex;
303   typedef struct {float real,imag;} petsc_mpiu_c_complex;
304 # endif
305 PETSC_EXTERN MPI_Datatype MPIU_C_DOUBLE_COMPLEX PetscAttrMPITypeTagLayoutCompatible(petsc_mpiu_c_double_complex);
306 PETSC_EXTERN MPI_Datatype MPIU_C_COMPLEX PetscAttrMPITypeTagLayoutCompatible(petsc_mpiu_c_complex);
307 #endif /* PETSC_HAVE_MPI_C_DOUBLE_COMPLEX */
308 #endif /* PETSC_HAVE_COMPLEX */
309 
310 #if defined(PETSC_HAVE_COMPLEX)
311 #  if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
312 #    define MPIU_COMPLEX MPIU_C_COMPLEX
313 #  elif defined(PETSC_USE_REAL_DOUBLE)
314 #    define MPIU_COMPLEX MPIU_C_DOUBLE_COMPLEX
315 #  elif defined(PETSC_USE_REAL___FLOAT128)
316 #    define MPIU_COMPLEX MPIU___COMPLEX128
317 #  endif /* PETSC_USE_REAL_* */
318 #endif
319 
320 #if (defined(PETSC_USE_COMPLEX) && !defined(PETSC_SKIP_COMPLEX))
321 typedef PetscComplex PetscScalar;
322 #define PetscRealPart(a)      PetscRealPartComplex(a)
323 #define PetscImaginaryPart(a) PetscImaginaryPartComplex(a)
324 #define PetscAbsScalar(a)     PetscAbsComplex(a)
325 #define PetscConj(a)          PetscConjComplex(a)
326 #define PetscSqrtScalar(a)    PetscSqrtComplex(a)
327 #define PetscPowScalar(a,b)   PetscPowComplex(a,b)
328 #define PetscExpScalar(a)     PetscExpComplex(a)
329 #define PetscLogScalar(a)     PetscLogComplex(a)
330 #define PetscSinScalar(a)     PetscSinComplex(a)
331 #define PetscCosScalar(a)     PetscCosComplex(a)
332 #define PetscAsinScalar(a)    PetscAsinComplex(a)
333 #define PetscAcosScalar(a)    PetscAcosComplex(a)
334 #define PetscTanScalar(a)     PetscTanComplex(a)
335 #define PetscSinhScalar(a)    PetscSinhComplex(a)
336 #define PetscCoshScalar(a)    PetscCoshComplex(a)
337 #define PetscTanhScalar(a)    PetscTanhComplex(a)
338 #define MPIU_SCALAR MPIU_COMPLEX
339 
340 /*
341     real number definitions
342  */
343 #else /* PETSC_USE_COMPLEX */
344 typedef PetscReal PetscScalar;
345 #define MPIU_SCALAR MPIU_REAL
346 
347 #define PetscRealPart(a)      (a)
348 #define PetscImaginaryPart(a) ((PetscReal)0.)
349 PETSC_STATIC_INLINE PetscReal PetscAbsScalar(PetscScalar a) {return a < 0.0 ? -a : a;}
350 #define PetscConj(a)          (a)
351 #if !defined(PETSC_USE_REAL___FLOAT128)  && !defined(PETSC_USE_REAL___FP16)
352 #define PetscSqrtScalar(a)    sqrt(a)
353 #define PetscPowScalar(a,b)   pow(a,b)
354 #define PetscExpScalar(a)     exp(a)
355 #define PetscLogScalar(a)     log(a)
356 #define PetscSinScalar(a)     sin(a)
357 #define PetscCosScalar(a)     cos(a)
358 #define PetscAsinScalar(a)    asin(a)
359 #define PetscAcosScalar(a)    acos(a)
360 #define PetscTanScalar(a)     tan(a)
361 #define PetscSinhScalar(a)    sinh(a)
362 #define PetscCoshScalar(a)    cosh(a)
363 #define PetscTanhScalar(a)    tanh(a)
364 #elif defined(PETSC_USE_REAL___FP16)
365 #define PetscSqrtScalar(a)    sqrtf(a)
366 #define PetscPowScalar(a,b)   powf(a,b)
367 #define PetscExpScalar(a)     expf(a)
368 #define PetscLogScalar(a)     logf(a)
369 #define PetscSinScalar(a)     sinf(a)
370 #define PetscCosScalar(a)     cosf(a)
371 #define PetscAsinScalar(a)    asinf(a)
372 #define PetscAcosScalar(a)    acosf(a)
373 #define PetscTanScalar(a)     tanf(a)
374 #define PetscSinhScalar(a)    sinhf(a)
375 #define PetscCoshScalar(a)    coshf(a)
376 #define PetscTanhScalar(a)    tanhf(a)
377 #else /* PETSC_USE_REAL___FLOAT128 */
378 #define PetscSqrtScalar(a)    sqrtq(a)
379 #define PetscPowScalar(a,b)   powq(a,b)
380 #define PetscExpScalar(a)     expq(a)
381 #define PetscLogScalar(a)     logq(a)
382 #define PetscSinScalar(a)     sinq(a)
383 #define PetscCosScalar(a)     cosq(a)
384 #define PetscAsinScalar(a)    asinq(a)
385 #define PetscAcosScalar(a)    acosq(a)
386 #define PetscTanScalar(a)     tanq(a)
387 #define PetscSinhScalar(a)    sinhq(a)
388 #define PetscCoshScalar(a)    coshq(a)
389 #define PetscTanhScalar(a)    tanhq(a)
390 #endif /* PETSC_USE_REAL___FLOAT128 */
391 
392 #endif /* PETSC_USE_COMPLEX */
393 
394 #define PetscSign(a) (((a) >= 0) ? ((a) == 0 ? 0 : 1) : -1)
395 #define PetscSignReal(a) (((a) >= 0.0) ? ((a) == 0.0 ? 0.0 : 1.0) : -1.0)
396 #define PetscAbs(a)  (((a) >= 0) ? (a) : (-(a)))
397 
398 /* --------------------------------------------------------------------------*/
399 
400 /*
401    Certain objects may be created using either single or double precision.
402    This is currently not used.
403 */
404 typedef enum { PETSC_SCALAR_DOUBLE,PETSC_SCALAR_SINGLE, PETSC_SCALAR_LONG_DOUBLE, PETSC_SCALAR_HALF } PetscScalarPrecision;
405 
406 #if defined(PETSC_HAVE_COMPLEX)
407 /* PETSC_i is the imaginary number, i */
408 PETSC_EXTERN PetscComplex PETSC_i;
409 
410 /* Try to do the right thing for complex number construction: see
411 
412   http://www.open-std.org/jtc1/sc22/wg14/www/docs/n1464.htm
413 
414   for details
415 */
416 PETSC_STATIC_INLINE PetscComplex PetscCMPLX(PetscReal x, PetscReal y)
417 {
418 #if   defined(__cplusplus) && !defined(PETSC_USE_REAL___FLOAT128)
419   return PetscComplex(x,y);
420 #elif defined(_Imaginary_I)
421   return x + y * _Imaginary_I;
422 #else
423   { /* In both C99 and C11 (ISO/IEC 9899, Section 6.2.5),
424 
425        "For each floating type there is a corresponding real type, which is always a real floating
426        type. For real floating types, it is the same type. For complex types, it is the type given
427        by deleting the keyword _Complex from the type name."
428 
429        So type punning should be portable. */
430     union { PetscComplex z; PetscReal f[2]; } uz;
431 
432     uz.f[0] = x;
433     uz.f[1] = y;
434     return uz.z;
435   }
436 #endif
437 }
438 #endif
439 
440 
441 /*MC
442    PetscMin - Returns minimum of two numbers
443 
444    Synopsis:
445    #include <petscmath.h>
446    type PetscMin(type v1,type v2)
447 
448    Not Collective
449 
450    Input Parameter:
451 +  v1 - first value to find minimum of
452 -  v2 - second value to find minimum of
453 
454    Notes: type can be integer or floating point value
455 
456    Level: beginner
457 
458 .seealso: PetscMax(), PetscClipInterval(), PetscAbsInt(), PetscAbsReal(), PetscSqr()
459 
460 M*/
461 #define PetscMin(a,b)   (((a)<(b)) ?  (a) : (b))
462 
463 /*MC
464    PetscMax - Returns maxium of two numbers
465 
466    Synopsis:
467    #include <petscmath.h>
468    type max PetscMax(type v1,type v2)
469 
470    Not Collective
471 
472    Input Parameter:
473 +  v1 - first value to find maximum of
474 -  v2 - second value to find maximum of
475 
476    Notes: type can be integer or floating point value
477 
478    Level: beginner
479 
480 .seealso: PetscMin(), PetscClipInterval(), PetscAbsInt(), PetscAbsReal(), PetscSqr()
481 
482 M*/
483 #define PetscMax(a,b)   (((a)<(b)) ?  (b) : (a))
484 
485 /*MC
486    PetscClipInterval - Returns a number clipped to be within an interval
487 
488    Synopsis:
489    #include <petscmath.h>
490    type clip PetscClipInterval(type x,type a,type b)
491 
492    Not Collective
493 
494    Input Parameter:
495 +  x - value to use if within interval (a,b)
496 .  a - lower end of interval
497 -  b - upper end of interval
498 
499    Notes: type can be integer or floating point value
500 
501    Level: beginner
502 
503 .seealso: PetscMin(), PetscMax(), PetscAbsInt(), PetscAbsReal(), PetscSqr()
504 
505 M*/
506 #define PetscClipInterval(x,a,b)   (PetscMax((a),PetscMin((x),(b))))
507 
508 /*MC
509    PetscAbsInt - Returns the absolute value of an integer
510 
511    Synopsis:
512    #include <petscmath.h>
513    int abs PetscAbsInt(int v1)
514 
515    Not Collective
516 
517    Input Parameter:
518 .   v1 - the integer
519 
520    Level: beginner
521 
522 .seealso: PetscMax(), PetscMin(), PetscAbsReal(), PetscSqr()
523 
524 M*/
525 #define PetscAbsInt(a)  (((a)<0)   ? (-(a)) : (a))
526 
527 /*MC
528    PetscAbsReal - Returns the absolute value of an real number
529 
530    Synopsis:
531    #include <petscmath.h>
532    Real abs PetscAbsReal(PetscReal v1)
533 
534    Not Collective
535 
536    Input Parameter:
537 .   v1 - the double
538 
539 
540    Level: beginner
541 
542 .seealso: PetscMax(), PetscMin(), PetscAbsInt(), PetscSqr()
543 
544 M*/
545 #if defined(PETSC_USE_REAL_SINGLE)
546 #define PetscAbsReal(a) fabsf(a)
547 #elif defined(PETSC_USE_REAL_DOUBLE)
548 #define PetscAbsReal(a) fabs(a)
549 #elif defined(PETSC_USE_REAL___FLOAT128)
550 #define PetscAbsReal(a) fabsq(a)
551 #elif defined(PETSC_USE_REAL___FP16)
552 #define PetscAbsReal(a) fabsf(a)
553 #endif
554 
555 /*MC
556    PetscSqr - Returns the square of a number
557 
558    Synopsis:
559    #include <petscmath.h>
560    type sqr PetscSqr(type v1)
561 
562    Not Collective
563 
564    Input Parameter:
565 .   v1 - the value
566 
567    Notes: type can be integer or floating point value
568 
569    Level: beginner
570 
571 .seealso: PetscMax(), PetscMin(), PetscAbsInt(), PetscAbsReal()
572 
573 M*/
574 #define PetscSqr(a)     ((a)*(a))
575 
576 /* ----------------------------------------------------------------------------*/
577 
578 #if defined(PETSC_USE_REAL_SINGLE)
579 #define PetscRealConstant(constant) constant##F
580 #elif defined(PETSC_USE_REAL___FLOAT128)
581 #define PetscRealConstant(constant) constant##Q
582 #else
583 #define PetscRealConstant(constant) constant
584 #endif
585 
586 /*
587      Basic constants
588 */
589 #define PETSC_PI   PetscRealConstant(3.1415926535897932384626433832795029)
590 #define PETSC_PHI  PetscRealConstant(1.6180339887498948482045868343656381)
591 
592 #if !defined(PETSC_USE_64BIT_INDICES)
593 #define PETSC_MAX_INT            2147483647
594 #define PETSC_MIN_INT            (-PETSC_MAX_INT - 1)
595 #else
596 #define PETSC_MAX_INT            9223372036854775807L
597 #define PETSC_MIN_INT            (-PETSC_MAX_INT - 1)
598 #endif
599 
600 #if defined(PETSC_USE_REAL_SINGLE)
601 #  define PETSC_MAX_REAL                3.40282346638528860e+38F
602 #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
603 #  define PETSC_MACHINE_EPSILON         1.19209290e-07F
604 #  define PETSC_SQRT_MACHINE_EPSILON    3.45266983e-04F
605 #  define PETSC_SMALL                   1.e-5F
606 #elif defined(PETSC_USE_REAL_DOUBLE)
607 #  define PETSC_MAX_REAL                1.7976931348623157e+308
608 #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
609 #  define PETSC_MACHINE_EPSILON         2.2204460492503131e-16
610 #  define PETSC_SQRT_MACHINE_EPSILON    1.490116119384766e-08
611 #  define PETSC_SMALL                   1.e-10
612 #elif defined(PETSC_USE_REAL___FLOAT128)
613 #  define PETSC_MAX_REAL                FLT128_MAX
614 #  define PETSC_MIN_REAL                (-FLT128_MAX)
615 #  define PETSC_MACHINE_EPSILON         FLT128_EPSILON
616 #  define PETSC_SQRT_MACHINE_EPSILON    1.38777878078144567552953958511352539e-17Q
617 #  define PETSC_SMALL                   1.e-20Q
618 #elif defined(PETSC_USE_REAL___FP16)  /* maybe should use single precision values for these? */
619 #  define PETSC_MAX_REAL                65504.
620 #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
621 #  define PETSC_MACHINE_EPSILON         .00097656
622 #  define PETSC_SQRT_MACHINE_EPSILON    .0312
623 #  define PETSC_SMALL                   5.e-3
624 #endif
625 
626 #define PETSC_INFINITY               (PETSC_MAX_REAL/4)
627 #define PETSC_NINFINITY              (-PETSC_INFINITY)
628 
629 PETSC_EXTERN PetscBool PetscIsInfReal(PetscReal);
630 PETSC_EXTERN PetscBool PetscIsNanReal(PetscReal);
631 PETSC_EXTERN PetscBool PetscIsNormalReal(PetscReal);
632 PETSC_STATIC_INLINE PetscBool PetscIsInfOrNanReal(PetscReal v) {return PetscIsInfReal(v) || PetscIsNanReal(v) ? PETSC_TRUE : PETSC_FALSE;}
633 PETSC_STATIC_INLINE PetscBool PetscIsInfScalar(PetscScalar v) {return PetscIsInfReal(PetscAbsScalar(v));}
634 PETSC_STATIC_INLINE PetscBool PetscIsNanScalar(PetscScalar v) {return PetscIsNanReal(PetscAbsScalar(v));}
635 PETSC_STATIC_INLINE PetscBool PetscIsInfOrNanScalar(PetscScalar v) {return PetscIsInfOrNanReal(PetscAbsScalar(v));}
636 PETSC_STATIC_INLINE PetscBool PetscIsNormalScalar(PetscScalar v) {return PetscIsNormalReal(PetscAbsScalar(v));}
637 
638 PETSC_EXTERN PetscBool PetscEqualReal(PetscReal,PetscReal);
639 PETSC_EXTERN PetscBool PetscEqualScalar(PetscScalar,PetscScalar);
640 
641 /*
642     These macros are currently hardwired to match the regular data types, so there is no support for a different
643     MatScalar from PetscScalar. We left the MatScalar in the source just in case we use it again.
644  */
645 #define MPIU_MATSCALAR MPIU_SCALAR
646 typedef PetscScalar MatScalar;
647 typedef PetscReal MatReal;
648 
649 struct petsc_mpiu_2scalar {PetscScalar a,b;};
650 PETSC_EXTERN MPI_Datatype MPIU_2SCALAR PetscAttrMPITypeTagLayoutCompatible(struct petsc_mpiu_2scalar);
651 #if defined(PETSC_USE_64BIT_INDICES) || !defined(MPI_2INT)
652 struct petsc_mpiu_2int {PetscInt a,b;};
653 PETSC_EXTERN MPI_Datatype MPIU_2INT PetscAttrMPITypeTagLayoutCompatible(struct petsc_mpiu_2int);
654 #else
655 #define MPIU_2INT MPI_2INT
656 #endif
657 
658 PETSC_STATIC_INLINE PetscInt PetscPowInt(PetscInt base,PetscInt power)
659 {
660   PetscInt result = 1;
661   while (power) {
662     if (power & 1) result *= base;
663     power >>= 1;
664     base *= base;
665   }
666   return result;
667 }
668 
669 PETSC_STATIC_INLINE PetscReal PetscPowRealInt(PetscReal base,PetscInt power)
670 {
671   PetscReal result = 1;
672   if (power < 0) {
673     power = -power;
674     base  = ((PetscReal)1)/base;
675   }
676   while (power) {
677     if (power & 1) result *= base;
678     power >>= 1;
679     base *= base;
680   }
681   return result;
682 }
683 
684 PETSC_STATIC_INLINE PetscScalar PetscPowScalarInt(PetscScalar base,PetscInt power)
685 {
686   PetscScalar result = 1;
687   if (power < 0) {
688     power = -power;
689     base  = ((PetscReal)1)/base;
690   }
691   while (power) {
692     if (power & 1) result *= base;
693     power >>= 1;
694     base *= base;
695   }
696   return result;
697 }
698 
699 PETSC_STATIC_INLINE PetscScalar PetscPowScalarReal(PetscScalar base,PetscReal power)
700 {
701   PetscScalar cpower = power;
702   return PetscPowScalar(base,cpower);
703 }
704 
705 #ifndef PETSC_HAVE_LOG2
706 PETSC_STATIC_INLINE PetscReal PetscLog2Real(PetscReal n)
707 {
708   return PetscLogReal(n)/PetscLogReal(2.0);
709 }
710 #endif
711 #endif
712