Actual source code: slepcmath.h

slepc-3.17.0 2022-03-31
Report Typos and Errors
  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: #if !defined(SLEPCMATH_H)
 16: #define SLEPCMATH_H

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

 31: #define SlepcDefaultTol(tol) (((tol)==PETSC_DEFAULT)?SLEPC_DEFAULT_TOL:(tol))

 33: /*@C
 34:    SlepcAbs - Returns sqrt(x**2+y**2), taking care not to cause unnecessary
 35:    overflow. It is based on LAPACK's DLAPY2.

 37:    Not Collective

 39:    Input parameters:
 40: .  x,y - the real numbers

 42:    Output parameter:
 43: .  return - the result

 45:    Note:
 46:    This function is not available from Fortran.

 48:    Level: developer
 49: @*/
 50: static inline PetscReal SlepcAbs(PetscReal x,PetscReal y)
 51: {
 52:   PetscReal w,z,t,xabs=PetscAbs(x),yabs=PetscAbs(y);

 54:   w = PetscMax(xabs,yabs);
 55:   z = PetscMin(xabs,yabs);
 56:   if (PetscUnlikely(z == 0.0)) return w;
 57:   t = z/w;
 58:   return w*PetscSqrtReal(1.0+t*t);
 59: }

 61: /*MC
 62:    SlepcAbsEigenvalue - Returns the absolute value of a complex number given
 63:    its real and imaginary parts.

 65:    Synopsis:
 66:    PetscReal SlepcAbsEigenvalue(PetscScalar x,PetscScalar y)

 68:    Not Collective

 70:    Input parameters:
 71: +  x  - the real part of the complex number
 72: -  y  - the imaginary part of the complex number

 74:    Notes:
 75:    This function computes sqrt(x**2+y**2), taking care not to cause unnecessary
 76:    overflow. It is based on LAPACK's DLAPY2.

 78:    In complex scalars, only the first argument is used.

 80:    This function is not available from Fortran.

 82:    Level: developer
 83: M*/
 84: #if !defined(PETSC_USE_COMPLEX)
 85: #define SlepcAbsEigenvalue(x,y) SlepcAbs(x,y)
 86: #else
 87: #define SlepcAbsEigenvalue(x,y) PetscAbsScalar(x)
 88: #endif

 90: #endif

 92: /*
 93:    SlepcSetFlushToZero - Set the FTZ flag in floating-point arithmetic.
 94: */
 95: static inline PetscErrorCode SlepcSetFlushToZero(unsigned int *state)
 96: {
 97: #if defined(PETSC_HAVE_XMMINTRIN_H) && defined(_MM_FLUSH_ZERO_ON) && defined(__SSE__)
 98:   *state = _MM_GET_FLUSH_ZERO_MODE();
 99:   _MM_SET_FLUSH_ZERO_MODE(_MM_FLUSH_ZERO_ON);
100: #else
101:   *state = 0;
102: #endif
103:   PetscFunctionReturn(0);
104: }

106: /*
107:    SlepcResetFlushToZero - Reset the FTZ flag in floating-point arithmetic.
108: */
109: static inline PetscErrorCode SlepcResetFlushToZero(unsigned int *state)
110: {
111: #if defined(PETSC_HAVE_XMMINTRIN_H) && defined(_MM_FLUSH_ZERO_MASK) && defined(__SSE__)
112:   _MM_SET_FLUSH_ZERO_MODE(*state & _MM_FLUSH_ZERO_MASK);
113: #else
114:   *state = 0;
115: #endif
116:   PetscFunctionReturn(0);
117: }