Actual source code: slepcmath.h

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc mathematics include file. Defines basic operations and functions.
 12:    This file is included by slepcsys.h and should not be used directly.
 13: */

 15: #pragma once

 17: /* SUBMANSEC = Sys */

 19: /*
 20:     Default tolerance for the different solvers, depending on the precision
 21: */
 22: #if defined(PETSC_USE_REAL_SINGLE)
 23: #  define SLEPC_DEFAULT_TOL   1e-5
 24: #elif defined(PETSC_USE_REAL_DOUBLE)
 25: #  define SLEPC_DEFAULT_TOL   1e-8
 26: #elif defined(PETSC_USE_REAL___FLOAT128)
 27: #  define SLEPC_DEFAULT_TOL   1e-16
 28: #elif defined(PETSC_USE_REAL___FP16)
 29: #  define SLEPC_DEFAULT_TOL   1e-2
 30: #endif

 32: static inline PetscReal SlepcDefaultTol(PetscReal tol)
 33: {
 34:   return tol == (PetscReal)PETSC_DETERMINE ? SLEPC_DEFAULT_TOL : tol;
 35: }

 37: /*@C
 38:    SlepcAbs - Returns $\sqrt{x^2+y^2}$, taking care not to cause unnecessary
 39:    overflow. It is based on LAPACK's `DLAPY2`.

 41:    Not Collective

 43:    Input parameters:
 44: +  x - the first real number
 45: -  y - the second real number

 47:    Output parameter:
 48: .  return - the result

 50:    Fortran Note:
 51:    This function is not available from Fortran.

 53:    Level: developer
 54: @*/
 55: static inline PetscReal SlepcAbs(PetscReal x,PetscReal y)
 56: {
 57:   PetscReal w,z,t,xabs=PetscAbs(x),yabs=PetscAbs(y);

 59:   w = PetscMax(xabs,yabs);
 60:   z = PetscMin(xabs,yabs);
 61:   if (PetscUnlikely(z == PetscRealConstant(0.0))) return w;
 62:   t = z/w;
 63:   return w*PetscSqrtReal(PetscRealConstant(1.0)+t*t);
 64: }

 66: /*MC
 67:    SlepcAbsEigenvalue - Returns the absolute value of a complex number given
 68:    its real and imaginary parts.

 70:    Synopsis:
 71:    PetscReal SlepcAbsEigenvalue(PetscScalar x,PetscScalar y)

 73:    Not Collective

 75:    Input parameters:
 76: +  x - the real part of the complex number
 77: -  y - the imaginary part of the complex number

 79:    Notes:
 80:    This function computes $\sqrt{x^2+y^2}$, taking care not to cause unnecessary
 81:    overflow. It is based on LAPACK's `DLAPY2`.

 83:    In complex scalars, only the first argument is used, i.e., the result is $|x|$.

 85:    Fortran Note:
 86:    This function is not available from Fortran.

 88:    Level: developer

 90: .seealso: `PetscAbsScalar()`
 91: M*/
 92: #if !defined(PETSC_USE_COMPLEX)
 93: #define SlepcAbsEigenvalue(x,y) SlepcAbs(x,y)
 94: #else
 95: #define SlepcAbsEigenvalue(x,y) PetscAbsScalar(x)
 96: #endif

 98: /*
 99:    SlepcSetFlushToZero - Set the FTZ flag in floating-point arithmetic.
100: */
101: static inline PetscErrorCode SlepcSetFlushToZero(unsigned int *state)
102: {
103:   PetscFunctionBegin;
104: #if defined(PETSC_HAVE_XMMINTRIN_H) && defined(_MM_FLUSH_ZERO_ON) && defined(__SSE__)
105:   *state = _MM_GET_FLUSH_ZERO_MODE();
106:   _MM_SET_FLUSH_ZERO_MODE(_MM_FLUSH_ZERO_ON);
107: #else
108:   *state = 0;
109: #endif
110:   PetscFunctionReturn(PETSC_SUCCESS);
111: }

113: /*
114:    SlepcResetFlushToZero - Reset the FTZ flag in floating-point arithmetic.
115: */
116: static inline PetscErrorCode SlepcResetFlushToZero(unsigned int *state)
117: {
118:   PetscFunctionBegin;
119: #if defined(PETSC_HAVE_XMMINTRIN_H) && defined(_MM_FLUSH_ZERO_MASK) && defined(__SSE__)
120:   _MM_SET_FLUSH_ZERO_MODE(*state & _MM_FLUSH_ZERO_MASK);
121: #else
122:   *state = 0;
123: #endif
124:   PetscFunctionReturn(PETSC_SUCCESS);
125: }