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 PetscErfReal(a)        erff(a)
 50: #define PetscCeilReal(a)       ceilf(a)
 51: #define PetscFloorReal(a)      floorf(a)
 52: #define PetscFmodReal(a,b)     fmodf(a,b)
 53: #define PetscCopysignReal(a,b) copysignf(a,b)
 54: #define PetscTGamma(a)         tgammaf(a)
 55: #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
 56: #define PetscLGamma(a)         gammaf(a)
 57: #else
 58: #define PetscLGamma(a)         lgammaf(a)
 59: #endif

 61: #elif defined(PETSC_USE_REAL_DOUBLE)
 62: #define PetscSqrtReal(a)       sqrt(a)
 63: #define PetscCbrtReal(a)       cbrt(a)
 64: #define PetscHypotReal(a,b)    hypot(a,b)
 65: #define PetscAtan2Real(a,b)    atan2(a,b)
 66: #define PetscPowReal(a,b)      pow(a,b)
 67: #define PetscExpReal(a)        exp(a)
 68: #define PetscLogReal(a)        log(a)
 69: #define PetscLog10Real(a)      log10(a)
 70: #define PetscLog2Real(a)       log2(a)
 71: #define PetscSinReal(a)        sin(a)
 72: #define PetscCosReal(a)        cos(a)
 73: #define PetscTanReal(a)        tan(a)
 74: #define PetscAsinReal(a)       asin(a)
 75: #define PetscAcosReal(a)       acos(a)
 76: #define PetscAtanReal(a)       atan(a)
 77: #define PetscSinhReal(a)       sinh(a)
 78: #define PetscCoshReal(a)       cosh(a)
 79: #define PetscTanhReal(a)       tanh(a)
 80: #define PetscAsinhReal(a)      asinh(a)
 81: #define PetscAcoshReal(a)      acosh(a)
 82: #define PetscAtanhReal(a)      atanh(a)
 83: #define PetscErfReal(a)        erf(a)
 84: #define PetscCeilReal(a)       ceil(a)
 85: #define PetscFloorReal(a)      floor(a)
 86: #define PetscFmodReal(a,b)     fmod(a,b)
 87: #define PetscCopysignReal(a,b) copysign(a,b)
 88: #define PetscTGamma(a)         tgamma(a)
 89: #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
 90: #define PetscLGamma(a)         gamma(a)
 91: #else
 92: #define PetscLGamma(a)         lgamma(a)
 93: #endif

 95: #elif defined(PETSC_USE_REAL___FLOAT128)
 96: #define PetscSqrtReal(a)       sqrtq(a)
 97: #define PetscCbrtReal(a)       cbrtq(a)
 98: #define PetscHypotReal(a,b)    hypotq(a,b)
 99: #define PetscAtan2Real(a,b)    atan2q(a,b)
100: #define PetscPowReal(a,b)      powq(a,b)
101: #define PetscExpReal(a)        expq(a)
102: #define PetscLogReal(a)        logq(a)
103: #define PetscLog10Real(a)      log10q(a)
104: #define PetscLog2Real(a)       log2q(a)
105: #define PetscSinReal(a)        sinq(a)
106: #define PetscCosReal(a)        cosq(a)
107: #define PetscTanReal(a)        tanq(a)
108: #define PetscAsinReal(a)       asinq(a)
109: #define PetscAcosReal(a)       acosq(a)
110: #define PetscAtanReal(a)       atanq(a)
111: #define PetscSinhReal(a)       sinhq(a)
112: #define PetscCoshReal(a)       coshq(a)
113: #define PetscTanhReal(a)       tanhq(a)
114: #define PetscAsinhReal(a)      asinhq(a)
115: #define PetscAcoshReal(a)      acoshq(a)
116: #define PetscAtanhReal(a)      atanhq(a)
117: #define PetscErfReal(a)        erfq(a)
118: #define PetscCeilReal(a)       ceilq(a)
119: #define PetscFloorReal(a)      floorq(a)
120: #define PetscFmodReal(a,b)     fmodq(a,b)
121: #define PetscCopysignReal(a,b) copysignq(a,b)
122: #define PetscTGamma(a)         tgammaq(a)
123: #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
124: #define PetscLGamma(a)         gammaq(a)
125: #else
126: #define PetscLGamma(a)         lgammaq(a)
127: #endif

