16#ifndef dealii_lac_utilities_h
17#define dealii_lac_utilities_h
63 template <
typename NumberType>
64 std::array<NumberType, 3>
93 template <
typename NumberType>
94 std::array<NumberType, 3>
132 template <
typename OperatorType,
typename VectorType>
135 const VectorType &
v0,
136 const unsigned int k,
187 template <
typename OperatorType,
typename VectorType>
190 const OperatorType & H,
191 const unsigned int n,
192 const std::pair<double, double> unwanted_spectrum,
207 namespace UtilitiesImplementation
214 template <
typename Number>
231 template <
typename NumberType>
232 std::array<std::complex<NumberType>, 3>
234 const std::complex<NumberType> & )
237 std::array<NumberType, 3> res;
243 template <
typename NumberType>
244 std::array<NumberType, 3>
248 const NumberType tau = g / f;
251 "real-valued Hyperbolic rotation does not exist for (" +
256 std::array<NumberType, 3> csr;
258 csr[1] = csr[0] * tau;
265 template <
typename NumberType>
266 std::array<std::complex<NumberType>, 3>
268 const std::complex<NumberType> & )
271 std::array<NumberType, 3> res;
277 template <
typename NumberType>
278 std::array<NumberType, 3>
281 std::array<NumberType, 3> res;
295 if (g == NumberType())
298 res[1] = NumberType();
301 else if (f == NumberType())
303 res[0] = NumberType();
309 const NumberType tau = g / f;
312 res[1] = res[0] * tau;
317 const NumberType tau = f / g;
320 res[0] = res[1] * tau;
329 template <
typename OperatorType,
typename VectorType>
332 const VectorType & v0_,
333 const unsigned int k,
350 std::vector<double> subdiagonal;
354 double a = v->l2_norm();
365 for (
unsigned int i = 1; i < k; ++i)
368 const double b = f->l2_norm();
382 subdiagonal.push_back(
b);
391 std::vector<double> Z;
393 std::vector<double> work;
417 return diagonal[k - 1] + f->l2_norm();
421 template <
typename OperatorType,
typename VectorType>
424 const OperatorType & op,
425 const unsigned int degree,
426 const std::pair<double, double> unwanted_spectrum,
430 const double a = unwanted_spectrum.first;
431 const double b = unwanted_spectrum.second;
434 const bool scale = (a_L < std::numeric_limits<double>::infinity());
438 "Lower bound of the unwanted spectrum should be smaller than the upper bound."));
442 "Scaling point should be outside of the unwanted spectrum."));
452 VectorType &y = *p_y;
453 VectorType &yn = *p_yn;
468 const double e = (
b - a) / 2.;
469 const double c = (a +
b) / 2.;
470 const double alpha = 1. /
e;
471 const double beta = -c /
e;
473 const double sigma1 =
475 double sigma =
scale ? sigma1 : 1.;
476 const double tau = 2. / sigma;
478 y.sadd(alpha * sigma, beta * sigma, x);
480 for (
unsigned int i = 2; i <= degree; ++i)
482 const double sigma_new =
scale ? 1. / (tau - sigma) : 1.;
484 yn.sadd(2. * alpha * sigma_new, 2. * beta * sigma_new, y);
485 yn.add(-sigma * sigma_new, x);
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcErrorCode(std::string arg1, types::blas_int arg2)
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDivideByZero()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Expression copysign(const Expression &value, const Expression &sign)
@ eigenvalues
Eigenvalue vector is filled.
@ diagonal
Matrix is diagonal.
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void chebyshev_filter(VectorType &x, const OperatorType &H, const unsigned int n, const std::pair< double, double > unwanted_spectrum, const double tau, VectorMemory< VectorType > &vector_memory)
std::array< NumberType, 3 > givens_rotation(const NumberType &x, const NumberType &y)
std::array< NumberType, 3 > hyperbolic_rotation(const NumberType &x, const NumberType &y)
double lanczos_largest_eigenvalue(const OperatorType &H, const VectorType &v0, const unsigned int k, VectorMemory< VectorType > &vector_memory, std::vector< double > *eigenvalues=nullptr)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
void call_stev(const char jobz, const types::blas_int n, Number *d, Number *e, Number *z, const types::blas_int ldz, Number *work, types::blas_int *info)
void copy(const T *begin, const T *end, U *dest)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)