.. _program_listing_file_src_utilities.tpp: Program Listing for File utilities.tpp ====================================== |exhale_lsh| :ref:`Return to documentation for file ` (``src/utilities.tpp``) .. |exhale_lsh| unicode:: U+021B0 .. UPWARDS ARROW WITH TIP LEFTWARDS .. code-block:: cpp /* This file is part of brille. Copyright © 2019,2020 Greg Tucker brille is free software: you can redistribute it and/or modify it under the terms of the GNU Affero General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. brille is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Affero General Public License for more details. You should have received a copy of the GNU Affero General Public License along with brille. If not, see . */ // trace of a square matrix template T trace(const T *M){ T out = T(0); for (int i=0; i void copy_array(T *D, const T *S){ for (int i=0; i void copy_matrix(T *M, const T *A){ copy_array(M,A); } template void copy_vector(T *V, const T *A){ copy_array(V,A); } // element-wise checking of approximate equality template bool equal_array(const T *A, const T *B, const T tol){ for (int i=0; i tol ) return false; return true; } template bool equal_matrix(const T *A, const T *B, const T tol){ return equal_array(A,B,tol); } template bool equal_vector(const T *A, const T *B, const T tol){ return equal_array(A,B,tol); } // array multiplication C = A * B -- where C is (N,M), A is (N,I) and B is (I,M) template void multiply_arrays(T *C, const R *A, const S *B){ for (size_t i=0;i void multiply_matrix_matrix(T *C, const R *A, const S *B){ multiply_arrays(C,A,B); } template void multiply_matrix_vector(T *C, const R *A, const S *b){ multiply_arrays(C,A,b); } template void multiply_vector_matrix(T *C, const R *a, const S *B){ multiply_arrays(C,a,B); } // array multiplication specialization for non-complex * complex arrays. template // std::enable_if_t> void mul_arrays(T* C, const I n, const I l, const I m, const R* A, const S* B){ for (I i=0;i(A[i*l+k]*B[k*m+j]); } template // std::enable_if_t> void mul_arrays(std::complex* C, const I n, const I l, const I m, const R* A, const std::complex* B){ for (I i=0;i(0); for (I i=0;i(A[i*l+k])*B[k*m+j]; } template // std::enable_if_t> void mul_arrays(std::complex* C, const I n, const I l, const I m, const std::complex* A, const S* B){ for (I i=0;i(0); for (I i=0;i(B[k*m+j]); } template // std::enable_if_t> void mul_arrays(std::complex* C, const I n, const I l, const I m, const std::complex* A, const std::complex* B){ for (I i=0;i(0); for (I i=0;i void add_arrays(T *C, const R *A, const S *B){ for (int i=0; i void add_matrix(T *C, const R *A, const S *B){ add_arrays(C,A,B); } // array element-wise subtraction template void subtract_arrays(T *C, const R *A, const S *B){ for (int i=0; i void subtract_matrix(T *C, const R *A, const S *B){ subtract_arrays(C,A,B); } // specialized casting of doubles to ints template T my_cast(const R a){ return T(a); } // inlined to avoid redefining this specialization each time the header is called template<> inline int my_cast(const double a){ return int( a<0 ? a-0.5 : a+0.5); } // element-wise re-(special)-casting of arrays template void cast_array(T *A, const R *B){ for (int i=0; i(B[i]); } template void cast_matrix(T *A, const R *B){ cast_array(A,B); } template void cast_vector(T *A, const R *B){ cast_array(A,B); } // The cofactor array C_ij of the N by M array A is a M-1 by N-1 array of transpose(A with row i and column j) missing: template void array_cofactor(R *C, const R *A, const int i, const int j, const int N, const int M){ int k = 0; for (int m=0; m void matrix_cofactor(R *C, const R *A, const int i, const int j, const int N){ array_cofactor(C,A,i,j,N,N); } // the determinant is only defined for square matrices template R matrix_determinant(const R *M, const int N){ if (1==N) return M[0]; // the determinant: R d{0}; // temporary cofactor storage R *cof = nullptr; cof = new R[(N-1)*(N-1)](); // loop over one row/column (let's go with row, for fun) R pm = 1; // +/-1 for alternating elements for (int i=0; i void matrix_adjoint(R *A, const R *M, const int N){ if (1==N){ A[0]=R(1); return; } R pm=1, *cof =nullptr; cof = new R[(N-1)*(N-1)](); for (int i=0; i R matrix_determinant_and_inverse(R *invM, const R *M, const R tol, const int N){ R d = matrix_determinant(M,N); if (std::abs(d) > tol ){ // the matrix is *not* singular and has an inverse matrix_adjoint(invM,M,N); // at this point invM is the Adjoint(M) // The inv(M) = adjoint(M)/det(M) int nn = N*N; for (int i=0; i bool matrix_inverse(R *invM, const R *M, R tol, const int N){ return ( std::abs( matrix_determinant_and_inverse(invM, M, tol, N) ) > tol) ; } template bool similar_matrix(T *M, const T *A, const T *B, const T tol){ T *C = nullptr; C = new T[N*N](); bool ok = matrix_inverse(C,B,tol,N); if ( ok ){ T *P = nullptr; P = new T[N*N](); multiply_matrix_matrix(P,A,B); multiply_matrix_matrix(M,C,P); delete[] P; } else { printf("spglib: No similar matrix due to 0 determinant.\n"); } delete[] C; return ok; } template void array_transpose(R *D, const R *S){ for (int i=0; i void matrix_transpose(R *D, const R *S){ array_transpose(D,S); } // template void complex_conj_array_transpose(std::complex *D, const std::complex *S){ for (int i=0; i void complex_conj_matrix_transpose(std::complex *D, const std::complex *S){ complex_conj_array_transpose(D,S); } // in place transpose: template void matrix_transpose(R *B){ R t; for (int i=0; i void matrix_metric(R *M, const R *L){ R *Lt = nullptr; Lt = new R[N*N](); matrix_transpose(Lt,L); multiply_matrix_matrix(M,Lt,L); delete[] Lt; } template R vector_norm_squared(const R *v){ R vv = 0; for (int i=0; i void vector_cross(T *c, const R *a, const S *b){ // if (3!=N) // throw std::domain_error("The cross product is only defined for 3-vectors"); // c[0] = static_cast(a[1])*static_cast(b[2]) - static_cast(a[2])*static_cast(b[1]); // c[1] = static_cast(a[2])*static_cast(b[0]) - static_cast(a[0])*static_cast(b[2]); // c[2] = static_cast(a[0])*static_cast(b[1]) - static_cast(a[1])*static_cast(b[0]); //} template S vector_dot(const size_t n, const T* a, const R* b){ S out = 0; for (size_t i=0; i(a[i])*static_cast(b[i]); return out; } template R vector_dot(const R *a, const R *b){ R out = 0; for (int i=0; i T mod1(const T a){ T b = a - my_cast(a); return ( (b < T(0) - 1e-10 ) ? b + T(1) : b ); } template bool is_int_matrix(const T * A, const R tol){ for (int i=0; i(A[i]) - A[i]) > tol ) return false; return true; } template bool is_int_matrix(const int *, const R){ return true; } template T frobenius_distance(const I n, const T* A, const T* B) { // sqrt(tr([A-B]*[A-B]')) == sqrt(A00^2 + A01^2 + A02^2 ... Ann^2) return vector_distance(n*n, A, B); } template T frobenius_distance(const I n, const std::complex* A, const std::complex* B) { // sqrt( tr( (A-B)x(A-B)') ) is identically the same as the vector distance of 9-vectors return vector_distance(n*n, A, B); } template T vector_angle(const I n, const T* A, const T* B){ T AA=0, BB=0, AB=0; for (I i=0; i1){ c_t /= act; act = std::abs(c_t); } if (act>1){ std::string msg = "vector angle cos(theta)=" + my_to_string(c_t) + " = " + my_to_string(AB) + "/(" + my_to_string(nA) + "×" + my_to_string(nB) + ")"; throw std::runtime_error(msg); } return std::acos(c_t); } template T vector_angle(const I n, const std::complex* A, const std::complex* B){ // return hermitian_angle(n,A,B); return euclidean_angle(n,A,B); } template T euclidean_angle(const I n, const std::complex* A, const std::complex* B){ T AB{0}, nA{0}, nB{0}, c_t; // Compute the products of complex n-vectors as if they were real 2n-vectors for (I i=0; i1){ c_t /= act; act = std::abs(c_t); } if (act>1){ std::string msg = "Complex Euclidean angle cos(theta)=" + my_to_string(c_t) + " = " + my_to_string(AB) + "/(" + my_to_string(nA) + "×" + my_to_string(nB) + ")"; throw std::runtime_error(msg); } return std::acos(c_t); } template T hermitian_product(const I n, const T* a, const T* b){ T h_dot{0}; for (I i=0; i std::complex hermitian_product(const I n, const T* a, const std::complex* b){ std::complex h_dot{0,0}; for (I i=0; i std::complex complex_prod(const std::complex& x, T y){ return std::conj(x)*y; } template std::complex complex_prod(const std::complex& x, const std::complex& y){ T xr{x.real()}, xi{x.imag()}, yr{y.real()}, yi{y.imag()}; // a*×b = (a'+ia")*×(b'+ib") = (a'-ia")×(b'+ib") = a'b' + a"b" + i(a'b" -a"b') return std::complex(xr*yr+xi*yi, xr*yi-xi*yr); } template std::complex hermitian_product(const I n, const std::complex* a, const T* b){ std::complex h_dot{0,0}; // for (size_t i=0; i std::complex hermitian_product(const I n, const std::complex* a, const std::complex* b){ T hr{0}, hi{0}; long long sn = u2s(n); #pragma omp parallel for shared(a,b) reduction(+: hr,hi) for (long long si=0; si h = std::conj(a[si])*b[si]; hr += h.real(); hi += h.imag(); } return std::complex(hr,hi); } template T hermitian_angle(const I n, const T* A, const T* B){ return vector_angle(n,A,B); } template T hermitian_angle(const I n, const std::complex* A, const std::complex* B){ T nAB, nA, nB, c_t; nAB = std::sqrt(vector_product(n,A,B)); nA = std::sqrt(std::real(hermitian_product(n,A,A))); nB = std::sqrt(std::real(hermitian_product(n,B,B))); if (nA && nB){ c_t = nAB/(nA*nB); } else { // if both are zero their angle is strictly undefined, but for our // purposes it would be better to set the angle to 0 [and cos(0)=1] c_t = (nA || nB) ? T(0) : T(1); } T act = std::abs(c_t); if (brille::approx::scalar(act,1.0) && act>1){ c_t/=act; // force close-to-one values to one, maintaining the sign act = std::abs(c_t); } if (act>1){ std::string msg = "Complex Hermitian angle cos(theta)=" + my_to_string(c_t) + " = " + my_to_string(nAB) + "/(" + my_to_string(nA) + "×" + my_to_string(nB) + ")"; throw std::runtime_error(msg); } return std::acos(c_t); } template T vector_distance(const I n, const T* a, const T* b){ T d, s=0; for (I i=0; i T vector_distance(const I n, const std::complex* a, const std::complex* b){ std::complex d; T s=0; for (I i=0; i T vector_product(const I n, const T* a, const T* b){ T h_dot{0}; for (size_t i=0; i T vector_product(const I n, const T* a, const std::complex* b){ std::complex h_dot = hermitian_product(n,a,b); return std::real(h_dot*std::conj(h_dot)); } template T vector_product(const I n, const std::complex* a, const T* b){ std::complex h_dot = hermitian_product(n,a,b); return std::real(h_dot*std::conj(h_dot)); } template T vector_product(const I n, const std::complex* a, const std::complex* b){ std::complex h_dot = hermitian_product(n,a,b); return std::real(h_dot*std::conj(h_dot)); } template T inner_product(const I n, const T* a, const T* b){ T h_dot{0}; for (I i=0; i T inner_product(const I n, const std::complex* a, const std::complex* b){ std::complex h_dot = hermitian_product(n,a,b); return std::real(h_dot); } template T squared_distance(const T&A, const T& B){ return (A-B)*(A-B); } template T squared_distance(const std::complex& A, const std::complex& B){ T r = std::real(A)-std::real(B); T i = std::imag(A)-std::imag(B); return r*r + i*i; } template T magnitude(const T a){ return std::abs(a); } template T magnitude(const std::complex a){ return std::sqrt(std::real(a*std::conj(a))); } /*! Determine an unsigned integer that encodes the signs of a complex vector. Each element of a complex valued vector has two signs associated with it, one for the real part and one for the imaginary part. Some properties of the vector, notably the squared modulus of its dot product with a real vector, depend on the relative signs of each element. This function encodes the four possible combinations into a quaternary based system, (++,+-,-+ ,--)→(0,1,2,3)₄ and combines the `n` values into a single unsigned integer using 2`n` bits. As an example, a vector (+a-𝑖b, -c+𝑖d, -e-𝑖f) has signs (++,-+,--) which give the quaternary numbers (0,2,3)₄ ⇒ 001011₂ = 9. */ template size_t encode_array_signs(const size_t n, const std::complex* a){ size_t signs=0; for (size_t i=0; i size_t encode_array_signs(const size_t n, const T* a){ size_t signs=0; for (size_t i=0; i int make_eigenvectors_equivalent(const size_t n, const T* v0, T* v1){ size_t s0, s1; s0 = encode_array_signs(n,v0); s1 = encode_array_signs(n,v1); if (s0 == s1) return 0; // the only valid permutation for real values is by 2. size_t onesign, s1mod{0}; for (size_t j=0; j> 2*(n-1-j)) - ((s1 >> 2*(n-j)) << 2*(n-j)); // permute it by 2 onesign = (onesign+2)%4; // and stick it in s1mod s1mod += onesign << 2*(n-1-j); } if (s0 == s1mod){ T m1{-1}; for (size_t j=0; j int make_eigenvectors_equivalent(const size_t n, const std::complex* v0, std::complex v1){ size_t s0, s1; s0 = encode_array_signs(n,v0); s1 = encode_array_signs(n,v1); if (s0 == s1) return 0; // The signs of each vector are not the same, so check for equivalence: size_t onesign; // we can permute each sign quaternary number up to three times for (size_t i=1; i<4; ++i){ size_t s1mod{0}; for (size_t j=0; j> 2*(n-1-j)) - ((s1 >> 2*(n-j)) << 2*(n-j)); // permute it onesign = (onesign+i)%4; // and stick it in s1mod s1mod += onesign << 2*(n-1-j); } if (s0 == s1mod){ // The two are equivalent for our purposes, but the signs of v1 need to be changed! T m1{-1}; switch (i){ case 1: // exchange the imaginary sign for (size_t j=0; j(std::real(v1[j]), m1*std::imag(v1[j])); break; case 2: // exchange both real and imaginary signs for (size_t j=0; j(m1*std::real(v1[j]), m1*std::imag(v1[j])); break; case 3: // exchange the real sign for (size_t j=0; j(m1*std::real(v1[j]), std::imag(v1[j])); break; } return 0; } } } template T antiphase(const T){ //return std::signbit(z) ? T(-1) : T(1); return T(1); } template std::complex antiphase(const std::complex z){ return std::polar(T(1),-std::arg(z)); } template T antiphase(const I, const T*, const T*){ return T(1); } template std::complex antiphase(const I n, const std::complex* a, const std::complex* b){ T real_dot{0}, imag_dot{0}; for (I i=0; i void inplace_antiphase(const I, const T*, const T*, T*) {} template void inplace_antiphase(const I n, const std::complex* a, const std::complex* b, std::complex* phased) { std::complex eith = antiphase(n, a, b); for (I i=0; i class A> T unsafe_antiphase(const A&, const A&){ return T(1); } template class A> std::complex unsafe_antiphase(const A>& a, const A>& b){ T real_dot{0}, imag_dot{0}; // a safe version of this would first ensure that a and b are the same shape for (auto [ox, ax, bx]: a.broadcastItr(b)){ std::complex av{a[ax]}, bv{b[bx]}; T ar{av.real()}, ai{av.imag()}, br{bv.real()}, bi{bv.imag()}; real_dot += ar*br + ai*bi; imag_dot += ar*bi - ai*br; } return std::polar(T(1), T(-1)*std::atan2(imag_dot, real_dot)); } template::type> S coth_over_en(const T en, const R beta){ S Sen = static_cast(en); S Sbeta = static_cast(beta); S den = std::tanh(Sen*Sbeta/S(2))*Sen; return S(1)/den; } template::type> S coth_over_en(const std::complex en, const R beta){ S den = static_cast(std::real(std::tanh(en*beta*0.5)*en)); return S(1)/den; } template std::enable_if_t::value, T> gcd(const T a, const T b){ // We include numeric for std::gcd in C++17 /* Pre-C++17 gcd was not part of the standard but may be present in an // experimental header -- which of course doesn't exist for MSVC */ #if defined(__cplusplus) && __cplusplus > 201700L return std::gcd(a,b); #else if (b==0) return a; return gcd(b, a % b); #endif } template std::enable_if_t::value&&std::is_unsigned::value, unsigned long long> binomial_coefficient(const T n, const R k){ unsigned long long ans{1}, num{1}, den{1}, comdiv; if (k>n){ std::string msg = "The binomial coefficient requires more choices than selections."; msg += "\n binomial_coefficient(" + std::to_string(n) + ":" + std::to_string(k) + ") is invalid"; throw std::domain_error(msg); } // the Binomial coefficient is symmetric due to the denominator k!(n-k)! R m = (n-k < k) ? static_cast(n-k) : k; bool overflow = false; for (R i=0; i(n-i); den *= static_cast(i)+1; if (lastnum > num || lastden > den){ comdiv = gcd(lastnum, lastden); if (comdiv > 1){ num = lastnum/comdiv; den = lastden/comdiv; --i; } else { overflow = true; } } if (overflow) break; } if (overflow){ long double dans{1}; for (T i=0; i(n-i)/(static_cast(i)+1); ans = std::llround(dans); } else { ans = num/den; } return ans; } template S u2s(const U u){ S smax = (std::numeric_limits::max)(); U usmax = static_cast(smax); if (u > usmax){ std::string msg = "unsigned_to_signed:: Value " + std::to_string(u) + " can not be stored in requested signed type" + " (max value "+ std::to_string(smax) + ")"; throw std::overflow_error(msg); } return static_cast(u); } template U s2u(const S s){ if (s < 0 //|| static_cast((std::numeric_limits::max)()) > (std::numeric_limits::max)() // the largest element of type S is too big for U //|| s > static_cast((std::numeric_limits::max)()) // s, specifically is too big to be expressed in type U ){ std::string msg = "signed_to_unsiged:: Value " + std::to_string(s) + " can not be stored in requested unsiged type" + " (max value "+ std::to_string((std::numeric_limits::max)()) + ")"; throw std::overflow_error(msg); } return static_cast(s); } // from https://en.cppreference.com/w/cpp/language/sizeof... template constexpr auto make_array(Ts&&... ts) -> std::array,sizeof...(Ts)> { return { std::forward(ts)... }; } template constexpr auto make_vector(Ts&&... ts) -> std::vector> { return { std::forward(ts)... }; } template using are_same = std::conjunction...>;