xref: /petsc/include/petscmath.h (revision 6ffe77eaecce1557e50d00ca5347a7f48e598865)
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 #if !defined(PETSCMATH_H)
11 #define PETSCMATH_H
12 
13 #include <math.h>
14 #include <petscsystypes.h>
15 
16 /* SUBMANSEC = Sys */
17 
18 /*
19 
20    Defines operations that are different for complex and real numbers.
21    All PETSc objects in one program are built around the object
22    PetscScalar which is either always a real or a complex.
23 
24 */
25 
26 /*
27     Real number definitions
28  */
29 #if defined(PETSC_USE_REAL_SINGLE)
30 #define PetscSqrtReal(a)       sqrtf(a)
31 #define PetscCbrtReal(a)       cbrtf(a)
32 #define PetscHypotReal(a,b)    hypotf(a,b)
33 #define PetscAtan2Real(a,b)    atan2f(a,b)
34 #define PetscPowReal(a,b)      powf(a,b)
35 #define PetscExpReal(a)        expf(a)
36 #define PetscLogReal(a)        logf(a)
37 #define PetscLog10Real(a)      log10f(a)
38 #define PetscLog2Real(a)       log2f(a)
39 #define PetscSinReal(a)        sinf(a)
40 #define PetscCosReal(a)        cosf(a)
41 #define PetscTanReal(a)        tanf(a)
42 #define PetscAsinReal(a)       asinf(a)
43 #define PetscAcosReal(a)       acosf(a)
44 #define PetscAtanReal(a)       atanf(a)
45 #define PetscSinhReal(a)       sinhf(a)
46 #define PetscCoshReal(a)       coshf(a)
47 #define PetscTanhReal(a)       tanhf(a)
48 #define PetscAsinhReal(a)      asinhf(a)
49 #define PetscAcoshReal(a)      acoshf(a)
50 #define PetscAtanhReal(a)      atanhf(a)
51 #define PetscErfReal(a)        erff(a)
52 #define PetscCeilReal(a)       ceilf(a)
53 #define PetscFloorReal(a)      floorf(a)
54 #define PetscFmodReal(a,b)     fmodf(a,b)
55 #define PetscCopysignReal(a,b) copysignf(a,b)
56 #define PetscTGamma(a)         tgammaf(a)
57 #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
58 #define PetscLGamma(a)         gammaf(a)
59 #else
60 #define PetscLGamma(a)         lgammaf(a)
61 #endif
62 
63 #elif defined(PETSC_USE_REAL_DOUBLE)
64 #define PetscSqrtReal(a)       sqrt(a)
65 #define PetscCbrtReal(a)       cbrt(a)
66 #define PetscHypotReal(a,b)    hypot(a,b)
67 #define PetscAtan2Real(a,b)    atan2(a,b)
68 #define PetscPowReal(a,b)      pow(a,b)
69 #define PetscExpReal(a)        exp(a)
70 #define PetscLogReal(a)        log(a)
71 #define PetscLog10Real(a)      log10(a)
72 #define PetscLog2Real(a)       log2(a)
73 #define PetscSinReal(a)        sin(a)
74 #define PetscCosReal(a)        cos(a)
75 #define PetscTanReal(a)        tan(a)
76 #define PetscAsinReal(a)       asin(a)
77 #define PetscAcosReal(a)       acos(a)
78 #define PetscAtanReal(a)       atan(a)
79 #define PetscSinhReal(a)       sinh(a)
80 #define PetscCoshReal(a)       cosh(a)
81 #define PetscTanhReal(a)       tanh(a)
82 #define PetscAsinhReal(a)      asinh(a)
83 #define PetscAcoshReal(a)      acosh(a)
84 #define PetscAtanhReal(a)      atanh(a)
85 #define PetscErfReal(a)        erf(a)
86 #define PetscCeilReal(a)       ceil(a)
87 #define PetscFloorReal(a)      floor(a)
88 #define PetscFmodReal(a,b)     fmod(a,b)
89 #define PetscCopysignReal(a,b) copysign(a,b)
90 #define PetscTGamma(a)         tgamma(a)
91 #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
92 #define PetscLGamma(a)         gamma(a)
93 #else
94 #define PetscLGamma(a)         lgamma(a)
95 #endif
96 
97 #elif defined(PETSC_USE_REAL___FLOAT128)
98 #define PetscSqrtReal(a)       sqrtq(a)
99 #define PetscCbrtReal(a)       cbrtq(a)
100 #define PetscHypotReal(a,b)    hypotq(a,b)
101 #define PetscAtan2Real(a,b)    atan2q(a,b)
102 #define PetscPowReal(a,b)      powq(a,b)
103 #define PetscExpReal(a)        expq(a)
104 #define PetscLogReal(a)        logq(a)
105 #define PetscLog10Real(a)      log10q(a)
106 #define PetscLog2Real(a)       log2q(a)
107 #define PetscSinReal(a)        sinq(a)
108 #define PetscCosReal(a)        cosq(a)
109 #define PetscTanReal(a)        tanq(a)
110 #define PetscAsinReal(a)       asinq(a)
111 #define PetscAcosReal(a)       acosq(a)
112 #define PetscAtanReal(a)       atanq(a)
113 #define PetscSinhReal(a)       sinhq(a)
114 #define PetscCoshReal(a)       coshq(a)
115 #define PetscTanhReal(a)       tanhq(a)
116 #define PetscAsinhReal(a)      asinhq(a)
117 #define PetscAcoshReal(a)      acoshq(a)
118 #define PetscAtanhReal(a)      atanhq(a)
119 #define PetscErfReal(a)        erfq(a)
120 #define PetscCeilReal(a)       ceilq(a)
121 #define PetscFloorReal(a)      floorq(a)
122 #define PetscFmodReal(a,b)     fmodq(a,b)
123 #define PetscCopysignReal(a,b) copysignq(a,b)
124 #define PetscTGamma(a)         tgammaq(a)
125 #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
126 #define PetscLGamma(a)         gammaq(a)
127 #else
128 #define PetscLGamma(a)         lgammaq(a)
129 #endif
130 
131 #elif defined(PETSC_USE_REAL___FP16)
132 #define PetscSqrtReal(a)       sqrtf(a)
133 #define PetscCbrtReal(a)       cbrtf(a)
134 #define PetscHypotReal(a,b)    hypotf(a,b)
135 #define PetscAtan2Real(a,b)    atan2f(a,b)
136 #define PetscPowReal(a,b)      powf(a,b)
137 #define PetscExpReal(a)        expf(a)
138 #define PetscLogReal(a)        logf(a)
139 #define PetscLog10Real(a)      log10f(a)
140 #define PetscLog2Real(a)       log2f(a)
141 #define PetscSinReal(a)        sinf(a)
142 #define PetscCosReal(a)        cosf(a)
143 #define PetscTanReal(a)        tanf(a)
144 #define PetscAsinReal(a)       asinf(a)
145 #define PetscAcosReal(a)       acosf(a)
146 #define PetscAtanReal(a)       atanf(a)
147 #define PetscSinhReal(a)       sinhf(a)
148 #define PetscCoshReal(a)       coshf(a)
149 #define PetscTanhReal(a)       tanhf(a)
150 #define PetscAsinhReal(a)      asinhf(a)
151 #define PetscAcoshReal(a)      acoshf(a)
152 #define PetscAtanhReal(a)      atanhf(a)
153 #define PetscErfReal(a)        erff(a)
154 #define PetscCeilReal(a)       ceilf(a)
155 #define PetscFloorReal(a)      floorf(a)
156 #define PetscFmodReal(a,b)     fmodf(a,b)
157 #define PetscCopySignReal(a,b) copysignf(a,b)
158 #define PetscTGamma(a)         tgammaf(a)
159 #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
160 #define PetscLGamma(a)         gammaf(a)
161 #else
162 #define PetscLGamma(a)         lgammaf(a)
163 #endif
164 
165 #endif /* PETSC_USE_REAL_* */
166 
167 static inline PetscReal PetscSignReal(PetscReal a)
168 {
169   return (PetscReal)((a < (PetscReal)0) ? -1 : ((a > (PetscReal)0) ? 1 : 0));
170 }
171 
172 #if !defined(PETSC_HAVE_LOG2)
173 #undef PetscLog2Real
174 static inline PetscReal PetscLog2Real(PetscReal a)
175 {
176   return PetscLogReal(a)/PetscLogReal((PetscReal)2);
177 }
178 #endif
179 
180 #if defined(PETSC_USE_REAL___FLOAT128)
181 PETSC_EXTERN MPI_Datatype MPIU___FLOAT128 PetscAttrMPITypeTag(__float128);
182 #endif
183 #if defined(PETSC_USE_REAL___FP16)
184 PETSC_EXTERN MPI_Datatype MPIU___FP16 PetscAttrMPITypeTag(__fp16);
185 #endif
186 
187 /*MC
188    MPIU_REAL - MPI datatype corresponding to PetscReal
189 
190    Notes:
191    In MPI calls that require an MPI datatype that matches a PetscReal or array of PetscReal values, pass this value.
192 
193    Level: beginner
194 
195 .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_SCALAR`, `MPIU_COMPLEX`, `MPIU_INT`
196 M*/
197 #if defined(PETSC_USE_REAL_SINGLE)
198 #  define MPIU_REAL MPI_FLOAT
199 #elif defined(PETSC_USE_REAL_DOUBLE)
200 #  define MPIU_REAL MPI_DOUBLE
201 #elif defined(PETSC_USE_REAL___FLOAT128)
202 #  define MPIU_REAL MPIU___FLOAT128
203 #elif defined(PETSC_USE_REAL___FP16)
204 #  define MPIU_REAL MPIU___FP16
205 #endif /* PETSC_USE_REAL_* */
206 
207 /*
208     Complex number definitions
209  */
210 #if defined(PETSC_HAVE_COMPLEX)
211 #if defined(__cplusplus) && !defined(PETSC_USE_REAL___FLOAT128)
212 /* C++ support of complex number */
213 
214 #define PetscRealPartComplex(a)      (a).real()
215 #define PetscImaginaryPartComplex(a) (a).imag()
216 #define PetscAbsComplex(a)           petsccomplexlib::abs(a)
217 #define PetscArgComplex(a)           petsccomplexlib::arg(a)
218 #define PetscConjComplex(a)          petsccomplexlib::conj(a)
219 #define PetscSqrtComplex(a)          petsccomplexlib::sqrt(a)
220 #define PetscPowComplex(a,b)         petsccomplexlib::pow(a,b)
221 #define PetscExpComplex(a)           petsccomplexlib::exp(a)
222 #define PetscLogComplex(a)           petsccomplexlib::log(a)
223 #define PetscSinComplex(a)           petsccomplexlib::sin(a)
224 #define PetscCosComplex(a)           petsccomplexlib::cos(a)
225 #define PetscTanComplex(a)           petsccomplexlib::tan(a)
226 #define PetscAsinComplex(a)          petsccomplexlib::asin(a)
227 #define PetscAcosComplex(a)          petsccomplexlib::acos(a)
228 #define PetscAtanComplex(a)          petsccomplexlib::atan(a)
229 #define PetscSinhComplex(a)          petsccomplexlib::sinh(a)
230 #define PetscCoshComplex(a)          petsccomplexlib::cosh(a)
231 #define PetscTanhComplex(a)          petsccomplexlib::tanh(a)
232 #define PetscAsinhComplex(a)         petsccomplexlib::asinh(a)
233 #define PetscAcoshComplex(a)         petsccomplexlib::acosh(a)
234 #define PetscAtanhComplex(a)         petsccomplexlib::atanh(a)
235 
236 /* TODO: Add configure tests
237 
238 #if !defined(PETSC_HAVE_CXX_TAN_COMPLEX)
239 #undef PetscTanComplex
240 static inline PetscComplex PetscTanComplex(PetscComplex z)
241 {
242   return PetscSinComplex(z)/PetscCosComplex(z);
243 }
244 #endif
245 
246 #if !defined(PETSC_HAVE_CXX_TANH_COMPLEX)
247 #undef PetscTanhComplex
248 static inline PetscComplex PetscTanhComplex(PetscComplex z)
249 {
250   return PetscSinhComplex(z)/PetscCoshComplex(z);
251 }
252 #endif
253 
254 #if !defined(PETSC_HAVE_CXX_ASIN_COMPLEX)
255 #undef PetscAsinComplex
256 static inline PetscComplex PetscAsinComplex(PetscComplex z)
257 {
258   const PetscComplex j(0,1);
259   return -j*PetscLogComplex(j*z+PetscSqrtComplex(1.0f-z*z));
260 }
261 #endif
262 
263 #if !defined(PETSC_HAVE_CXX_ACOS_COMPLEX)
264 #undef PetscAcosComplex
265 static inline PetscComplex PetscAcosComplex(PetscComplex z)
266 {
267   const PetscComplex j(0,1);
268   return j*PetscLogComplex(z-j*PetscSqrtComplex(1.0f-z*z));
269 }
270 #endif
271 
272 #if !defined(PETSC_HAVE_CXX_ATAN_COMPLEX)
273 #undef PetscAtanComplex
274 static inline PetscComplex PetscAtanComplex(PetscComplex z)
275 {
276   const PetscComplex j(0,1);
277   return 0.5f*j*PetscLogComplex((1.0f-j*z)/(1.0f+j*z));
278 }
279 #endif
280 
281 #if !defined(PETSC_HAVE_CXX_ASINH_COMPLEX)
282 #undef PetscAsinhComplex
283 static inline PetscComplex PetscAsinhComplex(PetscComplex z)
284 {
285   return PetscLogComplex(z+PetscSqrtComplex(z*z+1.0f));
286 }
287 #endif
288 
289 #if !defined(PETSC_HAVE_CXX_ACOSH_COMPLEX)
290 #undef PetscAcoshComplex
291 static inline PetscComplex PetscAcoshComplex(PetscComplex z)
292 {
293   return PetscLogComplex(z+PetscSqrtComplex(z*z-1.0f));
294 }
295 #endif
296 
297 #if !defined(PETSC_HAVE_CXX_ATANH_COMPLEX)
298 #undef PetscAtanhComplex
299 static inline PetscComplex PetscAtanhComplex(PetscComplex z)
300 {
301   return 0.5f*PetscLogComplex((1.0f+z)/(1.0f-z));
302 }
303 #endif
304 
305 */
306 
307 #else /* C99 support of complex number */
308 
309 #if defined(PETSC_USE_REAL_SINGLE)
310 #define PetscRealPartComplex(a)      crealf(a)
311 #define PetscImaginaryPartComplex(a) cimagf(a)
312 #define PetscAbsComplex(a)           cabsf(a)
313 #define PetscArgComplex(a)           cargf(a)
314 #define PetscConjComplex(a)          conjf(a)
315 #define PetscSqrtComplex(a)          csqrtf(a)
316 #define PetscPowComplex(a,b)         cpowf(a,b)
317 #define PetscExpComplex(a)           cexpf(a)
318 #define PetscLogComplex(a)           clogf(a)
319 #define PetscSinComplex(a)           csinf(a)
320 #define PetscCosComplex(a)           ccosf(a)
321 #define PetscTanComplex(a)           ctanf(a)
322 #define PetscAsinComplex(a)          casinf(a)
323 #define PetscAcosComplex(a)          cacosf(a)
324 #define PetscAtanComplex(a)          catanf(a)
325 #define PetscSinhComplex(a)          csinhf(a)
326 #define PetscCoshComplex(a)          ccoshf(a)
327 #define PetscTanhComplex(a)          ctanhf(a)
328 #define PetscAsinhComplex(a)         casinhf(a)
329 #define PetscAcoshComplex(a)         cacoshf(a)
330 #define PetscAtanhComplex(a)         catanhf(a)
331 
332 #elif defined(PETSC_USE_REAL_DOUBLE)
333 #define PetscRealPartComplex(a)      creal(a)
334 #define PetscImaginaryPartComplex(a) cimag(a)
335 #define PetscAbsComplex(a)           cabs(a)
336 #define PetscArgComplex(a)           carg(a)
337 #define PetscConjComplex(a)          conj(a)
338 #define PetscSqrtComplex(a)          csqrt(a)
339 #define PetscPowComplex(a,b)         cpow(a,b)
340 #define PetscExpComplex(a)           cexp(a)
341 #define PetscLogComplex(a)           clog(a)
342 #define PetscSinComplex(a)           csin(a)
343 #define PetscCosComplex(a)           ccos(a)
344 #define PetscTanComplex(a)           ctan(a)
345 #define PetscAsinComplex(a)          casin(a)
346 #define PetscAcosComplex(a)          cacos(a)
347 #define PetscAtanComplex(a)          catan(a)
348 #define PetscSinhComplex(a)          csinh(a)
349 #define PetscCoshComplex(a)          ccosh(a)
350 #define PetscTanhComplex(a)          ctanh(a)
351 #define PetscAsinhComplex(a)         casinh(a)
352 #define PetscAcoshComplex(a)         cacosh(a)
353 #define PetscAtanhComplex(a)         catanh(a)
354 
355 #elif defined(PETSC_USE_REAL___FLOAT128)
356 #define PetscRealPartComplex(a)      crealq(a)
357 #define PetscImaginaryPartComplex(a) cimagq(a)
358 #define PetscAbsComplex(a)           cabsq(a)
359 #define PetscArgComplex(a)           cargq(a)
360 #define PetscConjComplex(a)          conjq(a)
361 #define PetscSqrtComplex(a)          csqrtq(a)
362 #define PetscPowComplex(a,b)         cpowq(a,b)
363 #define PetscExpComplex(a)           cexpq(a)
364 #define PetscLogComplex(a)           clogq(a)
365 #define PetscSinComplex(a)           csinq(a)
366 #define PetscCosComplex(a)           ccosq(a)
367 #define PetscTanComplex(a)           ctanq(a)
368 #define PetscAsinComplex(a)          casinq(a)
369 #define PetscAcosComplex(a)          cacosq(a)
370 #define PetscAtanComplex(a)          catanq(a)
371 #define PetscSinhComplex(a)          csinhq(a)
372 #define PetscCoshComplex(a)          ccoshq(a)
373 #define PetscTanhComplex(a)          ctanhq(a)
374 #define PetscAsinhComplex(a)         casinhq(a)
375 #define PetscAcoshComplex(a)         cacoshq(a)
376 #define PetscAtanhComplex(a)         catanhq(a)
377 
378 #endif /* PETSC_USE_REAL_* */
379 #endif /* (__cplusplus) */
380 
381 /*
382    PETSC_i is the imaginary number, i
383 */
384 PETSC_EXTERN PetscComplex PETSC_i;
385 
386 /*
387    Try to do the right thing for complex number construction: see
388    http://www.open-std.org/jtc1/sc22/wg14/www/docs/n1464.htm
389    for details
390 */
391 static inline PetscComplex PetscCMPLX(PetscReal x, PetscReal y)
392 {
393 #if defined(__cplusplus) && !defined(PETSC_USE_REAL___FLOAT128)
394   return PetscComplex(x,y);
395 #elif defined(_Imaginary_I)
396   return x + y * _Imaginary_I;
397 #else
398   { /* In both C99 and C11 (ISO/IEC 9899, Section 6.2.5),
399 
400        "For each floating type there is a corresponding real type, which is always a real floating
401        type. For real floating types, it is the same type. For complex types, it is the type given
402        by deleting the keyword _Complex from the type name."
403 
404        So type punning should be portable. */
405     union { PetscComplex z; PetscReal f[2]; } uz;
406 
407     uz.f[0] = x;
408     uz.f[1] = y;
409     return uz.z;
410   }
411 #endif
412 }
413 
414 #define MPIU_C_COMPLEX MPI_C_COMPLEX PETSC_DEPRECATED_MACRO("GCC warning \"MPIU_C_COMPLEX macro is deprecated use MPI_C_COMPLEX (since version 3.15)\"")
415 #define MPIU_C_DOUBLE_COMPLEX MPI_C_DOUBLE_COMPLEX PETSC_DEPRECATED_MACRO("GCC warning \"MPIU_C_DOUBLE_COMPLEX macro is deprecated use MPI_C_DOUBLE_COMPLEX (since version 3.15)\"")
416 
417 #if defined(PETSC_USE_REAL___FLOAT128)
418 PETSC_EXTERN MPI_Datatype MPIU___COMPLEX128 PetscAttrMPITypeTag(__complex128);
419 #endif /* PETSC_USE_REAL___FLOAT128 */
420 
421 /*MC
422    MPIU_COMPLEX - MPI datatype corresponding to PetscComplex
423 
424    Notes:
425    In MPI calls that require an MPI datatype that matches a PetscComplex or array of PetscComplex values, pass this value.
426 
427    Level: beginner
428 
429 .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`, `MPIU_INT`, `PETSC_i`
430 M*/
431 #if defined(PETSC_USE_REAL_SINGLE)
432 #  define MPIU_COMPLEX MPI_C_COMPLEX
433 #elif defined(PETSC_USE_REAL_DOUBLE)
434 #  define MPIU_COMPLEX MPI_C_DOUBLE_COMPLEX
435 #elif defined(PETSC_USE_REAL___FLOAT128)
436 #  define MPIU_COMPLEX MPIU___COMPLEX128
437 #elif defined(PETSC_USE_REAL___FP16)
438 #  define MPIU_COMPLEX MPI_C_COMPLEX
439 #endif /* PETSC_USE_REAL_* */
440 
441 #endif /* PETSC_HAVE_COMPLEX */
442 
443 /*
444     Scalar number definitions
445  */
446 #if defined(PETSC_USE_COMPLEX) && defined(PETSC_HAVE_COMPLEX)
447 /*MC
448    MPIU_SCALAR - MPI datatype corresponding to PetscScalar
449 
450    Notes:
451    In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalar values, pass this value.
452 
453    Level: beginner
454 
455 .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_REAL`, `MPIU_COMPLEX`, `MPIU_INT`
456 M*/
457 #define MPIU_SCALAR MPIU_COMPLEX
458 
459 /*MC
460    PetscRealPart - Returns the real part of a PetscScalar
461 
462    Synopsis:
463    #include <petscmath.h>
464    PetscReal PetscRealPart(PetscScalar v)
465 
466    Not Collective
467 
468    Input Parameter:
469 .  v - value to find the real part of
470 
471    Level: beginner
472 
473 .seealso: `PetscScalar`, `PetscImaginaryPart()`, `PetscMax()`, `PetscClipInterval()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscSqr()`
474 
475 M*/
476 #define PetscRealPart(a)      PetscRealPartComplex(a)
477 
478 /*MC
479    PetscImaginaryPart - Returns the imaginary part of a PetscScalar
480 
481    Synopsis:
482    #include <petscmath.h>
483    PetscReal PetscImaginaryPart(PetscScalar v)
484 
485    Not Collective
486 
487    Input Parameter:
488 .  v - value to find the imaginary part of
489 
490    Level: beginner
491 
492    Notes:
493        If PETSc was configured for real numbers then this always returns the value 0
494 
495 .seealso: `PetscScalar`, `PetscRealPart()`, `PetscMax()`, `PetscClipInterval()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscSqr()`
496 
497 M*/
498 #define PetscImaginaryPart(a) PetscImaginaryPartComplex(a)
499 
500 #define PetscAbsScalar(a)     PetscAbsComplex(a)
501 #define PetscArgScalar(a)     PetscArgComplex(a)
502 #define PetscConj(a)          PetscConjComplex(a)
503 #define PetscSqrtScalar(a)    PetscSqrtComplex(a)
504 #define PetscPowScalar(a,b)   PetscPowComplex(a,b)
505 #define PetscExpScalar(a)     PetscExpComplex(a)
506 #define PetscLogScalar(a)     PetscLogComplex(a)
507 #define PetscSinScalar(a)     PetscSinComplex(a)
508 #define PetscCosScalar(a)     PetscCosComplex(a)
509 #define PetscTanScalar(a)     PetscTanComplex(a)
510 #define PetscAsinScalar(a)    PetscAsinComplex(a)
511 #define PetscAcosScalar(a)    PetscAcosComplex(a)
512 #define PetscAtanScalar(a)    PetscAtanComplex(a)
513 #define PetscSinhScalar(a)    PetscSinhComplex(a)
514 #define PetscCoshScalar(a)    PetscCoshComplex(a)
515 #define PetscTanhScalar(a)    PetscTanhComplex(a)
516 #define PetscAsinhScalar(a)   PetscAsinhComplex(a)
517 #define PetscAcoshScalar(a)   PetscAcoshComplex(a)
518 #define PetscAtanhScalar(a)   PetscAtanhComplex(a)
519 
520 #else /* PETSC_USE_COMPLEX */
521 #define MPIU_SCALAR MPIU_REAL
522 #define PetscRealPart(a)      (a)
523 #define PetscImaginaryPart(a) ((PetscReal)0)
524 #define PetscAbsScalar(a)     PetscAbsReal(a)
525 #define PetscArgScalar(a)     (((a) < (PetscReal)0) ? PETSC_PI : (PetscReal)0)
526 #define PetscConj(a)          (a)
527 #define PetscSqrtScalar(a)    PetscSqrtReal(a)
528 #define PetscPowScalar(a,b)   PetscPowReal(a,b)
529 #define PetscExpScalar(a)     PetscExpReal(a)
530 #define PetscLogScalar(a)     PetscLogReal(a)
531 #define PetscSinScalar(a)     PetscSinReal(a)
532 #define PetscCosScalar(a)     PetscCosReal(a)
533 #define PetscTanScalar(a)     PetscTanReal(a)
534 #define PetscAsinScalar(a)    PetscAsinReal(a)
535 #define PetscAcosScalar(a)    PetscAcosReal(a)
536 #define PetscAtanScalar(a)    PetscAtanReal(a)
537 #define PetscSinhScalar(a)    PetscSinhReal(a)
538 #define PetscCoshScalar(a)    PetscCoshReal(a)
539 #define PetscTanhScalar(a)    PetscTanhReal(a)
540 #define PetscAsinhScalar(a)   PetscAsinhReal(a)
541 #define PetscAcoshScalar(a)   PetscAcoshReal(a)
542 #define PetscAtanhScalar(a)   PetscAtanhReal(a)
543 
544 #endif /* PETSC_USE_COMPLEX */
545 
546 /*
547    Certain objects may be created using either single or double precision.
548    This is currently not used.
549 */
550 typedef enum { PETSC_SCALAR_DOUBLE, PETSC_SCALAR_SINGLE, PETSC_SCALAR_LONG_DOUBLE, PETSC_SCALAR_HALF } PetscScalarPrecision;
551 
552 /* --------------------------------------------------------------------------*/
553 
554 /*MC
555    PetscAbs - Returns the absolute value of a number
556 
557    Synopsis:
558    #include <petscmath.h>
559    type PetscAbs(type v)
560 
561    Not Collective
562 
563    Input Parameter:
564 .  v - the number
565 
566    Notes:
567     type can be integer or real floating point value
568 
569    Level: beginner
570 
571 .seealso: `PetscAbsInt()`, `PetscAbsReal()`, `PetscAbsScalar()`
572 
573 M*/
574 #define PetscAbs(a)  (((a) >= 0) ? (a) : (-(a)))
575 
576 /*MC
577    PetscSign - Returns the sign of a number as an integer
578 
579    Synopsis:
580    #include <petscmath.h>
581    int PetscSign(type v)
582 
583    Not Collective
584 
585    Input Parameter:
586 .  v - the number
587 
588    Notes:
589     type can be integer or real floating point value
590 
591    Level: beginner
592 
593 M*/
594 #define PetscSign(a) (((a) >= 0) ? ((a) == 0 ? 0 : 1) : -1)
595 
596 /*MC
597    PetscMin - Returns minimum of two numbers
598 
599    Synopsis:
600    #include <petscmath.h>
601    type PetscMin(type v1,type v2)
602 
603    Not Collective
604 
605    Input Parameters:
606 +  v1 - first value to find minimum of
607 -  v2 - second value to find minimum of
608 
609    Notes:
610     type can be integer or floating point value
611 
612    Level: beginner
613 
614 .seealso: `PetscMax()`, `PetscClipInterval()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscSqr()`
615 
616 M*/
617 #define PetscMin(a,b)   (((a)<(b)) ?  (a) : (b))
618 
619 /*MC
620    PetscMax - Returns maxium of two numbers
621 
622    Synopsis:
623    #include <petscmath.h>
624    type max PetscMax(type v1,type v2)
625 
626    Not Collective
627 
628    Input Parameters:
629 +  v1 - first value to find maximum of
630 -  v2 - second value to find maximum of
631 
632    Notes:
633     type can be integer or floating point value
634 
635    Level: beginner
636 
637 .seealso: `PetscMin()`, `PetscClipInterval()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscSqr()`
638 
639 M*/
640 #define PetscMax(a,b)   (((a)<(b)) ?  (b) : (a))
641 
642 /*MC
643    PetscClipInterval - Returns a number clipped to be within an interval
644 
645    Synopsis:
646    #include <petscmath.h>
647    type clip PetscClipInterval(type x,type a,type b)
648 
649    Not Collective
650 
651    Input Parameters:
652 +  x - value to use if within interval [a,b]
653 .  a - lower end of interval
654 -  b - upper end of interval
655 
656    Notes:
657     type can be integer or floating point value
658 
659    Level: beginner
660 
661 .seealso: `PetscMin()`, `PetscMax()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscSqr()`
662 
663 M*/
664 #define PetscClipInterval(x,a,b)   (PetscMax((a),PetscMin((x),(b))))
665 
666 /*MC
667    PetscAbsInt - Returns the absolute value of an integer
668 
669    Synopsis:
670    #include <petscmath.h>
671    int abs PetscAbsInt(int v1)
672 
673    Not Collective
674 
675    Input Parameter:
676 .   v1 - the integer
677 
678    Level: beginner
679 
680 .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsReal()`, `PetscSqr()`
681 
682 M*/
683 #define PetscAbsInt(a)  (((a)<0)   ? (-(a)) : (a))
684 
685 /*MC
686    PetscAbsReal - Returns the absolute value of an real number
687 
688    Synopsis:
689    #include <petscmath.h>
690    Real abs PetscAbsReal(PetscReal v1)
691 
692    Not Collective
693 
694    Input Parameter:
695 .   v1 - the double
696 
697    Level: beginner
698 
699 .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsInt()`, `PetscSqr()`
700 
701 M*/
702 #if defined(PETSC_USE_REAL_SINGLE)
703 #define PetscAbsReal(a) fabsf(a)
704 #elif defined(PETSC_USE_REAL_DOUBLE)
705 #define PetscAbsReal(a) fabs(a)
706 #elif defined(PETSC_USE_REAL___FLOAT128)
707 #define PetscAbsReal(a) fabsq(a)
708 #elif defined(PETSC_USE_REAL___FP16)
709 #define PetscAbsReal(a) fabsf(a)
710 #endif
711 
712 /*MC
713    PetscSqr - Returns the square of a number
714 
715    Synopsis:
716    #include <petscmath.h>
717    type sqr PetscSqr(type v1)
718 
719    Not Collective
720 
721    Input Parameter:
722 .   v1 - the value
723 
724    Notes:
725     type can be integer or floating point value
726 
727    Level: beginner
728 
729 .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsInt()`, `PetscAbsReal()`
730 
731 M*/
732 #define PetscSqr(a)     ((a)*(a))
733 
734 /* ----------------------------------------------------------------------------*/
735 
736 #if defined(PETSC_USE_REAL_SINGLE)
737 #define PetscRealConstant(constant) constant##F
738 #elif defined(PETSC_USE_REAL_DOUBLE)
739 #define PetscRealConstant(constant) constant
740 #elif defined(PETSC_USE_REAL___FLOAT128)
741 #define PetscRealConstant(constant) constant##Q
742 #elif defined(PETSC_USE_REAL___FP16)
743 #define PetscRealConstant(constant) constant##F
744 #endif
745 
746 /*
747      Basic constants
748 */
749 #define PETSC_PI    PetscRealConstant(3.1415926535897932384626433832795029)
750 #define PETSC_PHI   PetscRealConstant(1.6180339887498948482045868343656381)
751 #define PETSC_SQRT2 PetscRealConstant(1.4142135623730950488016887242096981)
752 
753 #if !defined(PETSC_USE_64BIT_INDICES)
754 #define PETSC_MAX_INT            2147483647
755 #define PETSC_MIN_INT            (-PETSC_MAX_INT - 1)
756 #else
757 #define PETSC_MAX_INT            9223372036854775807L
758 #define PETSC_MIN_INT            (-PETSC_MAX_INT - 1)
759 #endif
760 #define PETSC_MAX_UINT16         65535
761 
762 #if defined(PETSC_USE_REAL_SINGLE)
763 #  define PETSC_MAX_REAL                3.40282346638528860e+38F
764 #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
765 #  define PETSC_MACHINE_EPSILON         1.19209290e-07F
766 #  define PETSC_SQRT_MACHINE_EPSILON    3.45266983e-04F
767 #  define PETSC_SMALL                   1.e-5F
768 #elif defined(PETSC_USE_REAL_DOUBLE)
769 #  define PETSC_MAX_REAL                1.7976931348623157e+308
770 #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
771 #  define PETSC_MACHINE_EPSILON         2.2204460492503131e-16
772 #  define PETSC_SQRT_MACHINE_EPSILON    1.490116119384766e-08
773 #  define PETSC_SMALL                   1.e-10
774 #elif defined(PETSC_USE_REAL___FLOAT128)
775 #  define PETSC_MAX_REAL                FLT128_MAX
776 #  define PETSC_MIN_REAL                (-FLT128_MAX)
777 #  define PETSC_MACHINE_EPSILON         FLT128_EPSILON
778 #  define PETSC_SQRT_MACHINE_EPSILON    1.38777878078144567552953958511352539e-17Q
779 #  define PETSC_SMALL                   1.e-20Q
780 #elif defined(PETSC_USE_REAL___FP16)
781 #  define PETSC_MAX_REAL                65504.0F
782 #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
783 #  define PETSC_MACHINE_EPSILON         .0009765625F
784 #  define PETSC_SQRT_MACHINE_EPSILON    .03125F
785 #  define PETSC_SMALL                   5.e-3F
786 #endif
787 
788 #define PETSC_INFINITY               (PETSC_MAX_REAL/4)
789 #define PETSC_NINFINITY              (-PETSC_INFINITY)
790 
791 PETSC_EXTERN PetscBool PetscIsInfReal(PetscReal);
792 PETSC_EXTERN PetscBool PetscIsNanReal(PetscReal);
793 PETSC_EXTERN PetscBool PetscIsNormalReal(PetscReal);
794 static inline PetscBool PetscIsInfOrNanReal(PetscReal v) {return PetscIsInfReal(v) || PetscIsNanReal(v) ? PETSC_TRUE : PETSC_FALSE;}
795 static inline PetscBool PetscIsInfScalar(PetscScalar v) {return PetscIsInfReal(PetscAbsScalar(v));}
796 static inline PetscBool PetscIsNanScalar(PetscScalar v) {return PetscIsNanReal(PetscAbsScalar(v));}
797 static inline PetscBool PetscIsInfOrNanScalar(PetscScalar v) {return PetscIsInfOrNanReal(PetscAbsScalar(v));}
798 static inline PetscBool PetscIsNormalScalar(PetscScalar v) {return PetscIsNormalReal(PetscAbsScalar(v));}
799 
800 PETSC_EXTERN PetscBool PetscIsCloseAtTol(PetscReal,PetscReal,PetscReal,PetscReal);
801 PETSC_EXTERN PetscBool PetscEqualReal(PetscReal,PetscReal);
802 PETSC_EXTERN PetscBool PetscEqualScalar(PetscScalar,PetscScalar);
803 
804 /*
805     These macros are currently hardwired to match the regular data types, so there is no support for a different
806     MatScalar from PetscScalar. We left the MatScalar in the source just in case we use it again.
807  */
808 #define MPIU_MATSCALAR MPIU_SCALAR
809 typedef PetscScalar MatScalar;
810 typedef PetscReal MatReal;
811 
812 struct petsc_mpiu_2scalar {PetscScalar a,b;};
813 PETSC_EXTERN MPI_Datatype MPIU_2SCALAR PetscAttrMPITypeTagLayoutCompatible(struct petsc_mpiu_2scalar);
814 
815 /*
816    MPI Datatypes for composite reductions:
817    MPIU_REAL_INT -> struct { PetscReal; PetscInt; }
818    MPIU_SCALAR_INT -> struct { PetscScalar; PetscInt; }
819 */
820 PETSC_EXTERN MPI_Datatype MPIU_REAL_INT;
821 PETSC_EXTERN MPI_Datatype MPIU_SCALAR_INT;
822 
823 #if defined(PETSC_USE_64BIT_INDICES)
824 struct petsc_mpiu_2int {PetscInt a,b;};
825 PETSC_EXTERN MPI_Datatype MPIU_2INT PetscAttrMPITypeTagLayoutCompatible(struct petsc_mpiu_2int);
826 #else
827 #define MPIU_2INT MPI_2INT
828 #endif
829 PETSC_EXTERN MPI_Datatype MPI_4INT;
830 PETSC_EXTERN MPI_Datatype MPIU_4INT;
831 
832 static inline PetscInt PetscPowInt(PetscInt base,PetscInt power)
833 {
834   PetscInt result = 1;
835   while (power) {
836     if (power & 1) result *= base;
837     power >>= 1;
838     base *= base;
839   }
840   return result;
841 }
842 
843 static inline PetscInt64 PetscPowInt64(PetscInt base,PetscInt power)
844 {
845   PetscInt64 result = 1;
846   while (power) {
847     if (power & 1) result *= base;
848     power >>= 1;
849     base *= base;
850   }
851   return result;
852 }
853 
854 static inline PetscReal PetscPowRealInt(PetscReal base,PetscInt power)
855 {
856   PetscReal result = 1;
857   if (power < 0) {
858     power = -power;
859     base  = ((PetscReal)1)/base;
860   }
861   while (power) {
862     if (power & 1) result *= base;
863     power >>= 1;
864     base *= base;
865   }
866   return result;
867 }
868 
869 static inline PetscScalar PetscPowScalarInt(PetscScalar base,PetscInt power)
870 {
871   PetscScalar result = (PetscReal)1;
872   if (power < 0) {
873     power = -power;
874     base  = ((PetscReal)1)/base;
875   }
876   while (power) {
877     if (power & 1) result *= base;
878     power >>= 1;
879     base *= base;
880   }
881   return result;
882 }
883 
884 static inline PetscScalar PetscPowScalarReal(PetscScalar base,PetscReal power)
885 {
886   PetscScalar cpower = power;
887   return PetscPowScalar(base,cpower);
888 }
889 
890 /*MC
891     PetscApproximateLTE - Performs a less than or equal to on a given constant with a fudge for floating point numbers
892 
893    Synopsis:
894    #include <petscmath.h>
895    bool PetscApproximateLTE(PetscReal x,constant float)
896 
897    Not Collective
898 
899    Input Parameters:
900 +   x - the variable
901 -   b - the constant float it is checking if x is less than or equal to
902 
903    Notes:
904      The fudge factor is the value PETSC_SMALL
905 
906      The constant numerical value is automatically set to the appropriate precision of PETSc so can just be provided as, for example, 3.2
907 
908      This is used in several examples for setting initial conditions based on coordinate values that are computed with i*h that produces inexact
909      floating point results.
910 
911    Level: advanced
912 
913 .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscApproximateGTE()`
914 
915 M*/
916 #define PetscApproximateLTE(x,b)  ((x) <= (PetscRealConstant(b)+PETSC_SMALL))
917 
918 /*MC
919     PetscApproximateGTE - Performs a greater than or equal to on a given constant with a fudge for floating point numbers
920 
921    Synopsis:
922    #include <petscmath.h>
923    bool PetscApproximateGTE(PetscReal x,constant float)
924 
925    Not Collective
926 
927    Input Parameters:
928 +   x - the variable
929 -   b - the constant float it is checking if x is greater than or equal to
930 
931    Notes:
932      The fudge factor is the value PETSC_SMALL
933 
934      The constant numerical value is automatically set to the appropriate precision of PETSc so can just be provided as, for example, 3.2
935 
936      This is used in several examples for setting initial conditions based on coordinate values that are computed with i*h that produces inexact
937      floating point results.
938 
939    Level: advanced
940 
941 .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscApproximateLTE()`
942 
943 M*/
944 #define PetscApproximateGTE(x,b)  ((x) >= (PetscRealConstant(b)-PETSC_SMALL))
945 
946 /*MC
947     PetscCeilInt - Returns the ceiling of the quotation of two positive integers
948 
949    Synopsis:
950    #include <petscmath.h>
951    PetscInt PetscCeilInt(PetscInt x,PetscInt y)
952 
953    Not Collective
954 
955    Input Parameters:
956 +   x - the numerator
957 -   y - the denominator
958 
959    Level: advanced
960 
961 .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscApproximateLTE()`
962 
963 M*/
964 #define PetscCeilInt(x,y)  ((((PetscInt)(x))/((PetscInt)(y))) +  ((((PetscInt)(x)) % ((PetscInt)(y))) ? 1 : 0))
965 
966 #define PetscCeilInt64(x,y)  ((((PetscInt64)(x))/((PetscInt64)(y))) +  ((((PetscInt64)(x)) % ((PetscInt64)(y))) ? 1 : 0))
967 
968 PETSC_EXTERN PetscErrorCode PetscLinearRegression(PetscInt,const PetscReal[],const PetscReal[],PetscReal*,PetscReal*);
969 #endif
970