129: #elif defined(PETSC_USE_REAL___FP16)
130: #define PetscSqrtReal(a)       sqrtf(a)
131: #define PetscCbrtReal(a)       cbrtf(a)
132: #define PetscHypotReal(a,b)    hypotf(a,b)
133: #define PetscAtan2Real(a,b)    atan2f(a,b)
134: #define PetscPowReal(a,b)      powf(a,b)
135: #define PetscExpReal(a)        expf(a)
136: #define PetscLogReal(a)        logf(a)
137: #define PetscLog10Real(a)      log10f(a)
138: #define PetscLog2Real(a)       log2f(a)
139: #define PetscSinReal(a)        sinf(a)
140: #define PetscCosReal(a)        cosf(a)
141: #define PetscTanReal(a)        tanf(a)
142: #define PetscAsinReal(a)       asinf(a)
143: #define PetscAcosReal(a)       acosf(a)
144: #define PetscAtanReal(a)       atanf(a)
145: #define PetscSinhReal(a)       sinhf(a)
146: #define PetscCoshReal(a)       coshf(a)
147: #define PetscTanhReal(a)       tanhf(a)
148: #define PetscAsinhReal(a)      asinhf(a)
149: #define PetscAcoshReal(a)      acoshf(a)
150: #define PetscAtanhReal(a)      atanhf(a)
151: #define PetscErfReal(a)        erff(a)
152: #define PetscCeilReal(a)       ceilf(a)
153: #define PetscFloorReal(a)      floorf(a)
154: #define PetscFmodReal(a,b)     fmodf(a,b)
155: #define PetscCopySignReal(a,b) copysignf(a,b)
156: #define PetscTGamma(a)         tgammaf(a)
157: #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
158: #define PetscLGamma(a)         gammaf(a)
159: #else
160: #define PetscLGamma(a)         lgammaf(a)
161: #endif

163: #endif /* PETSC_USE_REAL_* */

165: static inline PetscReal PetscSignReal(PetscReal a)
166: {
167:   return (PetscReal)((a < (PetscReal)0) ? -1 : ((a > (PetscReal)0) ? 1 : 0));
168: }

170: #if !defined(PETSC_HAVE_LOG2)
171: #undef PetscLog2Real
172: static inline PetscReal PetscLog2Real(PetscReal a)
173: {
174:   return PetscLogReal(a)/PetscLogReal((PetscReal)2);
175: }
176: #endif

178: #if defined(PETSC_USE_REAL___FLOAT128)
179: PETSC_EXTERN MPI_Datatype MPIU___FLOAT128 PetscAttrMPITypeTag(__float128);
180: #endif
181: #if defined(PETSC_USE_REAL___FP16)
182: PETSC_EXTERN MPI_Datatype MPIU___FP16 PetscAttrMPITypeTag(__fp16);
183: #endif

185: /*MC
186:    MPIU_REAL - MPI datatype corresponding to PetscReal

188:    Notes:
189:    In MPI calls that require an MPI datatype that matches a PetscReal or array of PetscReal values, pass this value.

191:    Level: beginner

193: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_SCALAR, MPIU_COMPLEX, MPIU_INT
194: M*/
195: #if defined(PETSC_USE_REAL_SINGLE)
196: #  define MPIU_REAL MPI_FLOAT
197: #elif defined(PETSC_USE_REAL_DOUBLE)
198: #  define MPIU_REAL MPI_DOUBLE
199: #elif defined(PETSC_USE_REAL___FLOAT128)
200: #  define MPIU_REAL MPIU___FLOAT128
201: #elif defined(PETSC_USE_REAL___FP16)
202: #  define MPIU_REAL MPIU___FP16
203: #endif /* PETSC_USE_REAL_* */

205: /*
206:     Complex number definitions
207:  */
208: #if defined(PETSC_HAVE_COMPLEX)
209: #if defined(__cplusplus) && !defined(PETSC_USE_REAL___FLOAT128)
210: /* C++ support of complex number */

212: #define PetscRealPartComplex(a)      (a).real()
213: #define PetscImaginaryPartComplex(a) (a).imag()
214: #define PetscAbsComplex(a)           petsccomplexlib::abs(a)
215: #define PetscArgComplex(a)           petsccomplexlib::arg(a)
216: #define PetscConjComplex(a)          petsccomplexlib::conj(a)
217: #define PetscSqrtComplex(a)          petsccomplexlib::sqrt(a)
218: #define PetscPowComplex(a,b)         petsccomplexlib::pow(a,b)
219: #define PetscExpComplex(a)           petsccomplexlib::exp(a)
220: #define PetscLogComplex(a)           petsccomplexlib::log(a)
221: #define PetscSinComplex(a)           petsccomplexlib::sin(a)
222: #define PetscCosComplex(a)           petsccomplexlib::cos(a)
223: #define PetscTanComplex(a)           petsccomplexlib::tan(a)
224: #define PetscAsinComplex(a)          petsccomplexlib::asin(a)
225: #define PetscAcosComplex(a)          petsccomplexlib::acos(a)
226: #define PetscAtanComplex(a)          petsccomplexlib::atan(a)
227: #define PetscSinhComplex(a)          petsccomplexlib::sinh(a)
228: #define PetscCoshComplex(a)          petsccomplexlib::cosh(a)
229: #define PetscTanhComplex(a)          petsccomplexlib::tanh(a)
230: #define PetscAsinhComplex(a)         petsccomplexlib::asinh(a)
231: #define PetscAcoshComplex(a)         petsccomplexlib::acosh(a)
232: #define PetscAtanhComplex(a)         petsccomplexlib::atanh(a)

