| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | #include "zagryadskov_m_complex_spmm_ccs/stl/include/ops_stl.hpp" | ||
| 2 | |||
| 3 | #include <complex> | ||
| 4 | #include <functional> | ||
| 5 | #include <thread> | ||
| 6 | #include <vector> | ||
| 7 | |||
| 8 | #include "util/include/util.hpp" | ||
| 9 | #include "zagryadskov_m_complex_spmm_ccs/common/include/common.hpp" | ||
| 10 | |||
| 11 | namespace zagryadskov_m_complex_spmm_ccs { | ||
| 12 | |||
| 13 |
1/2✓ Branch 2 taken 24 times.
✗ Branch 3 not taken.
|
24 | ZagryadskovMComplexSpMMCCSSTL::ZagryadskovMComplexSpMMCCSSTL(const InType &in) { |
| 14 | SetTypeOfTask(GetStaticTypeOfTask()); | ||
| 15 | GetInput() = in; | ||
| 16 |
1/2✓ Branch 1 taken 24 times.
✗ Branch 2 not taken.
|
24 | GetOutput() = CCS(); |
| 17 | 24 | } | |
| 18 | |||
| 19 | 60 | void ZagryadskovMComplexSpMMCCSSTL::SpMMSymbolic(const CCS &a, const CCS &b, std::vector<int> &col_ptr, int jstart, | |
| 20 | int jend) { | ||
| 21 | 60 | std::vector<int> marker(a.m, -1); | |
| 22 | |||
| 23 |
2/2✓ Branch 0 taken 56 times.
✓ Branch 1 taken 60 times.
|
116 | for (int j = jstart; j < jend; ++j) { |
| 24 | int count = 0; | ||
| 25 | |||
| 26 |
2/2✓ Branch 0 taken 80 times.
✓ Branch 1 taken 56 times.
|
136 | for (int k = b.col_ptr[j]; k < b.col_ptr[j + 1]; ++k) { |
| 27 | 80 | int b_row = b.row_ind[k]; | |
| 28 |
2/2✓ Branch 0 taken 88 times.
✓ Branch 1 taken 80 times.
|
168 | for (int zp = a.col_ptr[b_row]; zp < a.col_ptr[b_row + 1]; ++zp) { |
| 29 |
2/2✓ Branch 0 taken 72 times.
✓ Branch 1 taken 16 times.
|
88 | int a_row = a.row_ind[zp]; |
| 30 |
2/2✓ Branch 0 taken 72 times.
✓ Branch 1 taken 16 times.
|
88 | if (marker[a_row] != j) { |
| 31 | 72 | marker[a_row] = j; | |
| 32 | 72 | ++count; | |
| 33 | } | ||
| 34 | } | ||
| 35 | } | ||
| 36 | 56 | col_ptr[j + 1] += count; | |
| 37 | } | ||
| 38 | 60 | } | |
| 39 | |||
| 40 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 46 times.
|
56 | void ZagryadskovMComplexSpMMCCSSTL::SpMMKernel(const CCS &a, const CCS &b, CCS &c, const std::complex<double> &zero, |
| 41 | std::vector<int> &rows, std::vector<std::complex<double>> &acc, | ||
| 42 | std::vector<int> &marker, int j) { | ||
| 43 | rows.clear(); | ||
| 44 | 56 | int write_ptr = c.col_ptr[j]; | |
| 45 | |||
| 46 |
2/2✓ Branch 0 taken 80 times.
✓ Branch 1 taken 56 times.
|
136 | for (int k = b.col_ptr[j]; k < b.col_ptr[j + 1]; ++k) { |
| 47 | 80 | std::complex<double> tmpval = b.values[k]; | |
| 48 | 80 | int b_row = b.row_ind[k]; | |
| 49 |
2/2✓ Branch 0 taken 88 times.
✓ Branch 1 taken 80 times.
|
168 | for (int zp = a.col_ptr[b_row]; zp < a.col_ptr[b_row + 1]; ++zp) { |
| 50 |
2/2✓ Branch 0 taken 72 times.
✓ Branch 1 taken 16 times.
|
88 | int a_row = a.row_ind[zp]; |
| 51 |
2/2✓ Branch 0 taken 72 times.
✓ Branch 1 taken 16 times.
|
88 | acc[a_row] += tmpval * a.values[zp]; |
| 52 |
2/2✓ Branch 0 taken 72 times.
✓ Branch 1 taken 16 times.
|
88 | if (marker[a_row] != j) { |
| 53 |
2/2✓ Branch 0 taken 12 times.
✓ Branch 1 taken 60 times.
|
72 | marker[a_row] = j; |
| 54 | rows.push_back(a_row); | ||
| 55 | } | ||
| 56 | } | ||
| 57 | } | ||
| 58 | |||
| 59 |
2/2✓ Branch 0 taken 72 times.
✓ Branch 1 taken 56 times.
|
128 | for (int r_idx : rows) { |
| 60 | 72 | c.row_ind[write_ptr] = r_idx; | |
| 61 | 72 | c.values[write_ptr] = acc[r_idx]; | |
| 62 | 72 | ++write_ptr; | |
| 63 | 72 | acc[r_idx] = zero; | |
| 64 | } | ||
| 65 | 56 | } | |
| 66 | |||
| 67 | 60 | void ZagryadskovMComplexSpMMCCSSTL::SpMMNumeric(const CCS &a, const CCS &b, CCS &c, const std::complex<double> &zero, | |
| 68 | int jstart, int jend) { | ||
| 69 | 60 | std::vector<int> marker(a.m, -1); | |
| 70 |
1/4✓ Branch 1 taken 60 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
60 | std::vector<std::complex<double>> acc(a.m, zero); |
| 71 | 60 | std::vector<int> rows; | |
| 72 | |||
| 73 |
2/2✓ Branch 0 taken 56 times.
✓ Branch 1 taken 60 times.
|
116 | for (int j = jstart; j < jend; ++j) { |
| 74 |
1/2✓ Branch 1 taken 56 times.
✗ Branch 2 not taken.
|
56 | SpMMKernel(a, b, c, zero, rows, acc, marker, j); |
| 75 | } | ||
| 76 | 60 | } | |
| 77 | |||
| 78 | 24 | void ZagryadskovMComplexSpMMCCSSTL::SpMM(const CCS &a, const CCS &b, CCS &c) { | |
| 79 | 24 | c.m = a.m; | |
| 80 | 24 | c.n = b.n; | |
| 81 | 24 | const int num_threads = ppc::util::GetNumThreads(); | |
| 82 | 24 | std::vector<std::thread> threads(num_threads); | |
| 83 | |||
| 84 | 24 | std::complex<double> zero(0.0, 0.0); | |
| 85 |
1/2✓ Branch 1 taken 24 times.
✗ Branch 2 not taken.
|
24 | c.col_ptr.assign(c.n + 1, 0); |
| 86 | |||
| 87 |
2/2✓ Branch 0 taken 60 times.
✓ Branch 1 taken 24 times.
|
84 | for (int tid = 0; tid < num_threads; ++tid) { |
| 88 | 60 | int jstart = (tid * b.n) / num_threads; | |
| 89 | 60 | int jend = ((tid + 1) * b.n) / num_threads; | |
| 90 |
2/4✓ Branch 1 taken 60 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 60 times.
|
60 | threads[tid] = std::thread(SpMMSymbolic, std::cref(a), std::cref(b), std::ref(c.col_ptr), jstart, jend); |
| 91 | } | ||
| 92 |
2/2✓ Branch 0 taken 60 times.
✓ Branch 1 taken 24 times.
|
84 | for (auto &th : threads) { |
| 93 |
1/2✓ Branch 1 taken 60 times.
✗ Branch 2 not taken.
|
60 | th.join(); |
| 94 | } | ||
| 95 | |||
| 96 |
2/2✓ Branch 0 taken 56 times.
✓ Branch 1 taken 24 times.
|
80 | for (int j = 0; j < c.n; ++j) { |
| 97 | 56 | c.col_ptr[j + 1] += c.col_ptr[j]; | |
| 98 | } | ||
| 99 |
1/2✓ Branch 1 taken 24 times.
✗ Branch 2 not taken.
|
24 | int nnz = c.col_ptr[b.n]; |
| 100 |
1/2✓ Branch 1 taken 24 times.
✗ Branch 2 not taken.
|
24 | c.row_ind.resize(nnz); |
| 101 |
1/2✓ Branch 1 taken 24 times.
✗ Branch 2 not taken.
|
24 | c.values.resize(nnz); |
| 102 | |||
| 103 |
2/2✓ Branch 0 taken 60 times.
✓ Branch 1 taken 24 times.
|
84 | for (int tid = 0; tid < num_threads; ++tid) { |
| 104 | 60 | int jstart = (tid * b.n) / num_threads; | |
| 105 | 60 | int jend = ((tid + 1) * b.n) / num_threads; | |
| 106 |
2/4✓ Branch 1 taken 60 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 60 times.
|
60 | threads[tid] = std::thread(SpMMNumeric, std::cref(a), std::cref(b), std::ref(c), std::cref(zero), jstart, jend); |
| 107 | } | ||
| 108 |
2/2✓ Branch 0 taken 60 times.
✓ Branch 1 taken 24 times.
|
84 | for (auto &th : threads) { |
| 109 |
1/2✓ Branch 1 taken 60 times.
✗ Branch 2 not taken.
|
60 | th.join(); |
| 110 | } | ||
| 111 | 24 | } | |
| 112 | |||
| 113 | 24 | bool ZagryadskovMComplexSpMMCCSSTL::ValidationImpl() { | |
| 114 | const CCS &a = std::get<0>(GetInput()); | ||
| 115 | const CCS &b = std::get<1>(GetInput()); | ||
| 116 | 24 | return a.n == b.m; | |
| 117 | } | ||
| 118 | |||
| 119 | 24 | bool ZagryadskovMComplexSpMMCCSSTL::PreProcessingImpl() { | |
| 120 | 24 | return true; | |
| 121 | } | ||
| 122 | |||
| 123 | 24 | bool ZagryadskovMComplexSpMMCCSSTL::RunImpl() { | |
| 124 | const CCS &a = std::get<0>(GetInput()); | ||
| 125 | const CCS &b = std::get<1>(GetInput()); | ||
| 126 | CCS &c = GetOutput(); | ||
| 127 | |||
| 128 | 24 | SpMM(a, b, c); | |
| 129 | |||
| 130 | 24 | return true; | |
| 131 | } | ||
| 132 | |||
| 133 | 24 | bool ZagryadskovMComplexSpMMCCSSTL::PostProcessingImpl() { | |
| 134 | 24 | return true; | |
| 135 | } | ||
| 136 | |||
| 137 | } // namespace zagryadskov_m_complex_spmm_ccs | ||
| 138 |