| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | /* | ||
| 2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
| 4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
| 5 | |||
| 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 | */ | ||
| 14 | |||
| 15 | #pragma once | ||
| 16 | |||
| 17 | /* SUBMANSEC = Sys */ | ||
| 18 | |||
| 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 | ||
| 31 | |||
| 32 | 31561 | static inline PetscReal SlepcDefaultTol(PetscReal tol) | |
| 33 | { | ||
| 34 |
5/6✓ Branch 0 taken 186 times.
✓ Branch 1 taken 189 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 16 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
|
31561 | return tol == (PetscReal)PETSC_DETERMINE ? SLEPC_DEFAULT_TOL : tol; |
| 35 | } | ||
| 36 | |||
| 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`. | ||
| 40 | |||
| 41 | Not Collective | ||
| 42 | |||
| 43 | Input parameters: | ||
| 44 | . x,y - the real numbers | ||
| 45 | |||
| 46 | Output parameter: | ||
| 47 | . return - the result | ||
| 48 | |||
| 49 | Fortran Note: | ||
| 50 | This function is not available from Fortran. | ||
| 51 | |||
| 52 | Level: developer | ||
| 53 | @*/ | ||
| 54 | 22239644 | static inline PetscReal SlepcAbs(PetscReal x,PetscReal y) | |
| 55 | { | ||
| 56 |
4/4✓ Branch 0 taken 220 times.
✓ Branch 1 taken 316 times.
✓ Branch 2 taken 150 times.
✓ Branch 3 taken 316 times.
|
22239644 | PetscReal w,z,t,xabs=PetscAbs(x),yabs=PetscAbs(y); |
| 57 | |||
| 58 |
2/2✓ Branch 0 taken 316 times.
✓ Branch 1 taken 246 times.
|
22239644 | w = PetscMax(xabs,yabs); |
| 59 |
2/2✓ Branch 0 taken 316 times.
✓ Branch 1 taken 246 times.
|
22239644 | z = PetscMin(xabs,yabs); |
| 60 |
2/2✓ Branch 0 taken 256 times.
✓ Branch 1 taken 223 times.
|
22239644 | if (PetscUnlikely(z == PetscRealConstant(0.0))) return w; |
| 61 | 1524768 | t = z/w; | |
| 62 | 1524768 | return w*PetscSqrtReal(PetscRealConstant(1.0)+t*t); | |
| 63 | ✗ | } | |
| 64 | |||
| 65 | /*MC | ||
| 66 | SlepcAbsEigenvalue - Returns the absolute value of a complex number given | ||
| 67 | its real and imaginary parts. | ||
| 68 | |||
| 69 | Synopsis: | ||
| 70 | PetscReal SlepcAbsEigenvalue(PetscScalar x,PetscScalar y) | ||
| 71 | |||
| 72 | Not Collective | ||
| 73 | |||
| 74 | Input parameters: | ||
| 75 | + x - the real part of the complex number | ||
| 76 | - y - the imaginary part of the complex number | ||
| 77 | |||
| 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`. | ||
| 81 | |||
| 82 | In complex scalars, only the first argument is used, i.e., the result is $|x|$. | ||
| 83 | |||
| 84 | Fortran Note: | ||
| 85 | This function is not available from Fortran. | ||
| 86 | |||
| 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 | ||
| 94 | |||
| 95 | /* | ||
| 96 | SlepcSetFlushToZero - Set the FTZ flag in floating-point arithmetic. | ||
| 97 | */ | ||
| 98 | 219 | static inline PetscErrorCode SlepcSetFlushToZero(unsigned int *state) | |
| 99 | { | ||
| 100 |
1/2✓ Branch 0 taken 4 times.
✗ Branch 1 not taken.
|
219 | 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 | 738 | *state = 0; | |
| 106 | #endif | ||
| 107 |
6/12✓ Branch 0 taken 4 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 4 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 4 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 4 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 4 times.
|
219 | PetscFunctionReturn(PETSC_SUCCESS); |
| 108 | } | ||
| 109 | |||
| 110 | /* | ||
| 111 | SlepcResetFlushToZero - Reset the FTZ flag in floating-point arithmetic. | ||
| 112 | */ | ||
| 113 | 144 | static inline PetscErrorCode SlepcResetFlushToZero(unsigned int *state) | |
| 114 | { | ||
| 115 |
1/2✓ Branch 0 taken 4 times.
✗ Branch 1 not taken.
|
144 | 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 | 738 | *state = 0; | |
| 120 | #endif | ||
| 121 |
6/12✓ Branch 0 taken 4 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 4 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 4 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 4 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 4 times.
|
144 | PetscFunctionReturn(PETSC_SUCCESS); |
| 122 | } | ||
| 123 |