234: /* TODO: Add configure tests

236: #if !defined(PETSC_HAVE_CXX_TAN_COMPLEX)
237: #undef PetscTanComplex
238: static inline PetscComplex PetscTanComplex(PetscComplex z)
239: {
240:   return PetscSinComplex(z)/PetscCosComplex(z);
241: }
242: #endif

244: #if !defined(PETSC_HAVE_CXX_TANH_COMPLEX)
245: #undef PetscTanhComplex
246: static inline PetscComplex PetscTanhComplex(PetscComplex z)
247: {
248:   return PetscSinhComplex(z)/PetscCoshComplex(z);
249: }
250: #endif

252: #if !defined(PETSC_HAVE_CXX_ASIN_COMPLEX)
253: #undef PetscAsinComplex
254: static inline PetscComplex PetscAsinComplex(PetscComplex z)
255: {
256:   const PetscComplex j(0,1);
257:   return -j*PetscLogComplex(j*z+PetscSqrtComplex(1.0f-z*z));
258: }
259: #endif

261: #if !defined(PETSC_HAVE_CXX_ACOS_COMPLEX)
262: #undef PetscAcosComplex
263: static inline PetscComplex PetscAcosComplex(PetscComplex z)
264: {
265:   const PetscComplex j(0,1);
266:   return j*PetscLogComplex(z-j*PetscSqrtComplex(1.0f-z*z));
267: }
268: #endif

270: #if !defined(PETSC_HAVE_CXX_ATAN_COMPLEX)
271: #undef PetscAtanComplex
272: static inline PetscComplex PetscAtanComplex(PetscComplex z)
273: {
274:   const PetscComplex j(0,1);
275:   return 0.5f*j*PetscLogComplex((1.0f-j*z)/(1.0f+j*z));
276: }
277: #endif

279: #if !defined(PETSC_HAVE_CXX_ASINH_COMPLEX)
280: #undef PetscAsinhComplex
281: static inline PetscComplex PetscAsinhComplex(PetscComplex z)
282: {
283:   return PetscLogComplex(z+PetscSqrtComplex(z*z+1.0f));
284: }
285: #endif

287: #if !defined(PETSC_HAVE_CXX_ACOSH_COMPLEX)
288: #undef PetscAcoshComplex
289: static inline PetscComplex PetscAcoshComplex(PetscComplex z)
290: {
291:   return PetscLogComplex(z+PetscSqrtComplex(z*z-1.0f));
292: }
293: #endif

295: #if !defined(PETSC_HAVE_CXX_ATANH_COMPLEX)
296: #undef PetscAtanhComplex
297: static inline PetscComplex PetscAtanhComplex(PetscComplex z)
298: {
299:   return 0.5f*PetscLogComplex((1.0f+z)/(1.0f-z));
300: }
301: #endif

303: */

305: #else /* C99 support of complex number */

307: #if defined(PETSC_USE_REAL_SINGLE)
308: #define PetscRealPartComplex(a)      crealf(a)
309: #define PetscImaginaryPartComplex(a) cimagf(a)
310: #define PetscAbsComplex(a)           cabsf(a)
311: #define PetscArgComplex(a)           cargf(a)
312: #define PetscConjComplex(a)          conjf(a)
313: #define PetscSqrtComplex(a)          csqrtf(a)
314: #define PetscPowComplex(a,b)         cpowf(a,b)
315: #define PetscExpComplex(a)           cexpf(a)
316: #define PetscLogComplex(a)           clogf(a)
317: #define PetscSinComplex(a)           csinf(a)
318: #define PetscCosComplex(a)           ccosf(a)
319: #define PetscTanComplex(a)           ctanf(a)
320: #define PetscAsinComplex(a)          casinf(a)
321: #define PetscAcosComplex(a)          cacosf(a)
322: #define PetscAtanComplex(a)          catanf(a)
323: #define PetscSinhComplex(a)          csinhf(a)
324: #define PetscCoshComplex(a)          ccoshf(a)
325: #define PetscTanhComplex(a)          ctanhf(a)
326: #define PetscAsinhComplex(a)         casinhf(a)
327: #define PetscAcoshComplex(a)         cacoshf(a)
328: #define PetscAtanhComplex(a)         catanhf(a)

330: #elif defined(PETSC_USE_REAL_DOUBLE)
331: #define PetscRealPartComplex(a)      creal(a)
332: #define PetscImaginaryPartComplex(a) cimag(a)
333: #define PetscAbsComplex(a)           cabs(a)
334: #define PetscArgComplex(a)           carg(a)
335: #define PetscConjComplex(a)          conj(a)
336: #define PetscSqrtComplex(a)          csqrt(a)
337: #define PetscPowComplex(a,b)         cpow(a,b)
338: #define PetscExpComplex(a)           cexp(a)
339: #define PetscLogComplex(a)           clog(a)
340: #define PetscSinComplex(a)           csin(a)
341: #define PetscCosComplex(a)           ccos(a)
342: #define PetscTanComplex(a)           ctan(a)
343: #define PetscAsinComplex(a)          casin(a)
344: #define PetscAcosComplex(a)          cacos(a)
345: #define PetscAtanComplex(a)          catan(a)
346: #define PetscSinhComplex(a)          csinh(a)
347: #define PetscCoshComplex(a)          ccosh(a)
348: #define PetscTanhComplex(a)          ctanh(a)
349: #define PetscAsinhComplex(a)         casinh(a)
350: #define PetscAcoshComplex(a)         cacosh(a)
351: #define PetscAtanhComplex(a)         catanh(a)

