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,y - the real numbers
 46:    Output parameter:
 47: .  return - the result
 49:    Fortran Note:
 50:    This function is not available from Fortran.
 52:    Level: developer
 53: @*/
 54: static inline PetscReal SlepcAbs(PetscReal x,PetscReal y)
 55: {
 56:   PetscReal w,z,t,xabs=PetscAbs(x),yabs=PetscAbs(y);
 58:   w = PetscMax(xabs,yabs);
 59:   z = PetscMin(xabs,yabs);
 60:   if (PetscUnlikely(z == PetscRealConstant(0.0))) return w;
 61:   t = z/w;
 62:   return w*PetscSqrtReal(PetscRealConstant(1.0)+t*t);
 63: }
 65: /*MC
 66:    SlepcAbsEigenvalue - Returns the absolute value of a complex number given
 67:    its real and imaginary parts.
 69:    Synopsis:
 70:    PetscReal SlepcAbsEigenvalue(PetscScalar x,PetscScalar y)
 72:    Not Collective
 74:    Input parameters:
 75: +  x  - the real part of the complex number
 76: -  y  - the imaginary part of the complex number
 78:    Notes:
 79:    This function computes $\sqrt{x^2+y^2}$, taking care not to cause unnecessary
 80:    overflow. It is based on LAPACK's `DLAPY2`.
 82:    In complex scalars, only the first argument is used, i.e., the result is $|x|$.
 84:    Fortran Note:
 85:    This function is not available from Fortran.
 87:    Level: developer
 88: M*/
 89: #if !defined(PETSC_USE_COMPLEX)
 90: #define SlepcAbsEigenvalue(x,y) SlepcAbs(x,y)
 91: #else
 92: #define SlepcAbsEigenvalue(x,y) PetscAbsScalar(x)
 93: #endif
 95: /*
 96:    SlepcSetFlushToZero - Set the FTZ flag in floating-point arithmetic.
 97: */
 98: static inline PetscErrorCode SlepcSetFlushToZero(unsigned int *state)
 99: {
100:   PetscFunctionBegin;
101: #if defined(PETSC_HAVE_XMMINTRIN_H) && defined(_MM_FLUSH_ZERO_ON) && defined(__SSE__)
102:   *state = _MM_GET_FLUSH_ZERO_MODE();
103:   _MM_SET_FLUSH_ZERO_MODE(_MM_FLUSH_ZERO_ON);
104: #else
105:   *state = 0;
106: #endif
107:   PetscFunctionReturn(PETSC_SUCCESS);
108: }
110: /*
111:    SlepcResetFlushToZero - Reset the FTZ flag in floating-point arithmetic.
112: */
113: static inline PetscErrorCode SlepcResetFlushToZero(unsigned int *state)
114: {
115:   PetscFunctionBegin;
116: #if defined(PETSC_HAVE_XMMINTRIN_H) && defined(_MM_FLUSH_ZERO_MASK) && defined(__SSE__)
117:   _MM_SET_FLUSH_ZERO_MODE(*state & _MM_FLUSH_ZERO_MASK);
118: #else
119:   *state = 0;
120: #endif
121:   PetscFunctionReturn(PETSC_SUCCESS);
122: }