bnmf-algs
util_details.hpp
Go to the documentation of this file.
1 #pragma once
2 
4 #include "defs.hpp"
5 #include <gsl/gsl_sf_psi.h>
6 #include <vector>
7 
8 namespace bnmf_algs {
9 namespace details {
19 template <typename T, typename Scalar>
20 void check_bld_params(const matrix_t<T>& X, size_t z,
21  const alloc_model::Params<Scalar>& model_params) {
22  BNMF_ASSERT((X.array() >= 0).all(), "X must be nonnegative");
23  BNMF_ASSERT(
24  model_params.alpha.size() == static_cast<size_t>(X.rows()),
25  "Number of alpha parameters must be equal to number of rows of X");
26  BNMF_ASSERT(model_params.beta.size() == static_cast<size_t>(z),
27  "Number of beta parameters must be equal to z");
28 }
29 
38 template <typename T, typename Scalar>
40  const matrix_t<T>& X,
41  const std::vector<alloc_model::Params<Scalar>>& param_vec) {
42 
43  for (long i = 0; i < X.rows(); ++i) {
44  for (long j = 0; j < X.cols(); ++j) {
45  BNMF_ASSERT(std::isnan(X(i, j)) || X(i, j) >= 0,
46  "X must contain nonnegative values or NaN");
47  }
48  }
49 
50  for (const auto& param : param_vec) {
51  BNMF_ASSERT(
52  param.alpha.size() == static_cast<size_t>(X.rows()),
53  "Number of alpha parameters must be equal to number of rows of X");
54  BNMF_ASSERT(param.beta.size() == static_cast<size_t>(X.cols()),
55  "Number of beta parameters must be equal to number of "
56  "columns of X");
57  }
58 }
59 
71 template <typename Real> Real gsl_psi_wrapper(Real x) {
72  return static_cast<Real>(gsl_sf_psi(static_cast<double>(x)));
73 }
74 } // namespace details
75 } // namespace bnmf_algs
Structure to hold the parameters for the Allocation Model .
Definition: alloc_model_params.hpp:25
void check_EM_params(const matrix_t< T > &X, const std::vector< alloc_model::Params< Scalar >> &param_vec)
Do parameter checks on EM algorithms for solving the BLD problem.
Definition: util_details.hpp:39
Real gsl_psi_wrapper(Real x)
A template wrapper function around gsl_sf_psi.
Definition: util_details.hpp:71
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > matrix_t
Matrix type used in the computations.
Definition: defs.hpp:41
std::vector< Scalar > alpha
Parameter vector of Dirichlet prior for matrix of size .
Definition: alloc_model_params.hpp:37
void check_bld_params(const matrix_t< T > &X, size_t z, const alloc_model::Params< Scalar > &model_params)
Do parameter checks on BLD computing function parameters and throw an assertion error if the paramete...
Definition: util_details.hpp:20
T size(T...args)
T isnan(T...args)
std::vector< Scalar > beta
Parameter vector of Dirichlet prior for matrix of size .
Definition: alloc_model_params.hpp:42
Main namespace for bnmf-algs library.
Definition: alloc_model_funcs.hpp:12