Skip to content

Commit

Permalink
added complex matrix apply template instantiation
Browse files Browse the repository at this point in the history
  • Loading branch information
mulaga committed Aug 26, 2024
1 parent 4933360 commit a48ebb9
Show file tree
Hide file tree
Showing 2 changed files with 14 additions and 0 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -30,6 +30,11 @@ void test_spinhalf_symmetric_apply(OpSum ops, PermutationGroup space_group,
arma::cx_vec w2(block.size(), arma::fill::zeros);
apply(ops, block, v, block, w2);
REQUIRE(close(w1, w2));
arma::cx_mat m(block.size(), 5, arma::fill::randn);
arma::cx_mat n1 = H * m;
arma::cx_mat n2(block.size(), 5, arma::fill::zeros);
apply(ops, block, m, block, n2);
REQUIRE(close(n1, n2));

// Compute eigenvalues and compare
arma::vec evals_mat;
Expand Down
9 changes: 9 additions & 0 deletions xdiag/algebra/apply.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -141,6 +141,9 @@ void apply(OpSum const &ops, Spinhalf const &block_in,
template void apply<double>(OpSum const &, Spinhalf const &,
arma::Mat<double> const &, Spinhalf const &,
arma::Mat<double> &, double);
template void apply<complex>(OpSum const &, Spinhalf const &,
arma::Mat<complex> const &, Spinhalf const &,
arma::Mat<complex> &, double);

template <typename coeff_t>
void apply(OpSum const &ops, tJ const &block_in,
Expand Down Expand Up @@ -177,6 +180,9 @@ void apply(OpSum const &ops, tJ const &block_in,
template void apply<double>(OpSum const &, tJ const &,
arma::Mat<double> const &, tJ const &,
arma::Mat<double> &, double);
template void apply<complex>(OpSum const &, tJ const &,
arma::Mat<complex> const &, tJ const &,
arma::Mat<complex> &, double);

template <typename coeff_t>
void apply(OpSum const &ops, Electron const &block_in,
Expand Down Expand Up @@ -213,6 +219,9 @@ void apply(OpSum const &ops, Electron const &block_in,
template void apply<double>(OpSum const &, Electron const &,
arma::Mat<double> const &, Electron const &,
arma::Mat<double> &, double);
template void apply<complex>(OpSum const &, Electron const &,
arma::Mat<complex> const &, Electron const &,
arma::Mat<complex> &, double);

#ifdef XDIAG_USE_MPI

Expand Down

0 comments on commit a48ebb9

Please sign in to comment.