xref: /petsc/include/petscmath.h (revision ac09b9214d23ea9ad238aa607de9fa447fd4e91b)
1e489efc1SBarry Smith /*
2314da920SBarry Smith 
3314da920SBarry Smith     PETSc mathematics include file. Defines certain basic mathematical
4a5057860SBarry Smith     constants and functions for working with single, double, and quad precision
5a5057860SBarry Smith     floating point numbers as well as complex single and double.
6314da920SBarry Smith 
7d382aafbSBarry Smith     This file is included by petscsys.h and should not be used directly.
8e7029fe1SSatish Balay 
9e489efc1SBarry Smith */
1026bd1501SBarry Smith #if !defined(PETSCMATH_H)
1126bd1501SBarry Smith #define PETSCMATH_H
12*ac09b921SBarry Smith 
130a5f7794SBarry Smith #include <math.h>
14df4397b0SStefano Zampini #include <petscsystypes.h>
15df4397b0SStefano Zampini 
16*ac09b921SBarry Smith /* SUBMANSEC = Sys */
17*ac09b921SBarry Smith 
185117d392SLisandro Dalcin /*
195117d392SLisandro Dalcin 
205117d392SLisandro Dalcin    Defines operations that are different for complex and real numbers.
215117d392SLisandro Dalcin    All PETSc objects in one program are built around the object
225117d392SLisandro Dalcin    PetscScalar which is either always a real or a complex.
235117d392SLisandro Dalcin 
245117d392SLisandro Dalcin */
255117d392SLisandro Dalcin 
265117d392SLisandro Dalcin /*
275117d392SLisandro Dalcin     Real number definitions
285117d392SLisandro Dalcin  */
295117d392SLisandro Dalcin #if defined(PETSC_USE_REAL_SINGLE)
305117d392SLisandro Dalcin #define PetscSqrtReal(a)       sqrtf(a)
315117d392SLisandro Dalcin #define PetscCbrtReal(a)       cbrtf(a)
325117d392SLisandro Dalcin #define PetscHypotReal(a,b)    hypotf(a,b)
335117d392SLisandro Dalcin #define PetscAtan2Real(a,b)    atan2f(a,b)
345117d392SLisandro Dalcin #define PetscPowReal(a,b)      powf(a,b)
355117d392SLisandro Dalcin #define PetscExpReal(a)        expf(a)
365117d392SLisandro Dalcin #define PetscLogReal(a)        logf(a)
375117d392SLisandro Dalcin #define PetscLog10Real(a)      log10f(a)
385117d392SLisandro Dalcin #define PetscLog2Real(a)       log2f(a)
395117d392SLisandro Dalcin #define PetscSinReal(a)        sinf(a)
405117d392SLisandro Dalcin #define PetscCosReal(a)        cosf(a)
415117d392SLisandro Dalcin #define PetscTanReal(a)        tanf(a)
425117d392SLisandro Dalcin #define PetscAsinReal(a)       asinf(a)
435117d392SLisandro Dalcin #define PetscAcosReal(a)       acosf(a)
445117d392SLisandro Dalcin #define PetscAtanReal(a)       atanf(a)
455117d392SLisandro Dalcin #define PetscSinhReal(a)       sinhf(a)
465117d392SLisandro Dalcin #define PetscCoshReal(a)       coshf(a)
475117d392SLisandro Dalcin #define PetscTanhReal(a)       tanhf(a)
485117d392SLisandro Dalcin #define PetscAsinhReal(a)      asinhf(a)
495117d392SLisandro Dalcin #define PetscAcoshReal(a)      acoshf(a)
505117d392SLisandro Dalcin #define PetscAtanhReal(a)      atanhf(a)
51d6685f55SMatthew G. Knepley #define PetscErfReal(a)        erff(a)
525117d392SLisandro Dalcin #define PetscCeilReal(a)       ceilf(a)
535117d392SLisandro Dalcin #define PetscFloorReal(a)      floorf(a)
545117d392SLisandro Dalcin #define PetscFmodReal(a,b)     fmodf(a,b)
559c3ee494SJed Brown #define PetscCopysignReal(a,b) copysignf(a,b)
565117d392SLisandro Dalcin #define PetscTGamma(a)         tgammaf(a)
571f17fa70SToby Isaac #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
581f17fa70SToby Isaac #define PetscLGamma(a)         gammaf(a)
591f17fa70SToby Isaac #else
601f17fa70SToby Isaac #define PetscLGamma(a)         lgammaf(a)
611f17fa70SToby Isaac #endif
625117d392SLisandro Dalcin 
635117d392SLisandro Dalcin #elif defined(PETSC_USE_REAL_DOUBLE)
645117d392SLisandro Dalcin #define PetscSqrtReal(a)       sqrt(a)
655117d392SLisandro Dalcin #define PetscCbrtReal(a)       cbrt(a)
665117d392SLisandro Dalcin #define PetscHypotReal(a,b)    hypot(a,b)
675117d392SLisandro Dalcin #define PetscAtan2Real(a,b)    atan2(a,b)
685117d392SLisandro Dalcin #define PetscPowReal(a,b)      pow(a,b)
695117d392SLisandro Dalcin #define PetscExpReal(a)        exp(a)
705117d392SLisandro Dalcin #define PetscLogReal(a)        log(a)
715117d392SLisandro Dalcin #define PetscLog10Real(a)      log10(a)
725117d392SLisandro Dalcin #define PetscLog2Real(a)       log2(a)
735117d392SLisandro Dalcin #define PetscSinReal(a)        sin(a)
745117d392SLisandro Dalcin #define PetscCosReal(a)        cos(a)
755117d392SLisandro Dalcin #define PetscTanReal(a)        tan(a)
765117d392SLisandro Dalcin #define PetscAsinReal(a)       asin(a)
775117d392SLisandro Dalcin #define PetscAcosReal(a)       acos(a)
785117d392SLisandro Dalcin #define PetscAtanReal(a)       atan(a)
795117d392SLisandro Dalcin #define PetscSinhReal(a)       sinh(a)
805117d392SLisandro Dalcin #define PetscCoshReal(a)       cosh(a)
815117d392SLisandro Dalcin #define PetscTanhReal(a)       tanh(a)
825117d392SLisandro Dalcin #define PetscAsinhReal(a)      asinh(a)
835117d392SLisandro Dalcin #define PetscAcoshReal(a)      acosh(a)
845117d392SLisandro Dalcin #define PetscAtanhReal(a)      atanh(a)
85d6685f55SMatthew G. Knepley #define PetscErfReal(a)        erf(a)
865117d392SLisandro Dalcin #define PetscCeilReal(a)       ceil(a)
875117d392SLisandro Dalcin #define PetscFloorReal(a)      floor(a)
885117d392SLisandro Dalcin #define PetscFmodReal(a,b)     fmod(a,b)
899c3ee494SJed Brown #define PetscCopysignReal(a,b) copysign(a,b)
905117d392SLisandro Dalcin #define PetscTGamma(a)         tgamma(a)
911f17fa70SToby Isaac #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
921f17fa70SToby Isaac #define PetscLGamma(a)         gamma(a)
931f17fa70SToby Isaac #else
941f17fa70SToby Isaac #define PetscLGamma(a)         lgamma(a)
951f17fa70SToby Isaac #endif
965117d392SLisandro Dalcin 
975117d392SLisandro Dalcin #elif defined(PETSC_USE_REAL___FLOAT128)
985117d392SLisandro Dalcin #define PetscSqrtReal(a)       sqrtq(a)
995117d392SLisandro Dalcin #define PetscCbrtReal(a)       cbrtq(a)
1005117d392SLisandro Dalcin #define PetscHypotReal(a,b)    hypotq(a,b)
1015117d392SLisandro Dalcin #define PetscAtan2Real(a,b)    atan2q(a,b)
1025117d392SLisandro Dalcin #define PetscPowReal(a,b)      powq(a,b)
1035117d392SLisandro Dalcin #define PetscExpReal(a)        expq(a)
1045117d392SLisandro Dalcin #define PetscLogReal(a)        logq(a)
1055117d392SLisandro Dalcin #define PetscLog10Real(a)      log10q(a)
1065117d392SLisandro Dalcin #define PetscLog2Real(a)       log2q(a)
1075117d392SLisandro Dalcin #define PetscSinReal(a)        sinq(a)
1085117d392SLisandro Dalcin #define PetscCosReal(a)        cosq(a)
1095117d392SLisandro Dalcin #define PetscTanReal(a)        tanq(a)
1105117d392SLisandro Dalcin #define PetscAsinReal(a)       asinq(a)
1115117d392SLisandro Dalcin #define PetscAcosReal(a)       acosq(a)
1125117d392SLisandro Dalcin #define PetscAtanReal(a)       atanq(a)
1135117d392SLisandro Dalcin #define PetscSinhReal(a)       sinhq(a)
1145117d392SLisandro Dalcin #define PetscCoshReal(a)       coshq(a)
1155117d392SLisandro Dalcin #define PetscTanhReal(a)       tanhq(a)
1165117d392SLisandro Dalcin #define PetscAsinhReal(a)      asinhq(a)
1175117d392SLisandro Dalcin #define PetscAcoshReal(a)      acoshq(a)
1185117d392SLisandro Dalcin #define PetscAtanhReal(a)      atanhq(a)
119d6685f55SMatthew G. Knepley #define PetscErfReal(a)        erfq(a)
1205117d392SLisandro Dalcin #define PetscCeilReal(a)       ceilq(a)
1215117d392SLisandro Dalcin #define PetscFloorReal(a)      floorq(a)
1225117d392SLisandro Dalcin #define PetscFmodReal(a,b)     fmodq(a,b)
1239c3ee494SJed Brown #define PetscCopysignReal(a,b) copysignq(a,b)
1245117d392SLisandro Dalcin #define PetscTGamma(a)         tgammaq(a)
1251f17fa70SToby Isaac #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
1261f17fa70SToby Isaac #define PetscLGamma(a)         gammaq(a)
1271f17fa70SToby Isaac #else
1281f17fa70SToby Isaac #define PetscLGamma(a)         lgammaq(a)
1291f17fa70SToby Isaac #endif
1305117d392SLisandro Dalcin 
1315117d392SLisandro Dalcin #elif defined(PETSC_USE_REAL___FP16)
1325117d392SLisandro Dalcin #define PetscSqrtReal(a)       sqrtf(a)
1335117d392SLisandro Dalcin #define PetscCbrtReal(a)       cbrtf(a)
1345117d392SLisandro Dalcin #define PetscHypotReal(a,b)    hypotf(a,b)
1355117d392SLisandro Dalcin #define PetscAtan2Real(a,b)    atan2f(a,b)
1365117d392SLisandro Dalcin #define PetscPowReal(a,b)      powf(a,b)
1375117d392SLisandro Dalcin #define PetscExpReal(a)        expf(a)
1385117d392SLisandro Dalcin #define PetscLogReal(a)        logf(a)
1395117d392SLisandro Dalcin #define PetscLog10Real(a)      log10f(a)
1405117d392SLisandro Dalcin #define PetscLog2Real(a)       log2f(a)
1415117d392SLisandro Dalcin #define PetscSinReal(a)        sinf(a)
1425117d392SLisandro Dalcin #define PetscCosReal(a)        cosf(a)
1435117d392SLisandro Dalcin #define PetscTanReal(a)        tanf(a)
1445117d392SLisandro Dalcin #define PetscAsinReal(a)       asinf(a)
1455117d392SLisandro Dalcin #define PetscAcosReal(a)       acosf(a)
1465117d392SLisandro Dalcin #define PetscAtanReal(a)       atanf(a)
1475117d392SLisandro Dalcin #define PetscSinhReal(a)       sinhf(a)
1485117d392SLisandro Dalcin #define PetscCoshReal(a)       coshf(a)
1495117d392SLisandro Dalcin #define PetscTanhReal(a)       tanhf(a)
1505117d392SLisandro Dalcin #define PetscAsinhReal(a)      asinhf(a)
1515117d392SLisandro Dalcin #define PetscAcoshReal(a)      acoshf(a)
1525117d392SLisandro Dalcin #define PetscAtanhReal(a)      atanhf(a)
153d6685f55SMatthew G. Knepley #define PetscErfReal(a)        erff(a)
1545117d392SLisandro Dalcin #define PetscCeilReal(a)       ceilf(a)
1555117d392SLisandro Dalcin #define PetscFloorReal(a)      floorf(a)
1565117d392SLisandro Dalcin #define PetscFmodReal(a,b)     fmodf(a,b)
1579c3ee494SJed Brown #define PetscCopySignReal(a,b) copysignf(a,b)
1585117d392SLisandro Dalcin #define PetscTGamma(a)         tgammaf(a)
1591f17fa70SToby Isaac #if defined(PETSC_HAVE_LGAMMA_IS_GAMMA)
1601f17fa70SToby Isaac #define PetscLGamma(a)         gammaf(a)
1611f17fa70SToby Isaac #else
1621f17fa70SToby Isaac #define PetscLGamma(a)         lgammaf(a)
1631f17fa70SToby Isaac #endif
1645117d392SLisandro Dalcin 
1655117d392SLisandro Dalcin #endif /* PETSC_USE_REAL_* */
1665117d392SLisandro Dalcin 
1679fbee547SJacob Faibussowitsch static inline PetscReal PetscSignReal(PetscReal a)
1685117d392SLisandro Dalcin {
1695117d392SLisandro Dalcin   return (PetscReal)((a < (PetscReal)0) ? -1 : ((a > (PetscReal)0) ? 1 : 0));
1705117d392SLisandro Dalcin }
1715117d392SLisandro Dalcin 
1725117d392SLisandro Dalcin #if !defined(PETSC_HAVE_LOG2)
1735117d392SLisandro Dalcin #undef PetscLog2Real
1749fbee547SJacob Faibussowitsch static inline PetscReal PetscLog2Real(PetscReal a)
1755117d392SLisandro Dalcin {
1765117d392SLisandro Dalcin   return PetscLogReal(a)/PetscLogReal((PetscReal)2);
1775117d392SLisandro Dalcin }
1785117d392SLisandro Dalcin #endif
1795117d392SLisandro Dalcin 
1805117d392SLisandro Dalcin #if defined(PETSC_USE_REAL___FLOAT128)
1815117d392SLisandro Dalcin PETSC_EXTERN MPI_Datatype MPIU___FLOAT128 PetscAttrMPITypeTag(__float128);
1825117d392SLisandro Dalcin #endif
1835117d392SLisandro Dalcin #if defined(PETSC_USE_REAL___FP16)
1845117d392SLisandro Dalcin PETSC_EXTERN MPI_Datatype MPIU___FP16 PetscAttrMPITypeTag(__fp16);
1855117d392SLisandro Dalcin #endif
1865117d392SLisandro Dalcin 
187df4397b0SStefano Zampini /*MC
188df4397b0SStefano Zampini    MPIU_REAL - MPI datatype corresponding to PetscReal
189df4397b0SStefano Zampini 
190df4397b0SStefano Zampini    Notes:
191df4397b0SStefano Zampini    In MPI calls that require an MPI datatype that matches a PetscReal or array of PetscReal values, pass this value.
192df4397b0SStefano Zampini 
193df4397b0SStefano Zampini    Level: beginner
194df4397b0SStefano Zampini 
195db781477SPatrick Sanan .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_SCALAR`, `MPIU_COMPLEX`, `MPIU_INT`
196df4397b0SStefano Zampini M*/
197c1d390e3SJed Brown #if defined(PETSC_USE_REAL_SINGLE)
198c1d390e3SJed Brown #  define MPIU_REAL MPI_FLOAT
199c1d390e3SJed Brown #elif defined(PETSC_USE_REAL_DOUBLE)
200c1d390e3SJed Brown #  define MPIU_REAL MPI_DOUBLE
201c1d390e3SJed Brown #elif defined(PETSC_USE_REAL___FLOAT128)
202c1d390e3SJed Brown #  define MPIU_REAL MPIU___FLOAT128
203570b7f6dSBarry Smith #elif defined(PETSC_USE_REAL___FP16)
204570b7f6dSBarry Smith #  define MPIU_REAL MPIU___FP16
205c1d390e3SJed Brown #endif /* PETSC_USE_REAL_* */
20659cb5930SBarry Smith 
2071093a601SBarry Smith /*
2081093a601SBarry Smith     Complex number definitions
2091093a601SBarry Smith  */
210df4397b0SStefano Zampini #if defined(PETSC_HAVE_COMPLEX)
211450fc7c9SSatish Balay #if defined(__cplusplus) && !defined(PETSC_USE_REAL___FLOAT128)
2121093a601SBarry Smith /* C++ support of complex number */
213b7940d39SSatish Balay 
21450f81f78SJed Brown #define PetscRealPartComplex(a)      (a).real()
21550f81f78SJed Brown #define PetscImaginaryPartComplex(a) (a).imag()
216df4397b0SStefano Zampini #define PetscAbsComplex(a)           petsccomplexlib::abs(a)
2175117d392SLisandro Dalcin #define PetscArgComplex(a)           petsccomplexlib::arg(a)
218df4397b0SStefano Zampini #define PetscConjComplex(a)          petsccomplexlib::conj(a)
219df4397b0SStefano Zampini #define PetscSqrtComplex(a)          petsccomplexlib::sqrt(a)
220df4397b0SStefano Zampini #define PetscPowComplex(a,b)         petsccomplexlib::pow(a,b)
221df4397b0SStefano Zampini #define PetscExpComplex(a)           petsccomplexlib::exp(a)
222df4397b0SStefano Zampini #define PetscLogComplex(a)           petsccomplexlib::log(a)
223df4397b0SStefano Zampini #define PetscSinComplex(a)           petsccomplexlib::sin(a)
224df4397b0SStefano Zampini #define PetscCosComplex(a)           petsccomplexlib::cos(a)
2255117d392SLisandro Dalcin #define PetscTanComplex(a)           petsccomplexlib::tan(a)
226df4397b0SStefano Zampini #define PetscAsinComplex(a)          petsccomplexlib::asin(a)
227df4397b0SStefano Zampini #define PetscAcosComplex(a)          petsccomplexlib::acos(a)
2285117d392SLisandro Dalcin #define PetscAtanComplex(a)          petsccomplexlib::atan(a)
229df4397b0SStefano Zampini #define PetscSinhComplex(a)          petsccomplexlib::sinh(a)
230df4397b0SStefano Zampini #define PetscCoshComplex(a)          petsccomplexlib::cosh(a)
231df4397b0SStefano Zampini #define PetscTanhComplex(a)          petsccomplexlib::tanh(a)
2325117d392SLisandro Dalcin #define PetscAsinhComplex(a)         petsccomplexlib::asinh(a)
2335117d392SLisandro Dalcin #define PetscAcoshComplex(a)         petsccomplexlib::acosh(a)
2345117d392SLisandro Dalcin #define PetscAtanhComplex(a)         petsccomplexlib::atanh(a)
2355117d392SLisandro Dalcin 
2365117d392SLisandro Dalcin /* TODO: Add configure tests
2375117d392SLisandro Dalcin 
2385117d392SLisandro Dalcin #if !defined(PETSC_HAVE_CXX_TAN_COMPLEX)
2395117d392SLisandro Dalcin #undef PetscTanComplex
2409fbee547SJacob Faibussowitsch static inline PetscComplex PetscTanComplex(PetscComplex z)
2415117d392SLisandro Dalcin {
2425117d392SLisandro Dalcin   return PetscSinComplex(z)/PetscCosComplex(z);
2435117d392SLisandro Dalcin }
244027d9794SBarry Smith #endif
245debe9ee2SPaul Mullowney 
2465117d392SLisandro Dalcin #if !defined(PETSC_HAVE_CXX_TANH_COMPLEX)
2475117d392SLisandro Dalcin #undef PetscTanhComplex
2489fbee547SJacob Faibussowitsch static inline PetscComplex PetscTanhComplex(PetscComplex z)
2495117d392SLisandro Dalcin {
2505117d392SLisandro Dalcin   return PetscSinhComplex(z)/PetscCoshComplex(z);
2515117d392SLisandro Dalcin }
2525117d392SLisandro Dalcin #endif
2535117d392SLisandro Dalcin 
2545117d392SLisandro Dalcin #if !defined(PETSC_HAVE_CXX_ASIN_COMPLEX)
2555117d392SLisandro Dalcin #undef PetscAsinComplex
2569fbee547SJacob Faibussowitsch static inline PetscComplex PetscAsinComplex(PetscComplex z)
2575117d392SLisandro Dalcin {
2585117d392SLisandro Dalcin   const PetscComplex j(0,1);
2595117d392SLisandro Dalcin   return -j*PetscLogComplex(j*z+PetscSqrtComplex(1.0f-z*z));
2605117d392SLisandro Dalcin }
2615117d392SLisandro Dalcin #endif
2625117d392SLisandro Dalcin 
2635117d392SLisandro Dalcin #if !defined(PETSC_HAVE_CXX_ACOS_COMPLEX)
2645117d392SLisandro Dalcin #undef PetscAcosComplex
2659fbee547SJacob Faibussowitsch static inline PetscComplex PetscAcosComplex(PetscComplex z)
2665117d392SLisandro Dalcin {
2675117d392SLisandro Dalcin   const PetscComplex j(0,1);
2685117d392SLisandro Dalcin   return j*PetscLogComplex(z-j*PetscSqrtComplex(1.0f-z*z));
2695117d392SLisandro Dalcin }
2705117d392SLisandro Dalcin #endif
2715117d392SLisandro Dalcin 
2725117d392SLisandro Dalcin #if !defined(PETSC_HAVE_CXX_ATAN_COMPLEX)
2735117d392SLisandro Dalcin #undef PetscAtanComplex
2749fbee547SJacob Faibussowitsch static inline PetscComplex PetscAtanComplex(PetscComplex z)
2755117d392SLisandro Dalcin {
2765117d392SLisandro Dalcin   const PetscComplex j(0,1);
2775117d392SLisandro Dalcin   return 0.5f*j*PetscLogComplex((1.0f-j*z)/(1.0f+j*z));
2785117d392SLisandro Dalcin }
2795117d392SLisandro Dalcin #endif
2805117d392SLisandro Dalcin 
2815117d392SLisandro Dalcin #if !defined(PETSC_HAVE_CXX_ASINH_COMPLEX)
2825117d392SLisandro Dalcin #undef PetscAsinhComplex
2839fbee547SJacob Faibussowitsch static inline PetscComplex PetscAsinhComplex(PetscComplex z)
2845117d392SLisandro Dalcin {
2855117d392SLisandro Dalcin   return PetscLogComplex(z+PetscSqrtComplex(z*z+1.0f));
2865117d392SLisandro Dalcin }
2875117d392SLisandro Dalcin #endif
2885117d392SLisandro Dalcin 
2895117d392SLisandro Dalcin #if !defined(PETSC_HAVE_CXX_ACOSH_COMPLEX)
2905117d392SLisandro Dalcin #undef PetscAcoshComplex
2919fbee547SJacob Faibussowitsch static inline PetscComplex PetscAcoshComplex(PetscComplex z)
2925117d392SLisandro Dalcin {
2935117d392SLisandro Dalcin   return PetscLogComplex(z+PetscSqrtComplex(z*z-1.0f));
2945117d392SLisandro Dalcin }
2955117d392SLisandro Dalcin #endif
2965117d392SLisandro Dalcin 
2975117d392SLisandro Dalcin #if !defined(PETSC_HAVE_CXX_ATANH_COMPLEX)
2985117d392SLisandro Dalcin #undef PetscAtanhComplex
2999fbee547SJacob Faibussowitsch static inline PetscComplex PetscAtanhComplex(PetscComplex z)
3005117d392SLisandro Dalcin {
3015117d392SLisandro Dalcin   return 0.5f*PetscLogComplex((1.0f+z)/(1.0f-z));
3025117d392SLisandro Dalcin }
3035117d392SLisandro Dalcin #endif
3045117d392SLisandro Dalcin 
3055117d392SLisandro Dalcin */
3065117d392SLisandro Dalcin 
3077a19d461SSatish Balay #else /* C99 support of complex number */
308519e2a1fSPaul Mullowney 
3097a19d461SSatish Balay #if defined(PETSC_USE_REAL_SINGLE)
31050f81f78SJed Brown #define PetscRealPartComplex(a)      crealf(a)
31150f81f78SJed Brown #define PetscImaginaryPartComplex(a) cimagf(a)
31250f81f78SJed Brown #define PetscAbsComplex(a)           cabsf(a)
3135117d392SLisandro Dalcin #define PetscArgComplex(a)           cargf(a)
31450f81f78SJed Brown #define PetscConjComplex(a)          conjf(a)
31550f81f78SJed Brown #define PetscSqrtComplex(a)          csqrtf(a)
31650f81f78SJed Brown #define PetscPowComplex(a,b)         cpowf(a,b)
31750f81f78SJed Brown #define PetscExpComplex(a)           cexpf(a)
31850f81f78SJed Brown #define PetscLogComplex(a)           clogf(a)
31950f81f78SJed Brown #define PetscSinComplex(a)           csinf(a)
32050f81f78SJed Brown #define PetscCosComplex(a)           ccosf(a)
3215117d392SLisandro Dalcin #define PetscTanComplex(a)           ctanf(a)
322255453a1SBarry Smith #define PetscAsinComplex(a)          casinf(a)
323255453a1SBarry Smith #define PetscAcosComplex(a)          cacosf(a)
3245117d392SLisandro Dalcin #define PetscAtanComplex(a)          catanf(a)
325a4bea5a6SPeter Brune #define PetscSinhComplex(a)          csinhf(a)
326a4bea5a6SPeter Brune #define PetscCoshComplex(a)          ccoshf(a)
327a4bea5a6SPeter Brune #define PetscTanhComplex(a)          ctanhf(a)
3285117d392SLisandro Dalcin #define PetscAsinhComplex(a)         casinhf(a)
3295117d392SLisandro Dalcin #define PetscAcoshComplex(a)         cacoshf(a)
3305117d392SLisandro Dalcin #define PetscAtanhComplex(a)         catanhf(a)
3311093a601SBarry Smith 
332ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
33350f81f78SJed Brown #define PetscRealPartComplex(a)      creal(a)
33450f81f78SJed Brown #define PetscImaginaryPartComplex(a) cimag(a)
33550f81f78SJed Brown #define PetscAbsComplex(a)           cabs(a)
3365117d392SLisandro Dalcin #define PetscArgComplex(a)           carg(a)
33750f81f78SJed Brown #define PetscConjComplex(a)          conj(a)
33850f81f78SJed Brown #define PetscSqrtComplex(a)          csqrt(a)
33950f81f78SJed Brown #define PetscPowComplex(a,b)         cpow(a,b)
34050f81f78SJed Brown #define PetscExpComplex(a)           cexp(a)
34150f81f78SJed Brown #define PetscLogComplex(a)           clog(a)
34250f81f78SJed Brown #define PetscSinComplex(a)           csin(a)
34350f81f78SJed Brown #define PetscCosComplex(a)           ccos(a)
3445117d392SLisandro Dalcin #define PetscTanComplex(a)           ctan(a)
345255453a1SBarry Smith #define PetscAsinComplex(a)          casin(a)
346255453a1SBarry Smith #define PetscAcosComplex(a)          cacos(a)
3475117d392SLisandro Dalcin #define PetscAtanComplex(a)          catan(a)
348a4bea5a6SPeter Brune #define PetscSinhComplex(a)          csinh(a)
349a4bea5a6SPeter Brune #define PetscCoshComplex(a)          ccosh(a)
350a4bea5a6SPeter Brune #define PetscTanhComplex(a)          ctanh(a)
3515117d392SLisandro Dalcin #define PetscAsinhComplex(a)         casinh(a)
3525117d392SLisandro Dalcin #define PetscAcoshComplex(a)         cacosh(a)
3535117d392SLisandro Dalcin #define PetscAtanhComplex(a)         catanh(a)
3541093a601SBarry Smith 
3558c764dc5SJose Roman #elif defined(PETSC_USE_REAL___FLOAT128)
35650f81f78SJed Brown #define PetscRealPartComplex(a)      crealq(a)
35750f81f78SJed Brown #define PetscImaginaryPartComplex(a) cimagq(a)
35850f81f78SJed Brown #define PetscAbsComplex(a)           cabsq(a)
3595117d392SLisandro Dalcin #define PetscArgComplex(a)           cargq(a)
36050f81f78SJed Brown #define PetscConjComplex(a)          conjq(a)
36150f81f78SJed Brown #define PetscSqrtComplex(a)          csqrtq(a)
36250f81f78SJed Brown #define PetscPowComplex(a,b)         cpowq(a,b)
36350f81f78SJed Brown #define PetscExpComplex(a)           cexpq(a)
36450f81f78SJed Brown #define PetscLogComplex(a)           clogq(a)
36550f81f78SJed Brown #define PetscSinComplex(a)           csinq(a)
36650f81f78SJed Brown #define PetscCosComplex(a)           ccosq(a)
3675117d392SLisandro Dalcin #define PetscTanComplex(a)           ctanq(a)
368255453a1SBarry Smith #define PetscAsinComplex(a)          casinq(a)
369255453a1SBarry Smith #define PetscAcosComplex(a)          cacosq(a)
3705117d392SLisandro Dalcin #define PetscAtanComplex(a)          catanq(a)
371a4bea5a6SPeter Brune #define PetscSinhComplex(a)          csinhq(a)
372a4bea5a6SPeter Brune #define PetscCoshComplex(a)          ccoshq(a)
373a4bea5a6SPeter Brune #define PetscTanhComplex(a)          ctanhq(a)
3745117d392SLisandro Dalcin #define PetscAsinhComplex(a)         casinhq(a)
3755117d392SLisandro Dalcin #define PetscAcoshComplex(a)         cacoshq(a)
3765117d392SLisandro Dalcin #define PetscAtanhComplex(a)         catanhq(a)
377a4bea5a6SPeter Brune 
378ce63c4c1SBarry Smith #endif /* PETSC_USE_REAL_* */
3797a19d461SSatish Balay #endif /* (__cplusplus) */
380e489efc1SBarry Smith 
3811093a601SBarry Smith /*
3825117d392SLisandro Dalcin    PETSC_i is the imaginary number, i
3831093a601SBarry Smith */
38450f81f78SJed Brown PETSC_EXTERN PetscComplex PETSC_i;
3858a351411SToby Isaac 
3865117d392SLisandro Dalcin /*
3875117d392SLisandro Dalcin    Try to do the right thing for complex number construction: see
3888a351411SToby Isaac    http://www.open-std.org/jtc1/sc22/wg14/www/docs/n1464.htm
3898a351411SToby Isaac    for details
3908a351411SToby Isaac */
3919fbee547SJacob Faibussowitsch static inline PetscComplex PetscCMPLX(PetscReal x, PetscReal y)
3928a351411SToby Isaac {
393450fc7c9SSatish Balay #if defined(__cplusplus) && !defined(PETSC_USE_REAL___FLOAT128)
3948a351411SToby Isaac   return PetscComplex(x,y);
3958a351411SToby Isaac #elif defined(_Imaginary_I)
3968a351411SToby Isaac   return x + y * _Imaginary_I;
3978a351411SToby Isaac #else
398616d7c5eSToby Isaac   { /* In both C99 and C11 (ISO/IEC 9899, Section 6.2.5),
399616d7c5eSToby Isaac 
400616d7c5eSToby Isaac        "For each floating type there is a corresponding real type, which is always a real floating
401616d7c5eSToby Isaac        type. For real floating types, it is the same type. For complex types, it is the type given
402616d7c5eSToby Isaac        by deleting the keyword _Complex from the type name."
403616d7c5eSToby Isaac 
404616d7c5eSToby Isaac        So type punning should be portable. */
405616d7c5eSToby Isaac     union { PetscComplex z; PetscReal f[2]; } uz;
406616d7c5eSToby Isaac 
407616d7c5eSToby Isaac     uz.f[0] = x;
408616d7c5eSToby Isaac     uz.f[1] = y;
409616d7c5eSToby Isaac     return uz.z;
410616d7c5eSToby Isaac   }
41150f81f78SJed Brown #endif
4128a351411SToby Isaac }
4138a351411SToby Isaac 
414de272c7aSSatish Balay #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)\"")
415de272c7aSSatish Balay #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)\"")
416de272c7aSSatish Balay 
4175117d392SLisandro Dalcin #if defined(PETSC_USE_REAL___FLOAT128)
4185117d392SLisandro Dalcin PETSC_EXTERN MPI_Datatype MPIU___COMPLEX128 PetscAttrMPITypeTag(__complex128);
4195117d392SLisandro Dalcin #endif /* PETSC_USE_REAL___FLOAT128 */
4205117d392SLisandro Dalcin 
4215117d392SLisandro Dalcin /*MC
4225117d392SLisandro Dalcin    MPIU_COMPLEX - MPI datatype corresponding to PetscComplex
4235117d392SLisandro Dalcin 
4245117d392SLisandro Dalcin    Notes:
4255117d392SLisandro Dalcin    In MPI calls that require an MPI datatype that matches a PetscComplex or array of PetscComplex values, pass this value.
4265117d392SLisandro Dalcin 
4275117d392SLisandro Dalcin    Level: beginner
4285117d392SLisandro Dalcin 
429db781477SPatrick Sanan .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_REAL`, `MPIU_SCALAR`, `MPIU_COMPLEX`, `MPIU_INT`, `PETSC_i`
4305117d392SLisandro Dalcin M*/
4315117d392SLisandro Dalcin #if defined(PETSC_USE_REAL_SINGLE)
432de272c7aSSatish Balay #  define MPIU_COMPLEX MPI_C_COMPLEX
4335117d392SLisandro Dalcin #elif defined(PETSC_USE_REAL_DOUBLE)
434de272c7aSSatish Balay #  define MPIU_COMPLEX MPI_C_DOUBLE_COMPLEX
4355117d392SLisandro Dalcin #elif defined(PETSC_USE_REAL___FLOAT128)
4365117d392SLisandro Dalcin #  define MPIU_COMPLEX MPIU___COMPLEX128
4375117d392SLisandro Dalcin #elif defined(PETSC_USE_REAL___FP16)
438de272c7aSSatish Balay #  define MPIU_COMPLEX MPI_C_COMPLEX
4395117d392SLisandro Dalcin #endif /* PETSC_USE_REAL_* */
4405117d392SLisandro Dalcin 
4415117d392SLisandro Dalcin #endif /* PETSC_HAVE_COMPLEX */
4425117d392SLisandro Dalcin 
4435117d392SLisandro Dalcin /*
4445117d392SLisandro Dalcin     Scalar number definitions
4455117d392SLisandro Dalcin  */
4467a19d461SSatish Balay #if defined(PETSC_USE_COMPLEX) && defined(PETSC_HAVE_COMPLEX)
4475117d392SLisandro Dalcin /*MC
4485117d392SLisandro Dalcin    MPIU_SCALAR - MPI datatype corresponding to PetscScalar
4495117d392SLisandro Dalcin 
4505117d392SLisandro Dalcin    Notes:
4515117d392SLisandro Dalcin    In MPI calls that require an MPI datatype that matches a PetscScalar or array of PetscScalar values, pass this value.
4525117d392SLisandro Dalcin 
4535117d392SLisandro Dalcin    Level: beginner
4545117d392SLisandro Dalcin 
455db781477SPatrick Sanan .seealso: `PetscReal`, `PetscScalar`, `PetscComplex`, `PetscInt`, `MPIU_REAL`, `MPIU_COMPLEX`, `MPIU_INT`
4565117d392SLisandro Dalcin M*/
4575117d392SLisandro Dalcin #define MPIU_SCALAR MPIU_COMPLEX
4585117d392SLisandro Dalcin 
4595117d392SLisandro Dalcin /*MC
4605117d392SLisandro Dalcin    PetscRealPart - Returns the real part of a PetscScalar
4615117d392SLisandro Dalcin 
4625117d392SLisandro Dalcin    Synopsis:
4635117d392SLisandro Dalcin    #include <petscmath.h>
4645117d392SLisandro Dalcin    PetscReal PetscRealPart(PetscScalar v)
4655117d392SLisandro Dalcin 
4665117d392SLisandro Dalcin    Not Collective
4675117d392SLisandro Dalcin 
4685117d392SLisandro Dalcin    Input Parameter:
4695117d392SLisandro Dalcin .  v - value to find the real part of
4705117d392SLisandro Dalcin 
4715117d392SLisandro Dalcin    Level: beginner
4725117d392SLisandro Dalcin 
473db781477SPatrick Sanan .seealso: `PetscScalar`, `PetscImaginaryPart()`, `PetscMax()`, `PetscClipInterval()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscSqr()`
4745117d392SLisandro Dalcin 
4755117d392SLisandro Dalcin M*/
4765117d392SLisandro Dalcin #define PetscRealPart(a)      PetscRealPartComplex(a)
4775117d392SLisandro Dalcin 
4785117d392SLisandro Dalcin /*MC
4795117d392SLisandro Dalcin    PetscImaginaryPart - Returns the imaginary part of a PetscScalar
4805117d392SLisandro Dalcin 
4815117d392SLisandro Dalcin    Synopsis:
4825117d392SLisandro Dalcin    #include <petscmath.h>
4835117d392SLisandro Dalcin    PetscReal PetscImaginaryPart(PetscScalar v)
4845117d392SLisandro Dalcin 
4855117d392SLisandro Dalcin    Not Collective
4865117d392SLisandro Dalcin 
4875117d392SLisandro Dalcin    Input Parameter:
4885117d392SLisandro Dalcin .  v - value to find the imaginary part of
4895117d392SLisandro Dalcin 
4905117d392SLisandro Dalcin    Level: beginner
4915117d392SLisandro Dalcin 
4925117d392SLisandro Dalcin    Notes:
4935117d392SLisandro Dalcin        If PETSc was configured for real numbers then this always returns the value 0
4945117d392SLisandro Dalcin 
495db781477SPatrick Sanan .seealso: `PetscScalar`, `PetscRealPart()`, `PetscMax()`, `PetscClipInterval()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscSqr()`
4965117d392SLisandro Dalcin 
4975117d392SLisandro Dalcin M*/
4985117d392SLisandro Dalcin #define PetscImaginaryPart(a) PetscImaginaryPartComplex(a)
4995117d392SLisandro Dalcin 
5005117d392SLisandro Dalcin #define PetscAbsScalar(a)     PetscAbsComplex(a)
5015117d392SLisandro Dalcin #define PetscArgScalar(a)     PetscArgComplex(a)
5025117d392SLisandro Dalcin #define PetscConj(a)          PetscConjComplex(a)
5035117d392SLisandro Dalcin #define PetscSqrtScalar(a)    PetscSqrtComplex(a)
5045117d392SLisandro Dalcin #define PetscPowScalar(a,b)   PetscPowComplex(a,b)
5055117d392SLisandro Dalcin #define PetscExpScalar(a)     PetscExpComplex(a)
5065117d392SLisandro Dalcin #define PetscLogScalar(a)     PetscLogComplex(a)
5075117d392SLisandro Dalcin #define PetscSinScalar(a)     PetscSinComplex(a)
5085117d392SLisandro Dalcin #define PetscCosScalar(a)     PetscCosComplex(a)
5095117d392SLisandro Dalcin #define PetscTanScalar(a)     PetscTanComplex(a)
5105117d392SLisandro Dalcin #define PetscAsinScalar(a)    PetscAsinComplex(a)
5115117d392SLisandro Dalcin #define PetscAcosScalar(a)    PetscAcosComplex(a)
5125117d392SLisandro Dalcin #define PetscAtanScalar(a)    PetscAtanComplex(a)
5135117d392SLisandro Dalcin #define PetscSinhScalar(a)    PetscSinhComplex(a)
5145117d392SLisandro Dalcin #define PetscCoshScalar(a)    PetscCoshComplex(a)
5155117d392SLisandro Dalcin #define PetscTanhScalar(a)    PetscTanhComplex(a)
5165117d392SLisandro Dalcin #define PetscAsinhScalar(a)   PetscAsinhComplex(a)
5175117d392SLisandro Dalcin #define PetscAcoshScalar(a)   PetscAcoshComplex(a)
5185117d392SLisandro Dalcin #define PetscAtanhScalar(a)   PetscAtanhComplex(a)
5195117d392SLisandro Dalcin 
5205117d392SLisandro Dalcin #else /* PETSC_USE_COMPLEX */
5215117d392SLisandro Dalcin #define MPIU_SCALAR MPIU_REAL
5225117d392SLisandro Dalcin #define PetscRealPart(a)      (a)
5235117d392SLisandro Dalcin #define PetscImaginaryPart(a) ((PetscReal)0)
5245117d392SLisandro Dalcin #define PetscAbsScalar(a)     PetscAbsReal(a)
5255117d392SLisandro Dalcin #define PetscArgScalar(a)     (((a) < (PetscReal)0) ? PETSC_PI : (PetscReal)0)
5265117d392SLisandro Dalcin #define PetscConj(a)          (a)
5275117d392SLisandro Dalcin #define PetscSqrtScalar(a)    PetscSqrtReal(a)
5285117d392SLisandro Dalcin #define PetscPowScalar(a,b)   PetscPowReal(a,b)
5295117d392SLisandro Dalcin #define PetscExpScalar(a)     PetscExpReal(a)
5305117d392SLisandro Dalcin #define PetscLogScalar(a)     PetscLogReal(a)
5315117d392SLisandro Dalcin #define PetscSinScalar(a)     PetscSinReal(a)
5325117d392SLisandro Dalcin #define PetscCosScalar(a)     PetscCosReal(a)
5335117d392SLisandro Dalcin #define PetscTanScalar(a)     PetscTanReal(a)
5345117d392SLisandro Dalcin #define PetscAsinScalar(a)    PetscAsinReal(a)
5355117d392SLisandro Dalcin #define PetscAcosScalar(a)    PetscAcosReal(a)
5365117d392SLisandro Dalcin #define PetscAtanScalar(a)    PetscAtanReal(a)
5375117d392SLisandro Dalcin #define PetscSinhScalar(a)    PetscSinhReal(a)
5385117d392SLisandro Dalcin #define PetscCoshScalar(a)    PetscCoshReal(a)
5395117d392SLisandro Dalcin #define PetscTanhScalar(a)    PetscTanhReal(a)
5405117d392SLisandro Dalcin #define PetscAsinhScalar(a)   PetscAsinhReal(a)
5415117d392SLisandro Dalcin #define PetscAcoshScalar(a)   PetscAcoshReal(a)
5425117d392SLisandro Dalcin #define PetscAtanhScalar(a)   PetscAtanhReal(a)
5435117d392SLisandro Dalcin 
5445117d392SLisandro Dalcin #endif /* PETSC_USE_COMPLEX */
5455117d392SLisandro Dalcin 
5465117d392SLisandro Dalcin /*
5475117d392SLisandro Dalcin    Certain objects may be created using either single or double precision.
5485117d392SLisandro Dalcin    This is currently not used.
5495117d392SLisandro Dalcin */
5505117d392SLisandro Dalcin typedef enum { PETSC_SCALAR_DOUBLE, PETSC_SCALAR_SINGLE, PETSC_SCALAR_LONG_DOUBLE, PETSC_SCALAR_HALF } PetscScalarPrecision;
5515117d392SLisandro Dalcin 
5525117d392SLisandro Dalcin /* --------------------------------------------------------------------------*/
5535117d392SLisandro Dalcin 
5545117d392SLisandro Dalcin /*MC
5555117d392SLisandro Dalcin    PetscAbs - Returns the absolute value of a number
5565117d392SLisandro Dalcin 
5575117d392SLisandro Dalcin    Synopsis:
5585117d392SLisandro Dalcin    #include <petscmath.h>
5595117d392SLisandro Dalcin    type PetscAbs(type v)
5605117d392SLisandro Dalcin 
5615117d392SLisandro Dalcin    Not Collective
5625117d392SLisandro Dalcin 
5635117d392SLisandro Dalcin    Input Parameter:
5645117d392SLisandro Dalcin .  v - the number
5655117d392SLisandro Dalcin 
5665117d392SLisandro Dalcin    Notes:
5675117d392SLisandro Dalcin     type can be integer or real floating point value
5685117d392SLisandro Dalcin 
5695117d392SLisandro Dalcin    Level: beginner
5705117d392SLisandro Dalcin 
571db781477SPatrick Sanan .seealso: `PetscAbsInt()`, `PetscAbsReal()`, `PetscAbsScalar()`
5725117d392SLisandro Dalcin 
5735117d392SLisandro Dalcin M*/
5745117d392SLisandro Dalcin #define PetscAbs(a)  (((a) >= 0) ? (a) : (-(a)))
5755117d392SLisandro Dalcin 
5765117d392SLisandro Dalcin /*MC
5775117d392SLisandro Dalcin    PetscSign - Returns the sign of a number as an integer
5785117d392SLisandro Dalcin 
5795117d392SLisandro Dalcin    Synopsis:
5805117d392SLisandro Dalcin    #include <petscmath.h>
5815117d392SLisandro Dalcin    int PetscSign(type v)
5825117d392SLisandro Dalcin 
5835117d392SLisandro Dalcin    Not Collective
5845117d392SLisandro Dalcin 
5855117d392SLisandro Dalcin    Input Parameter:
5865117d392SLisandro Dalcin .  v - the number
5875117d392SLisandro Dalcin 
5885117d392SLisandro Dalcin    Notes:
5895117d392SLisandro Dalcin     type can be integer or real floating point value
5905117d392SLisandro Dalcin 
5915117d392SLisandro Dalcin    Level: beginner
5925117d392SLisandro Dalcin 
5935117d392SLisandro Dalcin M*/
5945117d392SLisandro Dalcin #define PetscSign(a) (((a) >= 0) ? ((a) == 0 ? 0 : 1) : -1)
595e489efc1SBarry Smith 
596b6a5bde7SBarry Smith /*MC
597b6a5bde7SBarry Smith    PetscMin - Returns minimum of two numbers
598b6a5bde7SBarry Smith 
599eca87e8dSBarry Smith    Synopsis:
600aaa7dc30SBarry Smith    #include <petscmath.h>
601eca87e8dSBarry Smith    type PetscMin(type v1,type v2)
602eca87e8dSBarry Smith 
603eca87e8dSBarry Smith    Not Collective
604eca87e8dSBarry Smith 
605d8d19677SJose E. Roman    Input Parameters:
606b6a5bde7SBarry Smith +  v1 - first value to find minimum of
607b6a5bde7SBarry Smith -  v2 - second value to find minimum of
608b6a5bde7SBarry Smith 
60995452b02SPatrick Sanan    Notes:
61095452b02SPatrick Sanan     type can be integer or floating point value
611b6a5bde7SBarry Smith 
612b6a5bde7SBarry Smith    Level: beginner
613b6a5bde7SBarry Smith 
614db781477SPatrick Sanan .seealso: `PetscMax()`, `PetscClipInterval()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscSqr()`
615b6a5bde7SBarry Smith 
616b6a5bde7SBarry Smith M*/
617e489efc1SBarry Smith #define PetscMin(a,b)   (((a)<(b)) ?  (a) : (b))
618b6a5bde7SBarry Smith 
619b6a5bde7SBarry Smith /*MC
620b6a5bde7SBarry Smith    PetscMax - Returns maxium of two numbers
621b6a5bde7SBarry Smith 
622eca87e8dSBarry Smith    Synopsis:
623aaa7dc30SBarry Smith    #include <petscmath.h>
624eca87e8dSBarry Smith    type max PetscMax(type v1,type v2)
625eca87e8dSBarry Smith 
626eca87e8dSBarry Smith    Not Collective
627eca87e8dSBarry Smith 
628d8d19677SJose E. Roman    Input Parameters:
629b6a5bde7SBarry Smith +  v1 - first value to find maximum of
630b6a5bde7SBarry Smith -  v2 - second value to find maximum of
631b6a5bde7SBarry Smith 
63295452b02SPatrick Sanan    Notes:
63395452b02SPatrick Sanan     type can be integer or floating point value
634b6a5bde7SBarry Smith 
635b6a5bde7SBarry Smith    Level: beginner
636b6a5bde7SBarry Smith 
637db781477SPatrick Sanan .seealso: `PetscMin()`, `PetscClipInterval()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscSqr()`
638b6a5bde7SBarry Smith 
639b6a5bde7SBarry Smith M*/
640e489efc1SBarry Smith #define PetscMax(a,b)   (((a)<(b)) ?  (b) : (a))
641b6a5bde7SBarry Smith 
642b6a5bde7SBarry Smith /*MC
643d9a4bb16SJed Brown    PetscClipInterval - Returns a number clipped to be within an interval
644d9a4bb16SJed Brown 
645d9a4bb16SJed Brown    Synopsis:
646aaa7dc30SBarry Smith    #include <petscmath.h>
647d9a4bb16SJed Brown    type clip PetscClipInterval(type x,type a,type b)
648d9a4bb16SJed Brown 
649d9a4bb16SJed Brown    Not Collective
650d9a4bb16SJed Brown 
651d8d19677SJose E. Roman    Input Parameters:
6520d398bfeSStefano Zampini +  x - value to use if within interval [a,b]
653d9a4bb16SJed Brown .  a - lower end of interval
654d9a4bb16SJed Brown -  b - upper end of interval
655d9a4bb16SJed Brown 
65695452b02SPatrick Sanan    Notes:
65795452b02SPatrick Sanan     type can be integer or floating point value
658d9a4bb16SJed Brown 
659d9a4bb16SJed Brown    Level: beginner
660d9a4bb16SJed Brown 
661db781477SPatrick Sanan .seealso: `PetscMin()`, `PetscMax()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscSqr()`
662d9a4bb16SJed Brown 
663d9a4bb16SJed Brown M*/
664d9a4bb16SJed Brown #define PetscClipInterval(x,a,b)   (PetscMax((a),PetscMin((x),(b))))
665d9a4bb16SJed Brown 
666d9a4bb16SJed Brown /*MC
667b6a5bde7SBarry Smith    PetscAbsInt - Returns the absolute value of an integer
668b6a5bde7SBarry Smith 
669b6a5bde7SBarry Smith    Synopsis:
670aaa7dc30SBarry Smith    #include <petscmath.h>
671b6a5bde7SBarry Smith    int abs PetscAbsInt(int v1)
672b6a5bde7SBarry Smith 
673eca87e8dSBarry Smith    Not Collective
674eca87e8dSBarry Smith 
675eca87e8dSBarry Smith    Input Parameter:
676eca87e8dSBarry Smith .   v1 - the integer
677b6a5bde7SBarry Smith 
678b6a5bde7SBarry Smith    Level: beginner
679b6a5bde7SBarry Smith 
680db781477SPatrick Sanan .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsReal()`, `PetscSqr()`
681b6a5bde7SBarry Smith 
682b6a5bde7SBarry Smith M*/
6839fa7d148SSatish Balay #define PetscAbsInt(a)  (((a)<0)   ? (-(a)) : (a))
684b6a5bde7SBarry Smith 
685b6a5bde7SBarry Smith /*MC
686b6a5bde7SBarry Smith    PetscAbsReal - Returns the absolute value of an real number
687b6a5bde7SBarry Smith 
688eca87e8dSBarry Smith    Synopsis:
689aaa7dc30SBarry Smith    #include <petscmath.h>
690eca87e8dSBarry Smith    Real abs PetscAbsReal(PetscReal v1)
691eca87e8dSBarry Smith 
692eca87e8dSBarry Smith    Not Collective
693eca87e8dSBarry Smith 
694b6a5bde7SBarry Smith    Input Parameter:
695b6a5bde7SBarry Smith .   v1 - the double
696b6a5bde7SBarry Smith 
697b6a5bde7SBarry Smith    Level: beginner
698b6a5bde7SBarry Smith 
699db781477SPatrick Sanan .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsInt()`, `PetscSqr()`
700b6a5bde7SBarry Smith 
701b6a5bde7SBarry Smith M*/
7021118d4bcSLisandro Dalcin #if defined(PETSC_USE_REAL_SINGLE)
7031118d4bcSLisandro Dalcin #define PetscAbsReal(a) fabsf(a)
7041118d4bcSLisandro Dalcin #elif defined(PETSC_USE_REAL_DOUBLE)
7051118d4bcSLisandro Dalcin #define PetscAbsReal(a) fabs(a)
7061118d4bcSLisandro Dalcin #elif defined(PETSC_USE_REAL___FLOAT128)
7071118d4bcSLisandro Dalcin #define PetscAbsReal(a) fabsq(a)
7081118d4bcSLisandro Dalcin #elif defined(PETSC_USE_REAL___FP16)
7091118d4bcSLisandro Dalcin #define PetscAbsReal(a) fabsf(a)
7101118d4bcSLisandro Dalcin #endif
711b6a5bde7SBarry Smith 
712b6a5bde7SBarry Smith /*MC
713b6a5bde7SBarry Smith    PetscSqr - Returns the square of a number
714b6a5bde7SBarry Smith 
715b6a5bde7SBarry Smith    Synopsis:
716aaa7dc30SBarry Smith    #include <petscmath.h>
717b6a5bde7SBarry Smith    type sqr PetscSqr(type v1)
718b6a5bde7SBarry Smith 
719eca87e8dSBarry Smith    Not Collective
720eca87e8dSBarry Smith 
721eca87e8dSBarry Smith    Input Parameter:
722eca87e8dSBarry Smith .   v1 - the value
723eca87e8dSBarry Smith 
72495452b02SPatrick Sanan    Notes:
72595452b02SPatrick Sanan     type can be integer or floating point value
726b6a5bde7SBarry Smith 
727b6a5bde7SBarry Smith    Level: beginner
728b6a5bde7SBarry Smith 
729db781477SPatrick Sanan .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsInt()`, `PetscAbsReal()`
730b6a5bde7SBarry Smith 
731b6a5bde7SBarry Smith M*/
7324ebda54eSMatthew Knepley #define PetscSqr(a)     ((a)*(a))
733e489efc1SBarry Smith 
734314da920SBarry Smith /* ----------------------------------------------------------------------------*/
735ee223c85SLisandro Dalcin 
736ee223c85SLisandro Dalcin #if defined(PETSC_USE_REAL_SINGLE)
737ee223c85SLisandro Dalcin #define PetscRealConstant(constant) constant##F
7385117d392SLisandro Dalcin #elif defined(PETSC_USE_REAL_DOUBLE)
7395117d392SLisandro Dalcin #define PetscRealConstant(constant) constant
740ee223c85SLisandro Dalcin #elif defined(PETSC_USE_REAL___FLOAT128)
741ee223c85SLisandro Dalcin #define PetscRealConstant(constant) constant##Q
7425117d392SLisandro Dalcin #elif defined(PETSC_USE_REAL___FP16)
7435117d392SLisandro Dalcin #define PetscRealConstant(constant) constant##F
744ee223c85SLisandro Dalcin #endif
745ee223c85SLisandro Dalcin 
746314da920SBarry Smith /*
747d34fcf5fSBarry Smith      Basic constants
748314da920SBarry Smith */
7492fab75feSLisandro Dalcin #define PETSC_PI    PetscRealConstant(3.1415926535897932384626433832795029)
7502fab75feSLisandro Dalcin #define PETSC_PHI   PetscRealConstant(1.6180339887498948482045868343656381)
7517b156302SMatthew G. Knepley #define PETSC_SQRT2 PetscRealConstant(1.4142135623730950488016887242096981)
752d34fcf5fSBarry Smith 
753ab824b78SBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
75471fd2e92SBarry Smith #define PETSC_MAX_INT            2147483647
755ab824b78SBarry Smith #define PETSC_MIN_INT            (-PETSC_MAX_INT - 1)
756ab824b78SBarry Smith #else
757ab824b78SBarry Smith #define PETSC_MAX_INT            9223372036854775807L
758ab824b78SBarry Smith #define PETSC_MIN_INT            (-PETSC_MAX_INT - 1)
759ab824b78SBarry Smith #endif
760569ea7c4SPierre Jolivet #define PETSC_MAX_UINT16         65535
761e489efc1SBarry Smith 
762ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
763ab824b78SBarry Smith #  define PETSC_MAX_REAL                3.40282346638528860e+38F
7649fa7d148SSatish Balay #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
76582a7e548SBarry Smith #  define PETSC_MACHINE_EPSILON         1.19209290e-07F
76682a7e548SBarry Smith #  define PETSC_SQRT_MACHINE_EPSILON    3.45266983e-04F
767ee223c85SLisandro Dalcin #  define PETSC_SMALL                   1.e-5F
768ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL_DOUBLE)
769ab824b78SBarry Smith #  define PETSC_MAX_REAL                1.7976931348623157e+308
7709fa7d148SSatish Balay #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
77182a7e548SBarry Smith #  define PETSC_MACHINE_EPSILON         2.2204460492503131e-16
77282a7e548SBarry Smith #  define PETSC_SQRT_MACHINE_EPSILON    1.490116119384766e-08
773cf6e855fSSatish Balay #  define PETSC_SMALL                   1.e-10
774ce63c4c1SBarry Smith #elif defined(PETSC_USE_REAL___FLOAT128)
775ea345e14SBarry Smith #  define PETSC_MAX_REAL                FLT128_MAX
7769fa7d148SSatish Balay #  define PETSC_MIN_REAL                (-FLT128_MAX)
777d34fcf5fSBarry Smith #  define PETSC_MACHINE_EPSILON         FLT128_EPSILON
778ee223c85SLisandro Dalcin #  define PETSC_SQRT_MACHINE_EPSILON    1.38777878078144567552953958511352539e-17Q
779ee223c85SLisandro Dalcin #  define PETSC_SMALL                   1.e-20Q
7805117d392SLisandro Dalcin #elif defined(PETSC_USE_REAL___FP16)
7815117d392SLisandro Dalcin #  define PETSC_MAX_REAL                65504.0F
7829fa7d148SSatish Balay #  define PETSC_MIN_REAL                (-PETSC_MAX_REAL)
7835117d392SLisandro Dalcin #  define PETSC_MACHINE_EPSILON         .0009765625F
7845117d392SLisandro Dalcin #  define PETSC_SQRT_MACHINE_EPSILON    .03125F
7855117d392SLisandro Dalcin #  define PETSC_SMALL                   5.e-3F
7869cf09972SJed Brown #endif
7873e523bebSBarry Smith 
78825d0f998SSatish Balay #define PETSC_INFINITY               (PETSC_MAX_REAL/4)
7899fa7d148SSatish Balay #define PETSC_NINFINITY              (-PETSC_INFINITY)
790e270355aSBarry Smith 
7919f4f8022SLisandro Dalcin PETSC_EXTERN PetscBool PetscIsInfReal(PetscReal);
7923948c36eSLisandro Dalcin PETSC_EXTERN PetscBool PetscIsNanReal(PetscReal);
7938b49ba18SBarry Smith PETSC_EXTERN PetscBool PetscIsNormalReal(PetscReal);
7949fbee547SJacob Faibussowitsch static inline PetscBool PetscIsInfOrNanReal(PetscReal v) {return PetscIsInfReal(v) || PetscIsNanReal(v) ? PETSC_TRUE : PETSC_FALSE;}
7959fbee547SJacob Faibussowitsch static inline PetscBool PetscIsInfScalar(PetscScalar v) {return PetscIsInfReal(PetscAbsScalar(v));}
7969fbee547SJacob Faibussowitsch static inline PetscBool PetscIsNanScalar(PetscScalar v) {return PetscIsNanReal(PetscAbsScalar(v));}
7979fbee547SJacob Faibussowitsch static inline PetscBool PetscIsInfOrNanScalar(PetscScalar v) {return PetscIsInfOrNanReal(PetscAbsScalar(v));}
7989fbee547SJacob Faibussowitsch static inline PetscBool PetscIsNormalScalar(PetscScalar v) {return PetscIsNormalReal(PetscAbsScalar(v));}
7999a25a3ccSBarry Smith 
800b10005b4SLisandro Dalcin PETSC_EXTERN PetscBool PetscIsCloseAtTol(PetscReal,PetscReal,PetscReal,PetscReal);
801ce4818fdSLisandro Dalcin PETSC_EXTERN PetscBool PetscEqualReal(PetscReal,PetscReal);
802ce4818fdSLisandro Dalcin PETSC_EXTERN PetscBool PetscEqualScalar(PetscScalar,PetscScalar);
803ce4818fdSLisandro Dalcin 
80498725619SBarry Smith /*
80598725619SBarry Smith     These macros are currently hardwired to match the regular data types, so there is no support for a different
80698725619SBarry Smith     MatScalar from PetscScalar. We left the MatScalar in the source just in case we use it again.
80798725619SBarry Smith  */
80898725619SBarry Smith #define MPIU_MATSCALAR MPIU_SCALAR
80998725619SBarry Smith typedef PetscScalar MatScalar;
81098725619SBarry Smith typedef PetscReal MatReal;
81198725619SBarry Smith 
8128ad47952SJed Brown struct petsc_mpiu_2scalar {PetscScalar a,b;};
8138ad47952SJed Brown PETSC_EXTERN MPI_Datatype MPIU_2SCALAR PetscAttrMPITypeTagLayoutCompatible(struct petsc_mpiu_2scalar);
814df4397b0SStefano Zampini 
815092991acSStefano Zampini /*
816092991acSStefano Zampini    MPI Datatypes for composite reductions:
817092991acSStefano Zampini    MPIU_REAL_INT -> struct { PetscReal; PetscInt; }
818092991acSStefano Zampini    MPIU_SCALAR_INT -> struct { PetscScalar; PetscInt; }
819092991acSStefano Zampini */
820092991acSStefano Zampini PETSC_EXTERN MPI_Datatype MPIU_REAL_INT;
821092991acSStefano Zampini PETSC_EXTERN MPI_Datatype MPIU_SCALAR_INT;
822092991acSStefano Zampini 
823a616ada9SVaclav Hapla #if defined(PETSC_USE_64BIT_INDICES)
8248ad47952SJed Brown struct petsc_mpiu_2int {PetscInt a,b;};
8258ad47952SJed Brown PETSC_EXTERN MPI_Datatype MPIU_2INT PetscAttrMPITypeTagLayoutCompatible(struct petsc_mpiu_2int);
8268ad47952SJed Brown #else
8278ad47952SJed Brown #define MPIU_2INT MPI_2INT
8288ad47952SJed Brown #endif
829b5a892a1SMatthew G. Knepley PETSC_EXTERN MPI_Datatype MPI_4INT;
830b5a892a1SMatthew G. Knepley PETSC_EXTERN MPI_Datatype MPIU_4INT;
831e9fa29b7SSatish Balay 
8329fbee547SJacob Faibussowitsch static inline PetscInt PetscPowInt(PetscInt base,PetscInt power)
833b2fb0278SBarry Smith {
834fa711258SJed Brown   PetscInt result = 1;
835fa711258SJed Brown   while (power) {
836fa711258SJed Brown     if (power & 1) result *= base;
837fa711258SJed Brown     power >>= 1;
838fa711258SJed Brown     base *= base;
839fa711258SJed Brown   }
840fa711258SJed Brown   return result;
841fa711258SJed Brown }
842b2fb0278SBarry Smith 
8439fbee547SJacob Faibussowitsch static inline PetscInt64 PetscPowInt64(PetscInt base,PetscInt power)
844ad70a4c3SStefano Zampini {
845ad70a4c3SStefano Zampini   PetscInt64 result = 1;
846ad70a4c3SStefano Zampini   while (power) {
847ad70a4c3SStefano Zampini     if (power & 1) result *= base;
848ad70a4c3SStefano Zampini     power >>= 1;
849ad70a4c3SStefano Zampini     base *= base;
850ad70a4c3SStefano Zampini   }
851ad70a4c3SStefano Zampini   return result;
852ad70a4c3SStefano Zampini }
853ad70a4c3SStefano Zampini 
8549fbee547SJacob Faibussowitsch static inline PetscReal PetscPowRealInt(PetscReal base,PetscInt power)
855b2fb0278SBarry Smith {
856fa711258SJed Brown   PetscReal result = 1;
857d98d5da7SBarry Smith   if (power < 0) {
858d98d5da7SBarry Smith     power = -power;
85910d40e53SLisandro Dalcin     base  = ((PetscReal)1)/base;
860d98d5da7SBarry Smith   }
861fa711258SJed Brown   while (power) {
862fa711258SJed Brown     if (power & 1) result *= base;
863fa711258SJed Brown     power >>= 1;
864fa711258SJed Brown     base *= base;
865fa711258SJed Brown   }
866fa711258SJed Brown   return result;
867fa711258SJed Brown }
868fa711258SJed Brown 
8699fbee547SJacob Faibussowitsch static inline PetscScalar PetscPowScalarInt(PetscScalar base,PetscInt power)
870b2fb0278SBarry Smith {
8715117d392SLisandro Dalcin   PetscScalar result = (PetscReal)1;
8728b49ba18SBarry Smith   if (power < 0) {
8738b49ba18SBarry Smith     power = -power;
87410d40e53SLisandro Dalcin     base  = ((PetscReal)1)/base;
8758b49ba18SBarry Smith   }
8768b49ba18SBarry Smith   while (power) {
8778b49ba18SBarry Smith     if (power & 1) result *= base;
8788b49ba18SBarry Smith     power >>= 1;
8798b49ba18SBarry Smith     base *= base;
8808b49ba18SBarry Smith   }
8818b49ba18SBarry Smith   return result;
8828b49ba18SBarry Smith }
8838b49ba18SBarry Smith 
8849fbee547SJacob Faibussowitsch static inline PetscScalar PetscPowScalarReal(PetscScalar base,PetscReal power)
885b2fb0278SBarry Smith {
886b2fb0278SBarry Smith   PetscScalar cpower = power;
887b2fb0278SBarry Smith   return PetscPowScalar(base,cpower);
888b2fb0278SBarry Smith }
88978a59e97SMatthew G. Knepley 
890c803a25aSBarry Smith /*MC
89166baab88SBarry Smith     PetscApproximateLTE - Performs a less than or equal to on a given constant with a fudge for floating point numbers
892c803a25aSBarry Smith 
893c803a25aSBarry Smith    Synopsis:
894c803a25aSBarry Smith    #include <petscmath.h>
89566baab88SBarry Smith    bool PetscApproximateLTE(PetscReal x,constant float)
896c803a25aSBarry Smith 
897c803a25aSBarry Smith    Not Collective
898c803a25aSBarry Smith 
899c803a25aSBarry Smith    Input Parameters:
900c803a25aSBarry Smith +   x - the variable
901c803a25aSBarry Smith -   b - the constant float it is checking if x is less than or equal to
902c803a25aSBarry Smith 
903c803a25aSBarry Smith    Notes:
904c803a25aSBarry Smith      The fudge factor is the value PETSC_SMALL
905c803a25aSBarry Smith 
906c803a25aSBarry Smith      The constant numerical value is automatically set to the appropriate precision of PETSc so can just be provided as, for example, 3.2
907c803a25aSBarry Smith 
908c803a25aSBarry Smith      This is used in several examples for setting initial conditions based on coordinate values that are computed with i*h that produces inexact
909c803a25aSBarry Smith      floating point results.
910c803a25aSBarry Smith 
911c803a25aSBarry Smith    Level: advanced
912c803a25aSBarry Smith 
913db781477SPatrick Sanan .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscApproximateGTE()`
914c803a25aSBarry Smith 
915c803a25aSBarry Smith M*/
91666baab88SBarry Smith #define PetscApproximateLTE(x,b)  ((x) <= (PetscRealConstant(b)+PETSC_SMALL))
917c803a25aSBarry Smith 
918c803a25aSBarry Smith /*MC
91966baab88SBarry Smith     PetscApproximateGTE - Performs a greater than or equal to on a given constant with a fudge for floating point numbers
920c803a25aSBarry Smith 
921c803a25aSBarry Smith    Synopsis:
922c803a25aSBarry Smith    #include <petscmath.h>
92366baab88SBarry Smith    bool PetscApproximateGTE(PetscReal x,constant float)
924c803a25aSBarry Smith 
925c803a25aSBarry Smith    Not Collective
926c803a25aSBarry Smith 
927c803a25aSBarry Smith    Input Parameters:
928c803a25aSBarry Smith +   x - the variable
929c803a25aSBarry Smith -   b - the constant float it is checking if x is greater than or equal to
930c803a25aSBarry Smith 
931c803a25aSBarry Smith    Notes:
932c803a25aSBarry Smith      The fudge factor is the value PETSC_SMALL
933c803a25aSBarry Smith 
934c803a25aSBarry Smith      The constant numerical value is automatically set to the appropriate precision of PETSc so can just be provided as, for example, 3.2
935c803a25aSBarry Smith 
936c803a25aSBarry Smith      This is used in several examples for setting initial conditions based on coordinate values that are computed with i*h that produces inexact
937c803a25aSBarry Smith      floating point results.
938c803a25aSBarry Smith 
939c803a25aSBarry Smith    Level: advanced
940c803a25aSBarry Smith 
941db781477SPatrick Sanan .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscApproximateLTE()`
942c803a25aSBarry Smith 
943c803a25aSBarry Smith M*/
94466baab88SBarry Smith #define PetscApproximateGTE(x,b)  ((x) >= (PetscRealConstant(b)-PETSC_SMALL))
945c803a25aSBarry Smith 
946faa75363SBarry Smith /*MC
947faa75363SBarry Smith     PetscCeilInt - Returns the ceiling of the quotation of two positive integers
948faa75363SBarry Smith 
949faa75363SBarry Smith    Synopsis:
950faa75363SBarry Smith    #include <petscmath.h>
951faa75363SBarry Smith    PetscInt PetscCeilInt(PetscInt x,PetscInt y)
952faa75363SBarry Smith 
953faa75363SBarry Smith    Not Collective
954faa75363SBarry Smith 
955faa75363SBarry Smith    Input Parameters:
956faa75363SBarry Smith +   x - the numerator
957faa75363SBarry Smith -   y - the denominator
958faa75363SBarry Smith 
959faa75363SBarry Smith    Level: advanced
960faa75363SBarry Smith 
961db781477SPatrick Sanan .seealso: `PetscMax()`, `PetscMin()`, `PetscAbsInt()`, `PetscAbsReal()`, `PetscApproximateLTE()`
962faa75363SBarry Smith 
963faa75363SBarry Smith M*/
964faa75363SBarry Smith #define PetscCeilInt(x,y)  ((((PetscInt)(x))/((PetscInt)(y))) +  ((((PetscInt)(x)) % ((PetscInt)(y))) ? 1 : 0))
965faa75363SBarry Smith 
966faa75363SBarry Smith #define PetscCeilInt64(x,y)  ((((PetscInt64)(x))/((PetscInt64)(y))) +  ((((PetscInt64)(x)) % ((PetscInt64)(y))) ? 1 : 0))
967faa75363SBarry Smith 
968bebf13c0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLinearRegression(PetscInt,const PetscReal[],const PetscReal[],PetscReal*,PetscReal*);
969e489efc1SBarry Smith #endif
970