ViennaCL - The Vienna Computing Library  1.7.1
Free open-source GPU-accelerated linear algebra and solver library.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
nmf_operations.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_LINALG_OPENCL_NMF_OPERATIONS_HPP_
2 #define VIENNACL_LINALG_OPENCL_NMF_OPERATIONS_HPP_
3 
4 /* =========================================================================
5  Copyright (c) 2010-2016, Institute for Microelectronics,
6  Institute for Analysis and Scientific Computing,
7  TU Wien.
8  Portions of this software are copyright by UChicago Argonne, LLC.
9 
10  -----------------
11  ViennaCL - The Vienna Computing Library
12  -----------------
13 
14  Project Head: Karl Rupp rupp@iue.tuwien.ac.at
15 
16  (A list of authors and contributors can be found in the manual)
17 
18  License: MIT (X11), see file LICENSE in the base directory
19  ============================================================================= */
20 
27 
29 
30 namespace viennacl
31 {
32 namespace linalg
33 {
34 namespace opencl
35 {
36 
44 template<typename NumericT>
48  viennacl::linalg::nmf_config const & conf)
49 {
50  viennacl::ocl::context & ctx = const_cast<viennacl::ocl::context &>(viennacl::traits::opencl_handle(V).context());
51 
52  const std::string NMF_MUL_DIV_KERNEL = "el_wise_mul_div";
53 
55 
56  vcl_size_t k = W.size2();
57  conf.iters_ = 0;
58 
61 
64 
65  viennacl::matrix_base<NumericT> wn(V.size1(), k, W.row_major(), ctx);
66  viennacl::matrix_base<NumericT> wd(V.size1(), k, W.row_major(), ctx);
67  viennacl::matrix_base<NumericT> wtmp(V.size1(), V.size2(), W.row_major(), ctx);
68 
69  viennacl::matrix_base<NumericT> hn(k, V.size2(), H.row_major(), ctx);
70  viennacl::matrix_base<NumericT> hd(k, V.size2(), H.row_major(), ctx);
71  viennacl::matrix_base<NumericT> htmp(k, k, H.row_major(), ctx);
72 
73  viennacl::matrix_base<NumericT> appr(V.size1(), V.size2(), V.row_major(), ctx);
74 
75  NumericT last_diff = 0;
76  NumericT diff_init = 0;
77  bool stagnation_flag = false;
78 
79  for (vcl_size_t i = 0; i < conf.max_iterations(); i++)
80  {
81  conf.iters_ = i + 1;
82  {
83  hn = viennacl::linalg::prod(trans(W), V);
84  htmp = viennacl::linalg::prod(trans(W), W);
85  hd = viennacl::linalg::prod(htmp, H);
86 
88  viennacl::ocl::enqueue(mul_div_kernel(H, hn, hd, cl_uint(H.internal_size1() * H.internal_size2())));
89  }
90  {
91  wn = viennacl::linalg::prod(V, trans(H));
92  wtmp = viennacl::linalg::prod(W, H);
93  wd = viennacl::linalg::prod(wtmp, trans(H));
94 
96 
97  viennacl::ocl::enqueue(mul_div_kernel(W, wn, wd, cl_uint(W.internal_size1() * W.internal_size2())));
98  }
99 
100  if (i % conf.check_after_steps() == 0) //check for convergence
101  {
102  appr = viennacl::linalg::prod(W, H);
103 
104  appr -= V;
106 
107  if (i == 0)
108  diff_init = diff_val;
109 
110  if (conf.print_relative_error())
111  std::cout << diff_val / diff_init << std::endl;
112 
113  // Approximation check
114  if (diff_val / diff_init < conf.tolerance())
115  break;
116 
117  // Stagnation check
118  if (std::fabs(diff_val - last_diff) / (diff_val * NumericT(conf.check_after_steps())) < conf.stagnation_tolerance()) //avoid situations where convergence stagnates
119  {
120  if (stagnation_flag) // iteration stagnates (two iterates with no notable progress)
121  break;
122  else
123  // record stagnation in this iteration
124  stagnation_flag = true;
125  } else
126  // good progress in this iteration, so unset stagnation flag
127  stagnation_flag = false;
128 
129  // prepare for next iterate:
130  last_diff = diff_val;
131  }
132  }
133 }
134 
135 } //namespace opencl
136 } //namespace linalg
137 } //namespace viennacl
138 
139 #endif /* VIENNACL_LINALG_OPENCL_NMF_OPERATIONS_HPP_ */
void trans(const matrix_expression< const matrix_base< NumericT, SizeT, DistanceT >, const matrix_base< NumericT, SizeT, DistanceT >, op_trans > &proxy, matrix_base< NumericT > &temp_trans)
void nmf(viennacl::matrix_base< NumericT > const &V, viennacl::matrix_base< NumericT > &W, viennacl::matrix_base< NumericT > &H, viennacl::linalg::nmf_config const &conf)
The nonnegative matrix factorization (approximation) algorithm as suggested by Lee and Seung...
Represents an OpenCL kernel within ViennaCL.
Definition: kernel.hpp:58
vcl_size_t check_after_steps() const
Number of steps after which the convergence of NMF should be checked (again)
Manages an OpenCL context and provides the respective convenience functions for creating buffers...
Definition: context.hpp:55
Configuration class for the nonnegative-matrix-factorization algorithm. Specify tolerances, maximum iteration counts, etc., here.
float NumericT
Definition: bisect.cpp:40
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
Definition: prod.hpp:102
viennacl::ocl::kernel & get_kernel(std::string const &program_name, std::string const &kernel_name)
Convenience function for retrieving the kernel of a program directly from the context.
Definition: context.hpp:605
vcl_size_t max_iterations() const
Returns the maximum number of iterations for the NMF algorithm.
Represents a vector consisting of scalars 's' only, i.e. v[i] = s for all i. To be used as an initial...
Definition: matrix_def.hpp:93
double tolerance() const
Returns the relative tolerance for convergence.
Main kernel class for generating OpenCL kernels for nonnegative matrix factorization of a dense matri...
Definition: nmf.hpp:58
std::size_t vcl_size_t
Definition: forwards.h:75
size_type size2() const
Returns the number of columns.
Definition: matrix_def.hpp:226
size_type size1() const
Returns the number of rows.
Definition: matrix_def.hpp:224
double stagnation_tolerance() const
Relative tolerance for the stagnation check.
void enqueue(KernelType &k, viennacl::ocl::command_queue const &queue)
Enqueues a kernel in the provided queue.
Definition: enqueue.hpp:50
bool row_major() const
Definition: matrix_def.hpp:248
scalar_expression< const matrix_base< NumericT >, const matrix_base< NumericT >, op_norm_frobenius > norm_frobenius(const matrix_base< NumericT > &A)
OpenCL kernel file for nonnegative matrix factorization.
static void init(viennacl::ocl::context &ctx)
Definition: nmf.hpp:65
size_type internal_size2() const
Returns the internal number of columns. Usually required for launching OpenCL kernels only...
Definition: matrix_def.hpp:240
size_type internal_size1() const
Returns the internal number of rows. Usually required for launching OpenCL kernels only...
Definition: matrix_def.hpp:238
bool print_relative_error() const
Returns the flag specifying whether the relative tolerance should be printed in each iteration...