353: #elif defined(PETSC_USE_REAL___FLOAT128)
354: #define PetscRealPartComplex(a)      crealq(a)
355: #define PetscImaginaryPartComplex(a) cimagq(a)
356: #define PetscAbsComplex(a)           cabsq(a)
357: #define PetscArgComplex(a)           cargq(a)
358: #define PetscConjComplex(a)          conjq(a)
359: #define PetscSqrtComplex(a)          csqrtq(a)
360: #define PetscPowComplex(a,b)         cpowq(a,b)
361: #define PetscExpComplex(a)           cexpq(a)
362: #define PetscLogComplex(a)           clogq(a)
363: #define PetscSinComplex(a)           csinq(a)
364: #define PetscCosComplex(a)           ccosq(a)
365: #define PetscTanComplex(a)           ctanq(a)
366: #define PetscAsinComplex(a)          casinq(a)
367: #define PetscAcosComplex(a)          cacosq(a)
368: #define PetscAtanComplex(a)          catanq(a)
369: #define PetscSinhComplex(a)          csinhq(a)
370: #define PetscCoshComplex(a)          ccoshq(a)
371: #define PetscTanhComplex(a)          ctanhq(a)
372: #define PetscAsinhComplex(a)         casinhq(a)
373: #define PetscAcoshComplex(a)         cacoshq(a)
374: #define PetscAtanhComplex(a)         catanhq(a)

376: #endif /* PETSC_USE_REAL_* */
377: #endif /* (__cplusplus) */

379: /*
380:    PETSC_i is the imaginary number, i
381: */
382: PETSC_EXTERN PetscComplex PETSC_i;

384: /*
385:    Try to do the right thing for complex number construction: see
386:    http://www.open-std.org/jtc1/sc22/wg14/www/docs/n1464.htm
387:    for details
388: */
389: static inline PetscComplex PetscCMPLX(PetscReal x, PetscReal y)
390: {
391: #if defined(__cplusplus) && !defined(PETSC_USE_REAL___FLOAT128)
392:   return PetscComplex(x,y);
393: #elif defined(_Imaginary_I)
394:   return x + y * _Imaginary_I;
395: #else
396:   { /* In both C99 and C11 (ISO/IEC 9899, Section 6.2.5),

398:        "For each floating type there is a corresponding real type, which is always a real floating
399:        type. For real floating types, it is the same type. For complex types, it is the type given
400:        by deleting the keyword _Complex from the type name."

402:        So type punning should be portable. */
403:     union { PetscComplex z; PetscReal f[2]; } uz;

405:     uz.f[0] = x;
406:     uz.f[1] = y;
407:     return uz.z;
408:   }
409: #endif
410: }

412: #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)\"")
413: #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)\"")

415: #if defined(PETSC_USE_REAL___FLOAT128)
416: PETSC_EXTERN MPI_Datatype MPIU___COMPLEX128 PetscAttrMPITypeTag(__complex128);
417: #endif /* PETSC_USE_REAL___FLOAT128 */

419: /*MC
420:    MPIU_COMPLEX - MPI datatype corresponding to PetscComplex

422:    Notes:
423:    In MPI calls that require an MPI datatype that matches a PetscComplex or array of PetscComplex values, pass this value.

425:    Level: beginner

427: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_SCALAR, MPIU_COMPLEX, MPIU_INT, PETSC_i
428: M*/
429: #if defined(PETSC_USE_REAL_SINGLE)
430: #  define MPIU_COMPLEX MPI_C_COMPLEX
431: #elif defined(PETSC_USE_REAL_DOUBLE)
432: #  define MPIU_COMPLEX MPI_C_DOUBLE_COMPLEX
433: #elif defined(PETSC_USE_REAL___FLOAT128)
434: #  define MPIU_COMPLEX MPIU___COMPLEX128
435: #elif defined(PETSC_USE_REAL___FP16)
436: #  define MPIU_COMPLEX MPI_C_COMPLEX
437: #endif /* PETSC_USE_REAL_* */

439: #endif /* PETSC_HAVE_COMPLEX */

441: /*
442:     Scalar number definitions
443:  */
444: #if defined(PETSC_USE_COMPLEX) && defined(PETSC_HAVE_COMPLEX)
445: /*MC
446:    MPIU_SCALAR - MPI datatype corresponding to PetscScalar

448:    Notes:
449:    In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalar values, pass this value.

451:    Level: beginner

453: .seealso: PetscReal, PetscScalar, PetscComplex, PetscInt, MPIU_REAL, MPIU_COMPLEX, MPIU_INT
454: M*/
455: #define MPIU_SCALAR MPIU_COMPLEX

