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