planc
Parallel Lowrank Approximation with Non-negativity Constraints
distntf.cpp
Go to the documentation of this file.
1 /* Copyright 2016 Ramakrishnan Kannan */
2 
3 #include <string>
4 #include "common/distutils.hpp"
6 #include "common/tensor.hpp"
7 #include "common/utils.hpp"
8 #include "distntf/distauntf.hpp"
11 #include "distntf/distntfcpals.hpp"
12 #include "distntf/distntfhals.hpp"
13 #include "distntf/distntfio.hpp"
15 #include "distntf/distntfmu.hpp"
16 #include "distntf/distntfnes.hpp"
17 
18 namespace planc {
19 class DistNTF {
20  private:
21  int m_argc;
22  char **m_argv;
23  int m_k;
24  std::string m_Afile_name;
25  std::string m_outputfile_name;
26  int m_num_it;
27  UVEC m_proc_grids;
28  FVEC m_regs;
29  algotype m_ntfalgo;
30  double m_sparsity;
31  uint m_compute_error;
32  int m_num_k_blocks;
33  UVEC m_global_dims;
34  UVEC m_factor_local_dims;
35  UVEC m_nls_sizes;
36  UVEC m_nls_idxs;
37  bool m_enable_dim_tree;
38  static const int kprimeoffset = 17;
39 
40  void printConfig() {
41  std::cout << "a::" << this->m_ntfalgo << "::i::" << this->m_Afile_name
42  << "::k::" << this->m_k << "::dims::" << this->m_global_dims
43  << "::t::" << this->m_num_it
44  << "::proc_grids::" << this->m_proc_grids
45  << "::error::" << this->m_compute_error
46  << "::regs::" << this->m_regs
47  << "::num_k_blocks::" << m_num_k_blocks
48  << "::dim_tree::" << m_enable_dim_tree << std::endl;
49  }
50 
51  template <class NTFTYPE>
52  void callDistNTF() {
53  planc::Tensor A;
54  std::string rand_prefix("rand_");
55  planc::NTFMPICommunicator mpicomm(this->m_argc, this->m_argv,
56  this->m_proc_grids);
57  mpicomm.printConfig();
58  planc::DistNTFIO dio(mpicomm, A);
59  dio.readInput(m_Afile_name, this->m_global_dims, this->m_proc_grids,
60  this->m_k, this->m_sparsity);
61  this->m_global_dims = dio.global_dims();
62  memusage(mpicomm.rank(), "after input io");
63  INFO << mpicomm.rank()
64  << "::Completed generating tensor A=" << A.dimensions()
65  << "::start indices::" << A.global_idx()
66  << "::global dims::" << this->m_global_dims << std::endl;
67 #ifdef DISTNTF_VERBOSE
68  A.print();
69 #endif
70 #ifdef WRITE_RAND_INPUT
71  dio.writeRandInput();
72 #endif // ifdef WRITE_RAND_INPUT
73  // same matrix as only one of them will be used.
74 
75  // sometimes for really very large matrices starting w/
76  // rand initialization hurts ANLS BPP running time. For a better
77  // initializer we run couple of iterations of HALS.
78 #ifdef MPI_VERBOSE
79  INFO << mpicomm.rank() << "::" << __PRETTY_FUNCTION__
80  << "::" << factors.printinfo() << std::endl;
81 #endif // ifdef MPI_VERBOSE
82  MPI_Barrier(MPI_COMM_WORLD);
83  memusage(mpicomm.rank(), "b4 constructor ");
84  // TODO(ramkikannan): I was here. Need to modify the reallocations by
85  // using localOwnedRowCount instead of m_globalm.
86  // DistAUNTF(const Tensor &i_tensor, const int i_k, algotype i_algo,
87  // const UVEC &i_global_dims,
88  // const UVEC &i_local_dims,
89  // const NTFMPICommunicator &i_mpicomm)
90  // this->m_factor_local_dims = this->m_global_dims / mpicomm.size();
91  this->m_factor_local_dims = A.dimensions();
92  int num_modes = A.dimensions().n_rows;
93 
94  m_nls_sizes = arma::zeros<UVEC>(num_modes);
95  m_nls_idxs = arma::zeros<UVEC>(num_modes);
96  // Calculate NLS sizes
97  for (int i = 0; i < num_modes; i++) {
98  int slice_size = mpicomm.slice_size(i);
99  int slice_rank = mpicomm.slice_rank(i);
100  int num_rows = this->m_factor_local_dims[i];
101  m_nls_sizes[i] = itersplit(num_rows, slice_size, slice_rank);
102  m_nls_idxs[i] = startidx(num_rows, slice_size, slice_rank);
103  }
104 
105  MPI_Barrier(MPI_COMM_WORLD);
106 
107  NTFTYPE ntfsolver(A, this->m_k, this->m_ntfalgo, this->m_global_dims,
108  this->m_factor_local_dims, this->m_nls_sizes,
109  this->m_nls_idxs, mpicomm);
110  memusage(mpicomm.rank(), "after constructor ");
111  ntfsolver.num_iterations(this->m_num_it);
112  ntfsolver.compute_error(this->m_compute_error);
113  if (this->m_enable_dim_tree) {
114  ntfsolver.dim_tree(this->m_enable_dim_tree);
115  }
116  ntfsolver.regularizers(this->m_regs);
117  MPI_Barrier(MPI_COMM_WORLD);
118  // try {
119  mpitic();
120  ntfsolver.computeNTF();
121  double temp = mpitoc();
122  A.clear();
123  if (!this->m_outputfile_name.empty()) {
124  dio.write(this->m_outputfile_name, &ntfsolver);
125  }
126  if (mpicomm.rank() == 0) {
127  printf("NTF took %.3lf secs.\n", temp);
128  }
129  // } catch (std::exception& e) {
130  // printf("Failed rank %d: %s\n", mpicomm.rank(), e.what());
131  // MPI_Abort(MPI_COMM_WORLD, 1);
132  // }
133  }
134 
135  void parseCommandLine() {
136  planc::ParseCommandLine pc(this->m_argc, this->m_argv);
137  pc.parseplancopts();
138  this->m_ntfalgo = pc.lucalgo();
139  this->m_k = pc.lowrankk();
140  this->m_Afile_name = pc.input_file_name();
141  this->m_proc_grids = pc.processor_grids();
142  this->m_sparsity = pc.sparsity();
143  this->m_num_it = pc.iterations();
144  this->m_num_k_blocks = 1;
145  this->m_regs = pc.regularizers();
146  this->m_global_dims = pc.dimensions();
147  this->m_compute_error = pc.compute_error();
148  this->m_enable_dim_tree = pc.dim_tree();
149  this->m_outputfile_name = pc.output_file_name();
150  printConfig();
151  switch (this->m_ntfalgo) {
152  case MU:
153  callDistNTF<DistNTFMU>();
154  break;
155  case HALS:
156  callDistNTF<DistNTFHALS>();
157  break;
158  case ANLSBPP:
159  callDistNTF<DistNTFANLSBPP>();
160  break;
161  case AOADMM:
162  callDistNTF<DistNTFAOADMM>();
163  break;
164  case NESTEROV:
165  callDistNTF<DistNTFNES>();
166  break;
167  case CPALS:
168  callDistNTF<DistNTFCPALS>();
169  break;
170  default:
171  ERR << "Wrong algorithm choice. Quitting.." << this->m_ntfalgo
172  << std::endl;
173  }
174  }
175 
176  public:
178  DistNTF(int argc, char *argv[]) {
179  this->m_argc = argc;
180  this->m_argv = argv;
181  this->parseCommandLine();
182  }
183 };
184 
185 } // namespace planc
186 
187 int main(int argc, char *argv[]) {
188  planc::DistNTF dnd(argc, argv);
189  fflush(stdout);
190  MPI_Barrier(MPI_COMM_WORLD);
191  MPI_Finalize();
192 }
algotype lucalgo()
Returns the NMF algorithm to run. Passed as parameter –algo or -a.
int slice_size(int mode) const
std::string output_file_name()
Returns output file name.
float sparsity()
Input parameter for generating sparse matrix. Passed as -s or –sparsity.
Data is stored such that the unfolding is column major.
Definition: tensor.hpp:32
int rank(const int *i_coords) const
Returns the rank of current MPI process given the cartesian coordinates.
bool compute_error()
Returns whether to compute error not. Passed as parameter -e or –error.
void write(const std::string &output_file_name, DistAUNTF *ntfsolver)
Definition: distntfio.hpp:364
DistNTF(int argc, char *argv[])
Driver function for the distntf.
Definition: distntf.cpp:178
double mpitoc(int rank)
Definition: distutils.hpp:22
int startidx(int n, int p, int r)
Returns the start idx of the current rank r for a global dimension n across p processes.
Definition: distutils.hpp:90
void writeRandInput()
Definition: distntfio.hpp:391
int iterations()
Returns number of iterations. passed as -t or –iter.
#define FVEC
Definition: utils.h:55
Definition: utils.h:10
int main(int argc, char *argv[])
Definition: distntf.cpp:187
algotype
Definition: utils.h:10
#define UVEC
Definition: utils.h:58
int itersplit(int n, int p, int r)
The dimension a particular rank holds out of the global dimension n across p processes.
Definition: distutils.hpp:78
void print() const
prints the value of the tensor.
Definition: tensor.hpp:325
#define ERR
Definition: utils.h:28
void mpitic()
Definition: distutils.hpp:11
void parseplancopts()
parses the command line parameters
#define INFO
Definition: utils.h:36
void memusage(const int myrank, std::string event)
Captures the memory usage of every mpi process.
Definition: distutils.hpp:49
Definition: utils.h:10
UVEC global_idx() const
Definition: tensor.hpp:162
UWORD lowrankk()
returns the low rank. Passed as parameter –lowrank or -k
UVEC processor_grids()
Returns the process grid configuration.
Definition: utils.h:10
UVEC dimensions() const
Returns a vector of dimensions on every mode.
Definition: tensor.hpp:161
std::string input_file_name()
Returns input file name. Passed as -i or –input.
int slice_rank(int i) const
Returns the slice rank on a particular slice grid.
void readInput(const std::string file_name, UVEC i_global_dims, UVEC i_proc_grids, UWORD k=0, double sparsity=0)
Definition: distntfio.hpp:330
FVEC regularizers()
Returns the vector regularizers for all the modes.
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
UVEC global_dims() const
Definition: distntfio.hpp:394
void clear()
Clears the data and also destroys the storage.
Definition: tensor.hpp:153
Definition: utils.h:10
UVEC dimensions()
Returns vector of dimensions for every mode.
Definition: utils.h:10
bool dim_tree()
Enable dimension tree or not.
Definition: utils.h:10