planc
Parallel Lowrank Approximation with Non-negativity Constraints
distntfanlsbpp.hpp
Go to the documentation of this file.
1 /* Copyright Ramakrishnan Kannan 2018 */
2 
3 #ifndef DISTNTF_DISTNTFANLSBPP_HPP_
4 #define DISTNTF_DISTNTFANLSBPP_HPP_
5 
6 #include "distntf/distauntf.hpp"
7 #include "nnls/bppnnls.hpp"
8 
9 namespace planc {
10 
11 #define ONE_THREAD_MATRIX_SIZE 2000
12 
13 class DistNTFANLSBPP : public DistAUNTF {
14  protected:
22  MAT update(const int mode) {
23  MAT othermat(this->m_local_ncp_factors_t.factor(mode));
24  if (m_nls_sizes[mode] > 0) {
25  unsigned int numThreads =
26  (this->ncp_local_mttkrp_t[mode].n_cols / ONE_THREAD_MATRIX_SIZE) + 1;
27  #pragma omp parallel for schedule(dynamic)
28  for (UINT i = 0; i < numThreads; i++) {
29  UINT spanStart = i * ONE_THREAD_MATRIX_SIZE;
30  UINT spanEnd = (i + 1) * ONE_THREAD_MATRIX_SIZE - 1;
31  if (spanEnd > this->ncp_local_mttkrp_t[mode].n_cols - 1) {
32  spanEnd = this->ncp_local_mttkrp_t[mode].n_cols - 1;
33  }
34  // if it is exactly divisible, the last iteration is unnecessary.
35  BPPNNLS<MAT, VEC> *subProblem;
36  if (spanStart <= spanEnd) {
37  if (spanStart == spanEnd) {
38  subProblem = new BPPNNLS<MAT, VEC>(
39  this->global_gram,
40  (VEC)this->ncp_local_mttkrp_t[mode].col(spanStart), true);
41  } else { // if (spanStart < spanEnd)
42  subProblem = new BPPNNLS<MAT, VEC>(
43  this->global_gram,
44  (MAT)this->ncp_local_mttkrp_t[mode].cols(spanStart, spanEnd),
45  true);
46  }
47 #ifdef _VERBOSE
48  INFO << "Scheduling " << worh << " start=" << spanStart
49  << ", end=" << spanEnd << ", tid=" << omp_get_thread_num()
50  << std::endl;
51 #endif
52  // tic();
53  subProblem->solveNNLS();
54  // t2 = toc();
55 #ifdef _VERBOSE
56  INFO << "completed " << worh << " start=" << spanStart
57  << ", end=" << spanEnd << ", tid=" << omp_get_thread_num()
58  << " cpu=" << sched_getcpu() << " time taken=" << t2
59  << " num_iterations()=" << numIter << std::endl;
60 #endif
61  if (spanStart == spanEnd) {
62  VEC solVec = subProblem->getSolutionVector();
63  othermat.col(i) = solVec;
64  } else { // if (spanStart < spanEnd)
65  othermat.cols(spanStart, spanEnd) = subProblem->getSolutionMatrix();
66  }
67  subProblem->clear();
68  delete subProblem;
69  }
70  }
71  } else {
72  othermat.zeros();
73  }
74  return othermat;
75  }
76 
77  public:
78  DistNTFANLSBPP(const Tensor &i_tensor, const int i_k, algotype i_algo,
79  const UVEC &i_global_dims, const UVEC &i_local_dims,
80  const UVEC &i_nls_sizes, const UVEC &i_nls_idxs,
81  const NTFMPICommunicator &i_mpicomm)
82  : DistAUNTF(i_tensor, i_k, i_algo, i_global_dims, i_local_dims,
83  i_nls_sizes, i_nls_idxs, i_mpicomm) {}
84 }; // class DistNTFANLSBPP
85 
86 } // namespace planc
87 
88 #endif // DISTNTF_DISTNTFANLSBPP_HPP_
MATTYPE getSolutionMatrix()
Definition: nnls.hpp:79
Data is stored such that the unfolding is column major.
Definition: tensor.hpp:32
#define ONE_THREAD_MATRIX_SIZE
algotype
Definition: utils.h:10
#define UVEC
Definition: utils.h:58
#define INFO
Definition: utils.h:36
int solveNNLS()
Definition: bppnnls.hpp:30
DistNTFANLSBPP(const Tensor &i_tensor, const int i_k, algotype i_algo, const UVEC &i_global_dims, const UVEC &i_local_dims, const UVEC &i_nls_sizes, const UVEC &i_nls_idxs, const NTFMPICommunicator &i_mpicomm)
VECTYPE getSolutionVector()
Definition: nnls.hpp:76
unsigned int UINT
Definition: utils.h:68
#define MAT
Definition: utils.h:52
void clear()
Definition: nnls.hpp:82
ncp_factors contains the factors of the ncp every ith factor is of size n_i * k number of factors is ...
Definition: ncpfactors.hpp:20
#define VEC
Definition: utils.h:61