Hi, So I spent some time today investigating a 'bug'. The prototype for nmf is:
template <typename ScalarType> void nmf(viennacl::matrix<ScalarType> const & V, viennacl::matrix<ScalarType> & W, viennacl::matrix<ScalarType> & H, nmf_config const & conf) Now, if W or H is a zero matrix than prod(W, H) is also a zero matrix. So prod(W, H) - V = -V. And ||V|| = ||-V||, so the difference between the initial norm and the norm after the 1st iteration is 0, which is always less than conf.tolerance(). So the routine quits. I think this is a bug in the algorithm: why should the behaviour depend on the initial state of W and H, unless the initial W and H are supposed explicitly to be taken as approximations? But if so, I don't think a zero matrix approximation should cause the algorithm to quit! The following change is sufficient to fix this. I've committed it and a change ensuring any ViennaCL objects created during the run are in the right context. I don't know if it is ideal, since it means that in the (rare) case where W and H in their initial state are adequate, two sets of check_after_steps iterations are computed, where previously only one set would have been. There are alternatives: of course, we could compute norm_frobenius(prod(W, H)) and check that it is not 0 on the first iteration, but that's much more expensive in every other case. diff --git a/viennacl/linalg/nmf.hpp b/viennacl/linalg/nmf.hpp index e47712d..02e4007 100644 --- a/viennacl/linalg/nmf.hpp +++ b/viennacl/linalg/nmf.hpp @@ -173,7 +173,7 @@ namespace viennacl std::cout << diff_val / diff_init << std::endl; // Approximation check - if (diff_val / diff_init < conf.tolerance()) + if ((diff_val / diff_init < conf.tolerance()) && (conf.iters_ != conf.check_after_steps() + 1)) break; // Stagnation check Cheers, Toby -- Toby St Clere Smithe http://tsmithe.net ------------------------------------------------------------------------------ Infragistics Professional Build stunning WinForms apps today! Reboot your WinForms applications with our WinForms controls. Build a bridge from your legacy apps to the future. http://pubads.g.doubleclick.net/gampad/clk?id=153845071&iu=/4140/ostg.clktrk _______________________________________________ ViennaCL-devel mailing list ViennaCL-devel@lists.sourceforge.net https://lists.sourceforge.net/lists/listinfo/viennacl-devel