1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122
|
/*
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
SLEPc - Scalable Library for Eigenvalue Problem Computations
Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
This file is part of SLEPc.
SLEPc is distributed under a 2-clause BSD license (see LICENSE).
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
*/
/*
SLEPc mathematics include file. Defines basic operations and functions.
This file is included by slepcsys.h and should not be used directly.
*/
#pragma once
/* SUBMANSEC = Sys */
/*
Default tolerance for the different solvers, depending on the precision
*/
#if defined(PETSC_USE_REAL_SINGLE)
# define SLEPC_DEFAULT_TOL 1e-5
#elif defined(PETSC_USE_REAL_DOUBLE)
# define SLEPC_DEFAULT_TOL 1e-8
#elif defined(PETSC_USE_REAL___FLOAT128)
# define SLEPC_DEFAULT_TOL 1e-16
#elif defined(PETSC_USE_REAL___FP16)
# define SLEPC_DEFAULT_TOL 1e-2
#endif
static inline PetscReal SlepcDefaultTol(PetscReal tol)
{
return tol == (PetscReal)PETSC_DETERMINE ? SLEPC_DEFAULT_TOL : tol;
}
/*@C
SlepcAbs - Returns $\sqrt{x^2+y^2}$, taking care not to cause unnecessary
overflow. It is based on LAPACK's `DLAPY2`.
Not Collective
Input parameters:
. x,y - the real numbers
Output parameter:
. return - the result
Fortran Note:
This function is not available from Fortran.
Level: developer
@*/
static inline PetscReal SlepcAbs(PetscReal x,PetscReal y)
{
PetscReal w,z,t,xabs=PetscAbs(x),yabs=PetscAbs(y);
w = PetscMax(xabs,yabs);
z = PetscMin(xabs,yabs);
if (PetscUnlikely(z == PetscRealConstant(0.0))) return w;
t = z/w;
return w*PetscSqrtReal(PetscRealConstant(1.0)+t*t);
}
/*MC
SlepcAbsEigenvalue - Returns the absolute value of a complex number given
its real and imaginary parts.
Synopsis:
PetscReal SlepcAbsEigenvalue(PetscScalar x,PetscScalar y)
Not Collective
Input parameters:
+ x - the real part of the complex number
- y - the imaginary part of the complex number
Notes:
This function computes $\sqrt{x^2+y^2}$, taking care not to cause unnecessary
overflow. It is based on LAPACK's `DLAPY2`.
In complex scalars, only the first argument is used, i.e., the result is $|x|$.
Fortran Note:
This function is not available from Fortran.
Level: developer
M*/
#if !defined(PETSC_USE_COMPLEX)
#define SlepcAbsEigenvalue(x,y) SlepcAbs(x,y)
#else
#define SlepcAbsEigenvalue(x,y) PetscAbsScalar(x)
#endif
/*
SlepcSetFlushToZero - Set the FTZ flag in floating-point arithmetic.
*/
static inline PetscErrorCode SlepcSetFlushToZero(unsigned int *state)
{
PetscFunctionBegin;
#if defined(PETSC_HAVE_XMMINTRIN_H) && defined(_MM_FLUSH_ZERO_ON) && defined(__SSE__)
*state = _MM_GET_FLUSH_ZERO_MODE();
_MM_SET_FLUSH_ZERO_MODE(_MM_FLUSH_ZERO_ON);
#else
*state = 0;
#endif
PetscFunctionReturn(PETSC_SUCCESS);
}
/*
SlepcResetFlushToZero - Reset the FTZ flag in floating-point arithmetic.
*/
static inline PetscErrorCode SlepcResetFlushToZero(unsigned int *state)
{
PetscFunctionBegin;
#if defined(PETSC_HAVE_XMMINTRIN_H) && defined(_MM_FLUSH_ZERO_MASK) && defined(__SSE__)
_MM_SET_FLUSH_ZERO_MODE(*state & _MM_FLUSH_ZERO_MASK);
#else
*state = 0;
#endif
PetscFunctionReturn(PETSC_SUCCESS);
}
|