457: /*MC
458:    PetscRealPart - Returns the real part of a PetscScalar

460:    Synopsis:
461: #include <petscmath.h>
462:    PetscReal PetscRealPart(PetscScalar v)

464:    Not Collective

466:    Input Parameter:
467: .  v - value to find the real part of

469:    Level: beginner

471: .seealso: PetscScalar, PetscImaginaryPart(), PetscMax(), PetscClipInterval(), PetscAbsInt(), PetscAbsReal(), PetscSqr()

473: M*/
474: #define PetscRealPart(a)      PetscRealPartComplex(a)

476: /*MC
477:    PetscImaginaryPart - Returns the imaginary part of a PetscScalar

479:    Synopsis:
480: #include <petscmath.h>
481:    PetscReal PetscImaginaryPart(PetscScalar v)

483:    Not Collective

485:    Input Parameter:
486: .  v - value to find the imaginary part of

488:    Level: beginner

490:    Notes:
491:        If PETSc was configured for real numbers then this always returns the value 0

493: .seealso: PetscScalar, PetscRealPart(), PetscMax(), PetscClipInterval(), PetscAbsInt(), PetscAbsReal(), PetscSqr()

495: M*/
496: #define PetscImaginaryPart(a) PetscImaginaryPartComplex(a)

498: #define PetscAbsScalar(a)     PetscAbsComplex(a)
499: #define PetscArgScalar(a)     PetscArgComplex(a)
500: #define PetscConj(a)          PetscConjComplex(a)
501: #define PetscSqrtScalar(a)    PetscSqrtComplex(a)
502: #define PetscPowScalar(a,b)   PetscPowComplex(a,b)
503: #define PetscExpScalar(a)     PetscExpComplex(a)
504: #define PetscLogScalar(a)     PetscLogComplex(a)
505: #define PetscSinScalar(a)     PetscSinComplex(a)
506: #define PetscCosScalar(a)     PetscCosComplex(a)
507: #define PetscTanScalar(a)     PetscTanComplex(a)
508: #define PetscAsinScalar(a)    PetscAsinComplex(a)
509: #define PetscAcosScalar(a)    PetscAcosComplex(a)
510: #define PetscAtanScalar(a)    PetscAtanComplex(a)
511: #define PetscSinhScalar(a)    PetscSinhComplex(a)
512: #define PetscCoshScalar(a)    PetscCoshComplex(a)
513: #define PetscTanhScalar(a)    PetscTanhComplex(a)
514: #define PetscAsinhScalar(a)   PetscAsinhComplex(a)
515: #define PetscAcoshScalar(a)   PetscAcoshComplex(a)
516: #define PetscAtanhScalar(a)   PetscAtanhComplex(a)

518: #else /* PETSC_USE_COMPLEX */
519: #define MPIU_SCALAR MPIU_REAL
520: #define PetscRealPart(a)      (a)
521: #define PetscImaginaryPart(a) ((PetscReal)0)
522: #define PetscAbsScalar(a)     PetscAbsReal(a)
523: #define PetscArgScalar(a)     (((a) < (PetscReal)0) ? PETSC_PI : (PetscReal)0)
524: #define PetscConj(a)          (a)
525: #define PetscSqrtScalar(a)    PetscSqrtReal(a)
526: #define PetscPowScalar(a,b)   PetscPowReal(a,b)
527: #define PetscExpScalar(a)     PetscExpReal(a)
528: #define PetscLogScalar(a)     PetscLogReal(a)
529: #define PetscSinScalar(a)     PetscSinReal(a)
530: #define PetscCosScalar(a)     PetscCosReal(a)
531: #define PetscTanScalar(a)     PetscTanReal(a)
532: #define PetscAsinScalar(a)    PetscAsinReal(a)
533: #define PetscAcosScalar(a)    PetscAcosReal(a)
534: #define PetscAtanScalar(a)    PetscAtanReal(a)
535: #define PetscSinhScalar(a)    PetscSinhReal(a)
536: #define PetscCoshScalar(a)    PetscCoshReal(a)
537: #define PetscTanhScalar(a)    PetscTanhReal(a)
538: #define PetscAsinhScalar(a)   PetscAsinhReal(a)
539: #define PetscAcoshScalar(a)   PetscAcoshReal(a)
540: #define PetscAtanhScalar(a)   PetscAtanhReal(a)

542: #endif /* PETSC_USE_COMPLEX */

544: /*
545:    Certain objects may be created using either single or double precision.
546:    This is currently not used.
547: */
548: typedef enum { PETSC_SCALAR_DOUBLE, PETSC_SCALAR_SINGLE, PETSC_SCALAR_LONG_DOUBLE, PETSC_SCALAR_HALF } PetscScalarPrecision;

550: /* --------------------------------------------------------------------------*/

552: /*MC
553:    PetscAbs - Returns the absolute value of a number

555:    Synopsis:
556: #include <petscmath.h>
557:    type PetscAbs(type v)

559:    Not Collective

561:    Input Parameter:
562: .  v - the number

564:    Notes:
565:     type can be integer or real floating point value

567:    Level: beginner

569: .seealso: PetscAbsInt(), PetscAbsReal(), PetscAbsScalar()

571: M*/
572: #define PetscAbs(a)  (((a) >= 0) ? (a) : (-(a)))

