Actual source code: petscmath.h
1: /*
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.
7: This file is included by petscsys.h and should not be used directly.
9: */
11: #if !defined(PETSCMATH_H)
12: #define PETSCMATH_H
13: #include <math.h>
14: #include <petscsystypes.h>
16: /*
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.
22: */
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
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
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
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
155: #endif /* PETSC_USE_REAL_* */
157: PETSC_STATIC_INLINE PetscReal PetscSignReal(PetscReal a)
158: {
159: return (PetscReal)((a < (PetscReal)0) ? -1 : ((a > (PetscReal)0) ? 1 : 0));
160: }
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
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
177: /*MC
178: MPIU_REAL - MPI datatype corresponding to PetscReal
180: Notes:
181: In MPI calls that require an MPI datatype that matches a PetscReal or array of PetscReal values, pass this value.
183: Level: beginner
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_* */
197: /*
198: Complex number definitions
199: */
200: #if defined(PETSC_HAVE_COMPLEX)
201: #if defined(__cplusplus)
202: /* C++ support of complex number */
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)
226: /* TODO: Add configure tests
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
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
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
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
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
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
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
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
295: */
297: #else /* C99 support of complex number */
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)
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)
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)
368: #endif /* PETSC_USE_REAL_* */
369: #endif /* (__cplusplus) */
371: /*
372: PETSC_i is the imaginary number, i
373: */
374: PETSC_EXTERN PetscComplex PETSC_i;
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),
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."
394: So type punning should be portable. */
395: union { PetscComplex z; PetscReal f[2]; } uz;
397: uz.f[0] = x;
398: uz.f[1] = y;
399: return uz.z;
400: }
401: #endif
402: }
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)\"")
407: #if defined(PETSC_USE_REAL___FLOAT128)
408: PETSC_EXTERN MPI_Datatype MPIU___COMPLEX128 PetscAttrMPITypeTag(__complex128);
409: #endif /* PETSC_USE_REAL___FLOAT128 */
411: /*MC
412: MPIU_COMPLEX - MPI datatype corresponding to PetscComplex
414: Notes:
415: In MPI calls that require an MPI datatype that matches a PetscComplex or array of PetscComplex values, pass this value.
417: Level: beginner
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_* */
431: #endif /* PETSC_HAVE_COMPLEX */
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
440: Notes:
441: In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalar values, pass this value.
443: Level: beginner
445: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_COMPLEX, MPIU_INT
446: M*/
447: #define MPIU_SCALAR MPIU_COMPLEX
449: /*MC
450: PetscRealPart - Returns the real part of a PetscScalar
452: Synopsis:
453: #include <petscmath.h>
454: PetscReal PetscRealPart(PetscScalar v)
456: Not Collective
458: Input Parameter:
459: . v - value to find the real part of
461: Level: beginner
463: .seealso: PetscScalar, PetscImaginaryPart(), PetscMax(), PetscClipInterval(), PetscAbsInt(), PetscAbsReal(), PetscSqr()
465: M*/
466: #define PetscRealPart(a) PetscRealPartComplex(a)
468: /*MC
469: PetscImaginaryPart - Returns the imaginary part of a PetscScalar
471: Synopsis:
472: #include <petscmath.h>
473: PetscReal PetscImaginaryPart(PetscScalar v)
475: Not Collective
477: Input Parameter:
478: . v - value to find the imaginary part of
480: Level: beginner
482: Notes:
483: If PETSc was configured for real numbers then this always returns the value 0
485: .seealso: PetscScalar, PetscRealPart(), PetscMax(), PetscClipInterval(), PetscAbsInt(), PetscAbsReal(), PetscSqr()
487: M*/
488: #define PetscImaginaryPart(a) PetscImaginaryPartComplex(a)
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)
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)
534: #endif /* PETSC_USE_COMPLEX */
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;
542: /* --------------------------------------------------------------------------*/
544: /*MC
545: PetscAbs - Returns the absolute value of a number
547: Synopsis:
548: #include <petscmath.h>
549: type PetscAbs(type v)
551: Not Collective
553: Input Parameter:
554: . v - the number
556: Notes:
557: type can be integer or real floating point value
559: Level: beginner
561: .seealso: PetscAbsInt(), PetscAbsReal(), PetscAbsScalar()
563: M*/
564: #define PetscAbs(a) (((a) >= 0) ? (a) : (-(a)))
566: /*MC
567: PetscSign - Returns the sign of a number as an integer
569: Synopsis:
570: #include <petscmath.h>
571: int PetscSign(type v)
573: Not Collective
575: Input Parameter:
576: . v - the number
578: Notes:
579: type can be integer or real floating point value
581: Level: beginner
583: M*/
584: #define PetscSign(a) (((a) >= 0) ? ((a) == 0 ? 0 : 1) : -1)
586: /*MC
587: PetscMin - Returns minimum of two numbers
589: Synopsis:
590: #include <petscmath.h>
591: type PetscMin(type v1,type v2)
593: Not Collective
595: Input Parameter:
596: + v1 - first value to find minimum of
597: - v2 - second value to find minimum of
599: Notes:
600: type can be integer or floating point value
602: Level: beginner
604: .seealso: PetscMax(), PetscClipInterval(), PetscAbsInt(), PetscAbsReal(), PetscSqr()
606: M*/
607: #define PetscMin(a,b) (((a)<(b)) ? (a) : (b))
609: /*MC
610: PetscMax - Returns maxium of two numbers
612: Synopsis:
613: #include <petscmath.h>
614: type max PetscMax(type v1,type v2)
616: Not Collective
618: Input Parameter:
619: + v1 - first value to find maximum of
620: - v2 - second value to find maximum of
622: Notes:
623: type can be integer or floating point value
625: Level: beginner
627: .seealso: PetscMin(), PetscClipInterval(), PetscAbsInt(), PetscAbsReal(), PetscSqr()
629: M*/
630: #define PetscMax(a,b) (((a)<(b)) ? (b) : (a))
632: /*MC
633: PetscClipInterval - Returns a number clipped to be within an interval
635: Synopsis:
636: #include <petscmath.h>
637: type clip PetscClipInterval(type x,type a,type b)
639: Not Collective
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
646: Notes:
647: type can be integer or floating point value
649: Level: beginner
651: .seealso: PetscMin(), PetscMax(), PetscAbsInt(), PetscAbsReal(), PetscSqr()
653: M*/
654: #define PetscClipInterval(x,a,b) (PetscMax((a),PetscMin((x),(b))))
656: /*MC
657: PetscAbsInt - Returns the absolute value of an integer
659: Synopsis:
660: #include <petscmath.h>
661: int abs PetscAbsInt(int v1)
663: Not Collective
665: Input Parameter:
666: . v1 - the integer
668: Level: beginner
670: .seealso: PetscMax(), PetscMin(), PetscAbsReal(), PetscSqr()
672: M*/
673: #define PetscAbsInt(a) (((a)<0) ? (-(a)) : (a))
675: /*MC
676: PetscAbsReal - Returns the absolute value of an real number
678: Synopsis:
679: #include <petscmath.h>
680: Real abs PetscAbsReal(PetscReal v1)
682: Not Collective
684: Input Parameter:
685: . v1 - the double
688: Level: beginner
690: .seealso: PetscMax(), PetscMin(), PetscAbsInt(), PetscSqr()
692: M*/
693: #if defined(PETSC_USE_REAL_SINGLE)
694: #define PetscAbsReal(a) fabsf(a)
695: #elif defined(PETSC_USE_REAL_DOUBLE)
696: #define PetscAbsReal(a) fabs(a)
697: #elif defined(PETSC_USE_REAL___FLOAT128)
698: #define PetscAbsReal(a) fabsq(a)
699: #elif defined(PETSC_USE_REAL___FP16)
700: #define PetscAbsReal(a) fabsf(a)
701: #endif
703: /*MC
704: PetscSqr - Returns the square of a number
706: Synopsis:
707: #include <petscmath.h>
708: type sqr PetscSqr(type v1)
710: Not Collective
712: Input Parameter:
713: . v1 - the value
715: Notes:
716: type can be integer or floating point value
718: Level: beginner
720: .seealso: PetscMax(), PetscMin(), PetscAbsInt(), PetscAbsReal()
722: M*/
723: #define PetscSqr(a) ((a)*(a))
725: /* ----------------------------------------------------------------------------*/
727: #if defined(PETSC_USE_REAL_SINGLE)
728: #define PetscRealConstant(constant) constant##F
729: #elif defined(PETSC_USE_REAL_DOUBLE)
730: #define PetscRealConstant(constant) constant
731: #elif defined(PETSC_USE_REAL___FLOAT128)
732: #define PetscRealConstant(constant) constant##Q
733: #elif defined(PETSC_USE_REAL___FP16)
734: #define PetscRealConstant(constant) constant##F
735: #endif
737: /*
738: Basic constants
739: */
740: #define PETSC_PI PetscRealConstant(3.1415926535897932384626433832795029)
741: #define PETSC_PHI PetscRealConstant(1.6180339887498948482045868343656381)
742: #define PETSC_SQRT2 PetscRealConstant(1.4142135623730950488016887242096981)
744: #if !defined(PETSC_USE_64BIT_INDICES)
745: #define PETSC_MAX_INT 2147483647
746: #define PETSC_MIN_INT (-PETSC_MAX_INT - 1)
747: #else
748: #define PETSC_MAX_INT 9223372036854775807L
749: #define PETSC_MIN_INT (-PETSC_MAX_INT - 1)
750: #endif
751: #define PETSC_MAX_UINT16 65535
753: #if defined(PETSC_USE_REAL_SINGLE)
754: # define PETSC_MAX_REAL 3.40282346638528860e+38F
755: # define PETSC_MIN_REAL (-PETSC_MAX_REAL)
756: # define PETSC_MACHINE_EPSILON 1.19209290e-07F
757: # define PETSC_SQRT_MACHINE_EPSILON 3.45266983e-04F
758: # define PETSC_SMALL 1.e-5F
759: #elif defined(PETSC_USE_REAL_DOUBLE)
760: # define PETSC_MAX_REAL 1.7976931348623157e+308
761: # define PETSC_MIN_REAL (-PETSC_MAX_REAL)
762: # define PETSC_MACHINE_EPSILON 2.2204460492503131e-16
763: # define PETSC_SQRT_MACHINE_EPSILON 1.490116119384766e-08
764: # define PETSC_SMALL 1.e-10
765: #elif defined(PETSC_USE_REAL___FLOAT128)
766: # define PETSC_MAX_REAL FLT128_MAX
767: # define PETSC_MIN_REAL (-FLT128_MAX)
768: # define PETSC_MACHINE_EPSILON FLT128_EPSILON
769: # define PETSC_SQRT_MACHINE_EPSILON 1.38777878078144567552953958511352539e-17Q
770: # define PETSC_SMALL 1.e-20Q
771: #elif defined(PETSC_USE_REAL___FP16)
772: # define PETSC_MAX_REAL 65504.0F
773: # define PETSC_MIN_REAL (-PETSC_MAX_REAL)
774: # define PETSC_MACHINE_EPSILON .0009765625F
775: # define PETSC_SQRT_MACHINE_EPSILON .03125F
776: # define PETSC_SMALL 5.e-3F
777: #endif
779: #define PETSC_INFINITY (PETSC_MAX_REAL/4)
780: #define PETSC_NINFINITY (-PETSC_INFINITY)
782: PETSC_EXTERN PetscBool PetscIsInfReal(PetscReal);
783: PETSC_EXTERN PetscBool PetscIsNanReal(PetscReal);
784: PETSC_EXTERN PetscBool PetscIsNormalReal(PetscReal);
785: PETSC_STATIC_INLINE PetscBool PetscIsInfOrNanReal(PetscReal v) {return PetscIsInfReal(v) || PetscIsNanReal(v) ? PETSC_TRUE : PETSC_FALSE;}
786: PETSC_STATIC_INLINE PetscBool PetscIsInfScalar(PetscScalar v) {return PetscIsInfReal(PetscAbsScalar(v));}
787: PETSC_STATIC_INLINE PetscBool PetscIsNanScalar(PetscScalar v) {return PetscIsNanReal(PetscAbsScalar(v));}
788: PETSC_STATIC_INLINE PetscBool PetscIsInfOrNanScalar(PetscScalar v) {return PetscIsInfOrNanReal(PetscAbsScalar(v));}
789: PETSC_STATIC_INLINE PetscBool PetscIsNormalScalar(PetscScalar v) {return PetscIsNormalReal(PetscAbsScalar(v));}
791: PETSC_EXTERN PetscBool PetscIsCloseAtTol(PetscReal,PetscReal,PetscReal,PetscReal);
792: PETSC_EXTERN PetscBool PetscEqualReal(PetscReal,PetscReal);
793: PETSC_EXTERN PetscBool PetscEqualScalar(PetscScalar,PetscScalar);
795: /*
796: These macros are currently hardwired to match the regular data types, so there is no support for a different
797: MatScalar from PetscScalar. We left the MatScalar in the source just in case we use it again.
798: */
799: #define MPIU_MATSCALAR MPIU_SCALAR
800: typedef PetscScalar MatScalar;
801: typedef PetscReal MatReal;
803: struct petsc_mpiu_2scalar {PetscScalar a,b;};
804: PETSC_EXTERN MPI_Datatype MPIU_2SCALAR PetscAttrMPITypeTagLayoutCompatible(struct petsc_mpiu_2scalar);
806: #if defined(PETSC_USE_64BIT_INDICES)
807: struct petsc_mpiu_2int {PetscInt a,b;};
808: PETSC_EXTERN MPI_Datatype MPIU_2INT PetscAttrMPITypeTagLayoutCompatible(struct petsc_mpiu_2int);
809: #else
810: #define MPIU_2INT MPI_2INT
811: #endif
813: PETSC_STATIC_INLINE PetscInt PetscPowInt(PetscInt base,PetscInt power)
814: {
815: PetscInt result = 1;
816: while (power) {
817: if (power & 1) result *= base;
818: power >>= 1;
819: base *= base;
820: }
821: return result;
822: }
824: PETSC_STATIC_INLINE PetscInt64 PetscPowInt64(PetscInt base,PetscInt power)
825: {
826: PetscInt64 result = 1;
827: while (power) {
828: if (power & 1) result *= base;
829: power >>= 1;
830: base *= base;
831: }
832: return result;
833: }
835: PETSC_STATIC_INLINE PetscReal PetscPowRealInt(PetscReal base,PetscInt power)
836: {
837: PetscReal result = 1;
838: if (power < 0) {
839: power = -power;
840: base = ((PetscReal)1)/base;
841: }
842: while (power) {
843: if (power & 1) result *= base;
844: power >>= 1;
845: base *= base;
846: }
847: return result;
848: }
850: PETSC_STATIC_INLINE PetscScalar PetscPowScalarInt(PetscScalar base,PetscInt power)
851: {
852: PetscScalar result = (PetscReal)1;
853: if (power < 0) {
854: power = -power;
855: base = ((PetscReal)1)/base;
856: }
857: while (power) {
858: if (power & 1) result *= base;
859: power >>= 1;
860: base *= base;
861: }
862: return result;
863: }
865: PETSC_STATIC_INLINE PetscScalar PetscPowScalarReal(PetscScalar base,PetscReal power)
866: {
867: PetscScalar cpower = power;
868: return PetscPowScalar(base,cpower);
869: }
871: /*MC
872: PetscLTE - Performs a less than or equal to on a given constant with a fudge for floating point numbers
874: Synopsis:
875: #include <petscmath.h>
876: bool PetscLTE(PetscReal x,constant float)
878: Not Collective
880: Input Parameters:
881: + x - the variable
882: - b - the constant float it is checking if x is less than or equal to
884: Notes:
885: The fudge factor is the value PETSC_SMALL
887: The constant numerical value is automatically set to the appropriate precision of PETSc so can just be provided as, for example, 3.2
889: This is used in several examples for setting initial conditions based on coordinate values that are computed with i*h that produces inexact
890: floating point results.
892: Level: advanced
894: .seealso: PetscMax(), PetscMin(), PetscAbsInt(), PetscAbsReal(), PetscGTE()
896: M*/
897: #define PetscLTE(x,b) ((x) <= (PetscRealConstant(b)+PETSC_SMALL))
899: /*MC
900: PetscGTE - Performs a greater than or equal to on a given constant with a fudge for floating point numbers
902: Synopsis:
903: #include <petscmath.h>
904: bool PetscGTE(PetscReal x,constant float)
906: Not Collective
908: Input Parameters:
909: + x - the variable
910: - b - the constant float it is checking if x is greater than or equal to
912: Notes:
913: The fudge factor is the value PETSC_SMALL
915: The constant numerical value is automatically set to the appropriate precision of PETSc so can just be provided as, for example, 3.2
917: This is used in several examples for setting initial conditions based on coordinate values that are computed with i*h that produces inexact
918: floating point results.
920: Level: advanced
922: .seealso: PetscMax(), PetscMin(), PetscAbsInt(), PetscAbsReal(), PetscLTE()
924: M*/
925: #define PetscGTE(x,b) ((x) >= (PetscRealConstant(b)-PETSC_SMALL))
928: PETSC_EXTERN PetscErrorCode PetscLinearRegression(PetscInt,const PetscReal[],const PetscReal[],PetscReal*,PetscReal*);
929: #endif