| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | #include "maslova_u_mult_matr_crs/tbb/include/ops_tbb.hpp" | ||
| 2 | |||
| 3 | #include <tbb/tbb.h> | ||
| 4 | |||
| 5 | #include <algorithm> | ||
| 6 | #include <vector> | ||
| 7 | |||
| 8 | #include "maslova_u_mult_matr_crs/common/include/common.hpp" | ||
| 9 | |||
| 10 | namespace maslova_u_mult_matr_crs { | ||
| 11 | |||
| 12 |
1/2✓ Branch 2 taken 20 times.
✗ Branch 3 not taken.
|
20 | MaslovaUMultMatrTBB::MaslovaUMultMatrTBB(const InType &in) { |
| 13 | SetTypeOfTask(GetStaticTypeOfTask()); | ||
| 14 | GetInput() = in; | ||
| 15 | 20 | } | |
| 16 | |||
| 17 | 20 | bool MaslovaUMultMatrTBB::ValidationImpl() { | |
| 18 | const auto &input = GetInput(); | ||
| 19 | const auto &a = std::get<0>(input); | ||
| 20 | const auto &b = std::get<1>(input); | ||
| 21 |
3/6✓ Branch 0 taken 20 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 20 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 20 times.
|
20 | return (a.cols == b.rows && a.rows > 0 && b.cols > 0); |
| 22 | } | ||
| 23 | |||
| 24 | 20 | bool MaslovaUMultMatrTBB::PreProcessingImpl() { | |
| 25 | const auto &a = std::get<0>(GetInput()); | ||
| 26 | const auto &b = std::get<1>(GetInput()); | ||
| 27 | auto &c = GetOutput(); | ||
| 28 | 20 | c.rows = a.rows; | |
| 29 | 20 | c.cols = b.cols; | |
| 30 | 20 | return true; | |
| 31 | } | ||
| 32 | |||
| 33 | 32 | int MaslovaUMultMatrTBB::GetRowNNZ(int i, const CRSMatrix &a, const CRSMatrix &b, std::vector<int> &marker) { | |
| 34 | int row_nnz = 0; | ||
| 35 |
2/2✓ Branch 0 taken 48 times.
✓ Branch 1 taken 32 times.
|
80 | for (int j = a.row_ptr[i]; j < a.row_ptr[i + 1]; ++j) { |
| 36 | 48 | int col_a = a.col_ind[j]; | |
| 37 |
2/2✓ Branch 0 taken 72 times.
✓ Branch 1 taken 48 times.
|
120 | for (int k = b.row_ptr[col_a]; k < b.row_ptr[col_a + 1]; ++k) { |
| 38 |
2/2✓ Branch 0 taken 44 times.
✓ Branch 1 taken 28 times.
|
72 | int col_b = b.col_ind[k]; |
| 39 |
2/2✓ Branch 0 taken 44 times.
✓ Branch 1 taken 28 times.
|
72 | if (marker[col_b] != i) { |
| 40 | 44 | marker[col_b] = i; | |
| 41 | 44 | row_nnz++; | |
| 42 | } | ||
| 43 | } | ||
| 44 | } | ||
| 45 | 32 | return row_nnz; | |
| 46 | } | ||
| 47 | |||
| 48 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 32 times.
|
32 | void MaslovaUMultMatrTBB::FillRowValues(int i, const CRSMatrix &a, const CRSMatrix &b, CRSMatrix &c, |
| 49 | std::vector<double> &acc, std::vector<int> &marker, std::vector<int> &used) { | ||
| 50 | used.clear(); | ||
| 51 |
2/2✓ Branch 0 taken 48 times.
✓ Branch 1 taken 32 times.
|
80 | for (int j = a.row_ptr[i]; j < a.row_ptr[i + 1]; ++j) { |
| 52 | 48 | int col_a = a.col_ind[j]; | |
| 53 | 48 | double val_a = a.values[j]; | |
| 54 |
2/2✓ Branch 0 taken 72 times.
✓ Branch 1 taken 48 times.
|
120 | for (int k = b.row_ptr[col_a]; k < b.row_ptr[col_a + 1]; ++k) { |
| 55 |
2/2✓ Branch 0 taken 44 times.
✓ Branch 1 taken 28 times.
|
72 | int col_b = b.col_ind[k]; |
| 56 |
2/2✓ Branch 0 taken 44 times.
✓ Branch 1 taken 28 times.
|
72 | if (marker[col_b] != i) { |
| 57 |
1/2✓ Branch 0 taken 44 times.
✗ Branch 1 not taken.
|
44 | marker[col_b] = i; |
| 58 | used.push_back(col_b); | ||
| 59 | 44 | acc[col_b] = val_a * b.values[k]; | |
| 60 | } else { | ||
| 61 | 28 | acc[col_b] += val_a * b.values[k]; | |
| 62 | } | ||
| 63 | } | ||
| 64 | } | ||
| 65 | std::ranges::sort(used); | ||
| 66 | 32 | int write_pos = c.row_ptr[i]; | |
| 67 |
2/2✓ Branch 0 taken 44 times.
✓ Branch 1 taken 32 times.
|
76 | for (int col : used) { |
| 68 | 44 | c.values[write_pos] = acc[col]; | |
| 69 | 44 | c.col_ind[write_pos] = col; | |
| 70 | 44 | write_pos++; | |
| 71 | 44 | acc[col] = 0.0; | |
| 72 | } | ||
| 73 | 32 | } | |
| 74 | |||
| 75 | 20 | bool MaslovaUMultMatrTBB::RunImpl() { | |
| 76 | const auto &a = std::get<0>(GetInput()); | ||
| 77 | const auto &b = std::get<1>(GetInput()); | ||
| 78 | auto &c = GetOutput(); | ||
| 79 | |||
| 80 | 20 | c.row_ptr.assign(a.rows + 1, 0); | |
| 81 | |||
| 82 | 20 | tbb::parallel_for(tbb::blocked_range<int>(0, a.rows), [&](const tbb::blocked_range<int> &r) { | |
| 83 | 32 | std::vector<int> marker(b.cols, -1); | |
| 84 |
2/2✓ Branch 0 taken 32 times.
✓ Branch 1 taken 32 times.
|
64 | for (int i = r.begin(); i < r.end(); ++i) { |
| 85 | 32 | c.row_ptr[i + 1] = GetRowNNZ(i, a, b, marker); | |
| 86 | } | ||
| 87 | 32 | }); | |
| 88 | |||
| 89 |
2/2✓ Branch 0 taken 32 times.
✓ Branch 1 taken 20 times.
|
52 | for (int i = 0; i < a.rows; ++i) { |
| 90 | 32 | c.row_ptr[i + 1] += c.row_ptr[i]; | |
| 91 | } | ||
| 92 | |||
| 93 | 20 | c.values.resize(c.row_ptr[a.rows]); | |
| 94 | 20 | c.col_ind.resize(c.row_ptr[a.rows]); | |
| 95 | |||
| 96 | 20 | tbb::parallel_for(tbb::blocked_range<int>(0, a.rows), [&](const tbb::blocked_range<int> &r) { | |
| 97 | 32 | std::vector<double> acc(b.cols, 0.0); | |
| 98 |
1/4✓ Branch 1 taken 32 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
32 | std::vector<int> marker(b.cols, -1); |
| 99 | 32 | std::vector<int> used; | |
| 100 |
1/2✓ Branch 1 taken 32 times.
✗ Branch 2 not taken.
|
32 | used.reserve(b.cols); |
| 101 |
2/2✓ Branch 0 taken 32 times.
✓ Branch 1 taken 32 times.
|
64 | for (int i = r.begin(); i < r.end(); ++i) { |
| 102 |
1/2✓ Branch 1 taken 32 times.
✗ Branch 2 not taken.
|
32 | FillRowValues(i, a, b, c, acc, marker, used); |
| 103 | } | ||
| 104 | 32 | }); | |
| 105 | |||
| 106 | 20 | return true; | |
| 107 | } | ||
| 108 | |||
| 109 | 20 | bool MaslovaUMultMatrTBB::PostProcessingImpl() { | |
| 110 | 20 | return true; | |
| 111 | } | ||
| 112 | |||
| 113 | } // namespace maslova_u_mult_matr_crs | ||
| 114 |