574: /*MC
575:    PetscSign - Returns the sign of a number as an integer

577:    Synopsis:
578: #include <petscmath.h>
579:    int PetscSign(type v)

581:    Not Collective

583:    Input Parameter:
584: .  v - the number

586:    Notes:
587:     type can be integer or real floating point value

589:    Level: beginner

591: M*/
592: #define PetscSign(a) (((a) >= 0) ? ((a) == 0 ? 0 : 1) : -1)

594: /*MC
595:    PetscMin - Returns minimum of two numbers

597:    Synopsis:
598: #include <petscmath.h>
599:    type PetscMin(type v1,type v2)

601:    Not Collective

603:    Input Parameters:
604: +  v1 - first value to find minimum of
605: -  v2 - second value to find minimum of

607:    Notes:
608:     type can be integer or floating point value

610:    Level: beginner

612: .seealso: PetscMax(), PetscClipInterval(), PetscAbsInt(), PetscAbsReal(), PetscSqr()

614: M*/
615: #define PetscMin(a,b)   (((a)<(b)) ?  (a) : (b))

617: /*MC
618:    PetscMax - Returns maxium of two numbers

620:    Synopsis:
621: #include <petscmath.h>
622:    type max PetscMax(type v1,type v2)

624:    Not Collective

626:    Input Parameters:
627: +  v1 - first value to find maximum of
628: -  v2 - second value to find maximum of

630:    Notes:
631:     type can be integer or floating point value

633:    Level: beginner

635: .seealso: PetscMin(), PetscClipInterval(), PetscAbsInt(), PetscAbsReal(), PetscSqr()

637: M*/
638: #define PetscMax(a,b)   (((a)<(b)) ?  (b) : (a))

640: /*MC
641:    PetscClipInterval - Returns a number clipped to be within an interval

643:    Synopsis:
644: #include <petscmath.h>
645:    type clip PetscClipInterval(type x,type a,type b)

647:    Not Collective

649:    Input Parameters:
650: +  x - value to use if within interval [a,b]
651: .  a - lower end of interval
652: -  b - upper end of interval

654:    Notes:
655:     type can be integer or floating point value

657:    Level: beginner

659: .seealso: PetscMin(), PetscMax(), PetscAbsInt(), PetscAbsReal(), PetscSqr()

661: M*/
662: #define PetscClipInterval(x,a,b)   (PetscMax((a),PetscMin((x),(b))))

664: /*MC
665:    PetscAbsInt - Returns the absolute value of an integer

667:    Synopsis:
668: #include <petscmath.h>
669:    int abs PetscAbsInt(int v1)

671:    Not Collective

673:    Input Parameter:
674: .   v1 - the integer

676:    Level: beginner

678: .seealso: PetscMax(), PetscMin(), PetscAbsReal(), PetscSqr()

680: M*/
681: #define PetscAbsInt(a)  (((a)<0)   ? (-(a)) : (a))

683: /*MC
684:    PetscAbsReal - Returns the absolute value of an real number

686:    Synopsis:
687: #include <petscmath.h>
688:    Real abs PetscAbsReal(PetscReal v1)

690:    Not Collective

692:    Input Parameter:
693: .   v1 - the double

695:    Level: beginner

697: .seealso: PetscMax(), PetscMin(), PetscAbsInt(), PetscSqr()

699: M*/
700: #if defined(PETSC_USE_REAL_SINGLE)
701: #define PetscAbsReal(a) fabsf(a)
702: #elif defined(PETSC_USE_REAL_DOUBLE)
703: #define PetscAbsReal(a) fabs(a)
704: #elif defined(PETSC_USE_REAL___FLOAT128)
705: #define PetscAbsReal(a) fabsq(a)
706: #elif defined(PETSC_USE_REAL___FP16)
707: #define PetscAbsReal(a) fabsf(a)
708: #endif

710: /*MC
711:    PetscSqr - Returns the square of a number

713:    Synopsis:
714: #include <petscmath.h>
715:    type sqr PetscSqr(type v1)

717:    Not Collective

719:    Input Parameter:
720: .   v1 - the value

722:    Notes:
723:     type can be integer or floating point value

725:    Level: beginner

727: .seealso: PetscMax(), PetscMin(), PetscAbsInt(), PetscAbsReal()

729: M*/
730: #define PetscSqr(a)     ((a)*(a))

732: /* ----------------------------------------------------------------------------*/

734: #if defined(PETSC_USE_REAL_SINGLE)
735: #define PetscRealConstant(constant) constant##F
736: #elif defined(PETSC_USE_REAL_DOUBLE)
737: #define PetscRealConstant(constant) constant
738: #elif defined(PETSC_USE_REAL___FLOAT128)
739: #define PetscRealConstant(constant) constant##Q
740: #elif defined(PETSC_USE_REAL___FP16)
741: #define PetscRealConstant(constant) constant##F
742: #endif

744: /*
745:      Basic constants
746: */
747: #define PETSC_PI    PetscRealConstant(3.1415926535897932384626433832795029)
748: #define PETSC_PHI   PetscRealConstant(1.6180339887498948482045868343656381)
749: #define PETSC_SQRT2 PetscRealConstant(1.4142135623730950488016887242096981)

751: #if !defined(PETSC_USE_64BIT_INDICES)
752: #define PETSC_MAX_INT            2147483647
753: #define PETSC_MIN_INT            (-PETSC_MAX_INT - 1)
754: #else
755: #define PETSC_MAX_INT            9223372036854775807L
756: #define PETSC_MIN_INT            (-PETSC_MAX_INT - 1)
757: #endif
758: #define PETSC_MAX_UINT16         65535

760: #if defined(PETSC_USE_REAL_SINGLE)
761: #  define PETSC_MAX_REAL                3.40282346638528860e+38F
762: #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
763: #  define PETSC_MACHINE_EPSILON         1.19209290e-07F
764: #  define PETSC_SQRT_MACHINE_EPSILON    3.45266983e-04F
765: #  define PETSC_SMALL                   1.e-5F
766: #elif defined(PETSC_USE_REAL_DOUBLE)
767: #  define PETSC_MAX_REAL                1.7976931348623157e+308
768: #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
769: #  define PETSC_MACHINE_EPSILON         2.2204460492503131e-16
770: #  define PETSC_SQRT_MACHINE_EPSILON    1.490116119384766e-08
771: #  define PETSC_SMALL                   1.e-10
772: #elif defined(PETSC_USE_REAL___FLOAT128)
773: #  define PETSC_MAX_REAL                FLT128_MAX
774: #  define PETSC_MIN_REAL                (-FLT128_MAX)
775: #  define PETSC_MACHINE_EPSILON         FLT128_EPSILON
776: #  define PETSC_SQRT_MACHINE_EPSILON    1.38777878078144567552953958511352539e-17Q
777: #  define PETSC_SMALL                   1.e-20Q
778: #elif defined(PETSC_USE_REAL___FP16)
779: #  define PETSC_MAX_REAL                65504.0F
780: #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
781: #  define PETSC_MACHINE_EPSILON         .0009765625F
782: #  define PETSC_SQRT_MACHINE_EPSILON    .03125F
783: #  define PETSC_SMALL                   5.e-3F
784: #endif

786: #define PETSC_INFINITY               (PETSC_MAX_REAL/4)
787: #define PETSC_NINFINITY              (-PETSC_INFINITY)

789: PETSC_EXTERN PetscBool PetscIsInfReal(PetscReal);
790: PETSC_EXTERN PetscBool PetscIsNanReal(PetscReal);
791: PETSC_EXTERN PetscBool PetscIsNormalReal(PetscReal);
792: static inline PetscBool PetscIsInfOrNanReal(PetscReal v) {return PetscIsInfReal(v) || PetscIsNanReal(v) ? PETSC_TRUE : PETSC_FALSE;}
793: static inline PetscBool PetscIsInfScalar(PetscScalar v) {return PetscIsInfReal(PetscAbsScalar(v));}
794: static inline PetscBool PetscIsNanScalar(PetscScalar v) {return PetscIsNanReal(PetscAbsScalar(v));}
795: static inline PetscBool PetscIsInfOrNanScalar(PetscScalar v) {return PetscIsInfOrNanReal(PetscAbsScalar(v));}
796: static inline PetscBool PetscIsNormalScalar(PetscScalar v) {return PetscIsNormalReal(PetscAbsScalar(v));}

798: PETSC_EXTERN PetscBool PetscIsCloseAtTol(PetscReal,PetscReal,PetscReal,PetscReal);
799: PETSC_EXTERN PetscBool PetscEqualReal(PetscReal,PetscReal);
800: PETSC_EXTERN PetscBool PetscEqualScalar(PetscScalar,PetscScalar);

802: /*
803:     These macros are currently hardwired to match the regular data types, so there is no support for a different
804:     MatScalar from PetscScalar. We left the MatScalar in the source just in case we use it again.
805:  */
806: #define MPIU_MATSCALAR MPIU_SCALAR
807: typedef PetscScalar MatScalar;
808: typedef PetscReal MatReal;

810: struct petsc_mpiu_2scalar {PetscScalar a,b;};
811: PETSC_EXTERN MPI_Datatype MPIU_2SCALAR PetscAttrMPITypeTagLayoutCompatible(struct petsc_mpiu_2scalar);

813: /*
814:    MPI Datatypes for composite reductions:
815:    MPIU_REAL_INT -> struct { PetscReal; PetscInt; }
816:    MPIU_SCALAR_INT -> struct { PetscScalar; PetscInt; }
817: */
818: PETSC_EXTERN MPI_Datatype MPIU_REAL_INT;
819: PETSC_EXTERN MPI_Datatype MPIU_SCALAR_INT;

821: #if defined(PETSC_USE_64BIT_INDICES)
822: struct petsc_mpiu_2int {PetscInt a,b;};
823: PETSC_EXTERN MPI_Datatype MPIU_2INT PetscAttrMPITypeTagLayoutCompatible(struct petsc_mpiu_2int);
824: #else
825: #define MPIU_2INT MPI_2INT
826: #endif
827: PETSC_EXTERN MPI_Datatype MPI_4INT;
828: PETSC_EXTERN MPI_Datatype MPIU_4INT;

830: static inline PetscInt PetscPowInt(PetscInt base,PetscInt power)
831: {
832:   PetscInt result = 1;
833:   while (power) {
834:     if (power & 1) result *= base;
835:     power >>= 1;
836:     base *= base;
837:   }
838:   return result;
839: }

841: static inline PetscInt64 PetscPowInt64(PetscInt base,PetscInt power)
842: {
843:   PetscInt64 result = 1;
844:   while (power) {
845:     if (power & 1) result *= base;
846:     power >>= 1;
847:     base *= base;
848:   }
849:   return result;
850: }

852: static inline PetscReal PetscPowRealInt(PetscReal base,PetscInt power)
853: {
854:   PetscReal result = 1;
855:   if (power < 0) {
856:     power = -power;
857:     base  = ((PetscReal)1)/base;
858:   }
859:   while (power) {
860:     if (power & 1) result *= base;
861:     power >>= 1;
862:     base *= base;
863:   }
864:   return result;
865: }

867: static inline PetscScalar PetscPowScalarInt(PetscScalar base,PetscInt power)
868: {
869:   PetscScalar result = (PetscReal)1;
870:   if (power < 0) {
871:     power = -power;
872:     base  = ((PetscReal)1)/base;
873:   }
874:   while (power) {
875:     if (power & 1) result *= base;
876:     power >>= 1;
877:     base *= base;
878:   }
879:   return result;
880: }

882: static inline PetscScalar PetscPowScalarReal(PetscScalar base,PetscReal power)
883: {
884:   PetscScalar cpower = power;
885:   return PetscPowScalar(base,cpower);
886: }

888: /*MC
889:     PetscApproximateLTE - Performs a less than or equal to on a given constant with a fudge for floating point numbers

891:    Synopsis:
892: #include <petscmath.h>
893:    bool PetscApproximateLTE(PetscReal x,constant float)

895:    Not Collective

897:    Input Parameters:
898: +   x - the variable
899: -   b - the constant float it is checking if x is less than or equal to

901:    Notes:
902:      The fudge factor is the value PETSC_SMALL

904:      The constant numerical value is automatically set to the appropriate precision of PETSc so can just be provided as, for example, 3.2

906:      This is used in several examples for setting initial conditions based on coordinate values that are computed with i*h that produces inexact
907:      floating point results.

909:    Level: advanced

911: .seealso: PetscMax(), PetscMin(), PetscAbsInt(), PetscAbsReal(), PetscApproximateGTE()

913: M*/
914: #define PetscApproximateLTE(x,b)  ((x) <= (PetscRealConstant(b)+PETSC_SMALL))

916: /*MC
917:     PetscApproximateGTE - Performs a greater than or equal to on a given constant with a fudge for floating point numbers

919:    Synopsis:
920: #include <petscmath.h>
921:    bool PetscApproximateGTE(PetscReal x,constant float)

923:    Not Collective

925:    Input Parameters:
926: +   x - the variable
927: -   b - the constant float it is checking if x is greater than or equal to

929:    Notes:
930:      The fudge factor is the value PETSC_SMALL

932:      The constant numerical value is automatically set to the appropriate precision of PETSc so can just be provided as, for example, 3.2

934:      This is used in several examples for setting initial conditions based on coordinate values that are computed with i*h that produces inexact
935:      floating point results.

937:    Level: advanced

939: .seealso: PetscMax(), PetscMin(), PetscAbsInt(), PetscAbsReal(), PetscApproximateLTE()

941: M*/
942: #define PetscApproximateGTE(x,b)  ((x) >= (PetscRealConstant(b)-PETSC_SMALL))

944: /*MC
945:     PetscCeilInt - Returns the ceiling of the quotation of two positive integers

947:    Synopsis:
948: #include <petscmath.h>
949:    PetscInt PetscCeilInt(PetscInt x,PetscInt y)

951:    Not Collective

953:    Input Parameters:
954: +   x - the numerator
955: -   y - the denominator

957:    Level: advanced

959: .seealso: PetscMax(), PetscMin(), PetscAbsInt(), PetscAbsReal(), PetscApproximateLTE()

961: M*/
962: #define PetscCeilInt(x,y)  ((((PetscInt)(x))/((PetscInt)(y))) +  ((((PetscInt)(x)) % ((PetscInt)(y))) ? 1 : 0))

964: #define PetscCeilInt64(x,y)  ((((PetscInt64)(x))/((PetscInt64)(y))) +  ((((PetscInt64)(x)) % ((PetscInt64)(y))) ? 1 : 0))

966: PETSC_EXTERN PetscErrorCode PetscLinearRegression(PetscInt,const PetscReal[],const PetscReal[],PetscReal*,PetscReal*);
967: #endif