inline
void
op_reshape::apply(Mat<typename T1::elem_type>& out, const Op<T1,op_reshape>& in)
  {
  arma_extra_debug_sigprint();
  
  typedef typename T1::elem_type eT;
  
  const Proxy<T1> P(in.m);
  
  const uword in_n_rows = in.aux_uword_a;
  const uword in_n_cols = in.aux_uword_b;
  
  if( (is_Mat<typename Proxy<T1>::stored_type>::value == true) && (Proxy<T1>::fake_mat == false) )
    {
    // not checking for aliasing here, as this might be an inplace reshape
    
    const unwrap<typename Proxy<T1>::stored_type> tmp(P.Q);
    
    op_reshape::apply_unwrap(out, tmp.M, in_n_rows, in_n_cols, uword(0));
    }
  else
    {
    if(P.is_alias(out))
      {
      Mat<eT> tmp;
      
      op_reshape::apply_proxy(tmp, P, in_n_rows, in_n_cols);
      
      out.steal_mem(tmp);
      }
    else
      {
      op_reshape::apply_proxy(out, P, in_n_rows, in_n_cols);
      }
    }
  }
arma_hot
inline
void
spop_strans::apply_spmat(SpMat<eT>& out, const SpMat<eT>& X)
  {
  arma_extra_debug_sigprint();
  
  typedef typename umat::elem_type ueT;
  
  const uword N = X.n_nonzero;
  
  if(N == uword(0))
    {
    out.zeros(X.n_cols, X.n_rows);
    return;
    }
  
  umat locs(2, N);
  
  typename SpMat<eT>::const_iterator it = X.begin();
  
  for(uword count = 0; count < N; ++count)
    {
    ueT* locs_ptr = locs.colptr(count);
    
    locs_ptr[0] = it.col();
    locs_ptr[1] = it.row();
    
    ++it;
    }
  
  const Col<eT> vals(const_cast<eT*>(X.values), N, false);
  
  SpMat<eT> tmp(locs, vals, X.n_cols, X.n_rows);
  
  out.steal_mem(tmp);
  }
Exemplo n.º 3
0
inline 
void
op_median::direct_cx_median_index
  (
  uword& out_index1, 
  uword& out_index2, 
  std::vector< arma_cx_median_packet<T> >& X
  )
  {
  arma_extra_debug_sigprint();
  
  typedef arma_cx_median_packet<T> eT;
  
  const uword n_elem = uword(X.size());
  const uword half   = n_elem/2;
  
  typename std::vector<eT>::iterator first    = X.begin();
  typename std::vector<eT>::iterator nth      = first + half;
  typename std::vector<eT>::iterator pastlast = X.end();
  
  std::nth_element(first, nth, pastlast);
  
  out_index1 = (*nth).index;
  
  if((n_elem % 2) == 0)  // even number of elements
    {
    typename std::vector<eT>::iterator start   = X.begin();
    typename std::vector<eT>::iterator pastend = start + half;
    
    out_index2 = (*(std::max_element(start, pastend))).index;
    }
  else  // odd number of elements
    {
    out_index2 = out_index1;
    }
  }
Exemplo n.º 4
0
inline
arma_warn_unused
typename enable_if2< is_arma_sparse_type<T1>::value, typename T1::pod_type >::result
norm
  (
  const T1&   X,
  const uword k = uword(2),
  const typename arma_real_or_cx_only<typename T1::elem_type>::result* junk = 0
  )
  {
  arma_extra_debug_sigprint();
  arma_ignore(junk);
  
  typedef typename T1::elem_type eT;
  typedef typename T1::pod_type   T;
  
  const SpProxy<T1> P(X);
  
  if(P.get_n_nonzero() == 0)
    {
    return T(0);
    }
  
  const bool is_vec = (P.get_n_rows() == 1) || (P.get_n_cols() == 1);
  
  if(is_vec == true)
    {
    const unwrap_spmat<typename SpProxy<T1>::stored_type> tmp(P.Q);
    const SpMat<eT>& A = tmp.M;
    
    // create a fake dense vector to allow reuse of code for dense vectors
    Col<eT> fake_vector( access::rwp(A.values), A.n_nonzero, false );
    
    const Proxy< Col<eT> > P_fake_vector(fake_vector);
    
    switch(k)
      {
      case 1:
        return op_norm::vec_norm_1(P_fake_vector);
        break;
      
      case 2:
        return op_norm::vec_norm_2(P_fake_vector);
        break;
      
      default:
        {
        arma_debug_check( (k == 0), "norm(): k must be greater than zero"   );
        return op_norm::vec_norm_k(P_fake_vector, int(k));
        }
      }
    }
  else
    {
    switch(k)
      {
      case 1:
        return op_norm::mat_norm_1(P);
        break;
      
      case 2:
        return op_norm::mat_norm_2(P);
        break;
      
      default:
        arma_stop_logic_error("norm(): unsupported or unimplemented norm type for sparse matrices");
        return T(0);
      }
    }
  }
inline
void
op_repmat::apply_noalias(Mat<typename obj::elem_type>& out, const obj& X, const uword copies_per_row, const uword copies_per_col)
  {
  arma_extra_debug_sigprint();
  
  typedef typename obj::elem_type eT;
  
  const uword X_n_rows = obj::is_row ? uword(1) : X.n_rows;
  const uword X_n_cols = obj::is_col ? uword(1) : X.n_cols;
  
  out.set_size(X_n_rows * copies_per_row, X_n_cols * copies_per_col);
  
  const uword out_n_rows = out.n_rows;
  const uword out_n_cols = out.n_cols;
  
  // if( (out_n_rows > 0) && (out_n_cols > 0) )
  //   {
  //   for(uword col = 0; col < out_n_cols; col += X_n_cols)
  //   for(uword row = 0; row < out_n_rows; row += X_n_rows)
  //     {
  //     out.submat(row, col, row+X_n_rows-1, col+X_n_cols-1) = X;
  //     }
  //   }
  
  if( (out_n_rows > 0) && (out_n_cols > 0) )
    {
    if(copies_per_row != 1)
      {
      for(uword col_copy=0; col_copy < copies_per_col; ++col_copy)
        {
        const uword out_col_offset = X_n_cols * col_copy;
        
        for(uword col=0; col < X_n_cols; ++col)
          {
                eT* out_colptr = out.colptr(col + out_col_offset);
          const eT* X_colptr   = X.colptr(col);
          
          for(uword row_copy=0; row_copy < copies_per_row; ++row_copy)
            {
            const uword out_row_offset = X_n_rows * row_copy;
            
            arrayops::copy( &out_colptr[out_row_offset], X_colptr, X_n_rows );
            }
          }
        }
      }
    else
      {
      for(uword col_copy=0; col_copy < copies_per_col; ++col_copy)
        {
        const uword out_col_offset = X_n_cols * col_copy;
        
        for(uword col=0; col < X_n_cols; ++col)
          {
                eT* out_colptr = out.colptr(col + out_col_offset);
          const eT* X_colptr   = X.colptr(col);
          
          arrayops::copy( out_colptr, X_colptr, X_n_rows );
          }
        }
      }
    }
  
  }
Exemplo n.º 6
0
arma_hot
inline
void
spglue_minus::apply_noalias(SpMat<eT>& out, const SpProxy<T1>& pa, const SpProxy<T2>& pb)
  {
  arma_extra_debug_sigprint();
  
  arma_debug_assert_same_size(pa.get_n_rows(), pa.get_n_cols(), pb.get_n_rows(), pb.get_n_cols(), "subtraction");
  
  if(pa.get_n_nonzero() == 0)  { out = pb.Q; out *= eT(-1); return; }
  if(pb.get_n_nonzero() == 0)  { out = pa.Q;                return; }
  
  const uword max_n_nonzero = spglue_elem_helper::max_n_nonzero_plus(pa, pb);
  
  // Resize memory to upper bound
  out.reserve(pa.get_n_rows(), pa.get_n_cols(), max_n_nonzero);
  
  // Now iterate across both matrices.
  typename SpProxy<T1>::const_iterator_type x_it  = pa.begin();
  typename SpProxy<T1>::const_iterator_type x_end = pa.end();
  
  typename SpProxy<T2>::const_iterator_type y_it  = pb.begin();
  typename SpProxy<T2>::const_iterator_type y_end = pb.end();
  
  uword count = 0;
  
  while( (x_it != x_end) || (y_it != y_end) )
    {
    eT out_val;
    
    const uword x_it_row = x_it.row();
    const uword x_it_col = x_it.col();
    
    const uword y_it_row = y_it.row();
    const uword y_it_col = y_it.col();
    
    bool use_y_loc = false;
    
    if(x_it == y_it)
      {
      out_val = (*x_it) - (*y_it);
      
      ++x_it;
      ++y_it;
      }
    else
      {
      if((x_it_col < y_it_col) || ((x_it_col == y_it_col) && (x_it_row < y_it_row))) // if y is closer to the end
        {
        out_val = (*x_it);
        
        ++x_it;
        }
      else
        {
        out_val = -(*y_it);  // take the negative
        
        ++y_it;
        
        use_y_loc = true;
        }
      }
    
    if(out_val != eT(0))
      {
      access::rw(out.values[count]) = out_val;
      
      const uword out_row = (use_y_loc == false) ? x_it_row : y_it_row;
      const uword out_col = (use_y_loc == false) ? x_it_col : y_it_col;
      
      access::rw(out.row_indices[count]) = out_row;
      access::rw(out.col_ptrs[out_col + 1])++;
      ++count;
      }
    }
  
  const uword out_n_cols = out.n_cols;
  
  uword* col_ptrs = access::rwp(out.col_ptrs);
  
  // Fix column pointers to be cumulative.
  for(uword c = 1; c <= out_n_cols; ++c)
    {
    col_ptrs[c] += col_ptrs[c - 1];
    }
  
  if(count < max_n_nonzero)
    {
    if(count <= (max_n_nonzero/2))
      {
      out.mem_resize(count);
      }
    else
      {
      // quick resize without reallocating memory and copying data
      access::rw(         out.n_nonzero) = count;
      access::rw(     out.values[count]) = eT(0);
      access::rw(out.row_indices[count]) = uword(0);
      }
    }
  }
Exemplo n.º 7
0
inline
void
op_pinv::apply(Mat<typename T1::elem_type>& out, const Op<T1,op_pinv>& in)
  {
  arma_extra_debug_sigprint();
  
  typedef typename T1::elem_type            eT;
  typedef typename get_pod_type<eT>::result  T;
  
  const bool use_divide_and_conquer = (in.aux_uword_a == 1);
  
  T tol = access::tmp_real(in.aux);
  
  arma_debug_check((tol < T(0)), "pinv(): tolerance must be >= 0");
  
  const Proxy<T1> P(in.m);
  
  const uword n_rows = P.get_n_rows();
  const uword n_cols = P.get_n_cols();
  
  if( (n_rows*n_cols) == 0 )
    {
    out.set_size(n_cols,n_rows);
    return;
    }
  
  
  // economical SVD decomposition 
  Mat<eT> U;
  Col< T> s;
  Mat<eT> V;
  
  bool status = false;
  
  if(use_divide_and_conquer)
    {
    status = (n_cols > n_rows) ? auxlib::svd_dc_econ(U, s, V, trans(P.Q)) : auxlib::svd_dc_econ(U, s, V, P.Q);
    }
  else
    {
    status = (n_cols > n_rows) ? auxlib::svd_econ(U, s, V, trans(P.Q), 'b') : auxlib::svd_econ(U, s, V, P.Q, 'b');
    }
  
  if(status == false)
    {
    out.reset();
    arma_bad("pinv(): svd failed");
    return;
    }
  
  const uword s_n_elem = s.n_elem;
  const T*    s_mem    = s.memptr();
  
  // set tolerance to default if it hasn't been specified
  if( (tol == T(0)) && (s_n_elem > 0) )
    {
    tol = (std::max)(n_rows, n_cols) * s_mem[0] * std::numeric_limits<T>::epsilon();
    }
  
  
  uword count = 0;
  
  for(uword i = 0; i < s_n_elem; ++i)
    {
    count += (s_mem[i] >= tol) ? uword(1) : uword(0);
    }
  
  
  if(count > 0)
    {
    Col<T> s2(count);
    
    T* s2_mem = s2.memptr();
    
    uword count2 = 0;
    
    for(uword i=0; i < s_n_elem; ++i)
      {
      const T val = s_mem[i];
      
      if(val >= tol)  {  s2_mem[count2] = T(1) / val;  ++count2; }
      }
    
    
    if(n_rows >= n_cols)
      {
      out = ( (V.n_cols > count) ? V.cols(0,count-1) : V ) * diagmat(s2) * trans( (U.n_cols > count) ? U.cols(0,count-1) : U );
      }
    else
      {
      out = ( (U.n_cols > count) ? U.cols(0,count-1) : U ) * diagmat(s2) * trans( (V.n_cols > count) ? V.cols(0,count-1) : V );
      }
    }
  else
    {
    out.zeros(n_cols, n_rows);
    }
  }
Exemplo n.º 8
0
inline
bool
op_all::all_vec_helper
  (
  const mtOp<uword, T1, op_type>& X,
  const typename arma_op_rel_only<op_type>::result           junk1,
  const typename arma_not_cx<typename T1::elem_type>::result junk2
  )
  {
  arma_extra_debug_sigprint();
  arma_ignore(junk1);
  arma_ignore(junk2);
  
  typedef typename T1::elem_type eT;
  
  const eT val = X.aux;
  
  const Proxy<T1> P(X.m);
  
  const uword n_elem = P.get_n_elem();
  
  uword count = 0;
  
  if(Proxy<T1>::prefer_at_accessor == false)
    {
    typename Proxy<T1>::ea_type Pea = P.get_ea();
    
    for(uword i=0; i < n_elem; ++i)
      {
      const eT tmp = Pea[i];
      
           if(is_same_type<op_type, op_rel_lt_pre   >::yes)  { count += (val <  tmp) ? uword(1) : uword(0); }
      else if(is_same_type<op_type, op_rel_lt_post  >::yes)  { count += (tmp <  val) ? uword(1) : uword(0); }
      else if(is_same_type<op_type, op_rel_gt_pre   >::yes)  { count += (val >  tmp) ? uword(1) : uword(0); }
      else if(is_same_type<op_type, op_rel_gt_post  >::yes)  { count += (tmp >  val) ? uword(1) : uword(0); }
      else if(is_same_type<op_type, op_rel_lteq_pre >::yes)  { count += (val <= tmp) ? uword(1) : uword(0); }
      else if(is_same_type<op_type, op_rel_lteq_post>::yes)  { count += (tmp <= val) ? uword(1) : uword(0); }
      else if(is_same_type<op_type, op_rel_gteq_pre >::yes)  { count += (val >= tmp) ? uword(1) : uword(0); }
      else if(is_same_type<op_type, op_rel_gteq_post>::yes)  { count += (tmp >= val) ? uword(1) : uword(0); }
      else if(is_same_type<op_type, op_rel_eq       >::yes)  { count += (tmp == val) ? uword(1) : uword(0); }
      else if(is_same_type<op_type, op_rel_noteq    >::yes)  { count += (tmp != val) ? uword(1) : uword(0); }
      }
    }
  else
    {
    const uword n_rows = P.get_n_rows();
    const uword n_cols = P.get_n_cols();
    
    for(uword col=0; col < n_cols; ++col)
    for(uword row=0; row < n_rows; ++row)
      {
      const eT tmp = P.at(row,col);
      
           if(is_same_type<op_type, op_rel_lt_pre   >::yes)  { if(val <  tmp) { ++count; } }
      else if(is_same_type<op_type, op_rel_lt_post  >::yes)  { if(tmp <  val) { ++count; } }
      else if(is_same_type<op_type, op_rel_gt_pre   >::yes)  { if(val >  tmp) { ++count; } }
      else if(is_same_type<op_type, op_rel_gt_post  >::yes)  { if(tmp >  val) { ++count; } }
      else if(is_same_type<op_type, op_rel_lteq_pre >::yes)  { if(val <= tmp) { ++count; } }
      else if(is_same_type<op_type, op_rel_lteq_post>::yes)  { if(tmp <= val) { ++count; } }
      else if(is_same_type<op_type, op_rel_gteq_pre >::yes)  { if(val >= tmp) { ++count; } }
      else if(is_same_type<op_type, op_rel_gteq_post>::yes)  { if(tmp >= val) { ++count; } }
      else if(is_same_type<op_type, op_rel_eq       >::yes)  { if(tmp == val) { ++count; } }
      else if(is_same_type<op_type, op_rel_noteq    >::yes)  { if(tmp != val) { ++count; } }
      }
    }
  
  return (n_elem == count);
  }
Exemplo n.º 9
0
void DHT::temperature (uint16_t * ipart, uint16_t * dpart) const {
	decode_dht22 (uword (data_[kTl] & 0x7F, data_[kTh]), ipart, dpart);
}
Exemplo n.º 10
0
arma_inline
typename enable_if2< (is_supported_blas_type<typename T1::elem_type>::value && is_cx<typename T1::elem_type>::no), const mtOp<std::complex<typename T1::elem_type>, T1, op_logmat> >::result
logmat(const Base<typename T1::elem_type,T1>& X, const uword n_iters = 100u)
  {
  arma_extra_debug_sigprint();
  
  return mtOp<std::complex<typename T1::elem_type>, T1, op_logmat>(X.get_ref(), n_iters, uword(0));
  }
Exemplo n.º 11
0
inline
void
op_any::apply_helper(Mat<uword>& out, const Proxy<T1>& P, const uword dim)
  {
  arma_extra_debug_sigprint();
  
  const uword n_rows = P.get_n_rows();
  const uword n_cols = P.get_n_cols();
  
  typedef typename Proxy<T1>::elem_type eT;
  
  if(dim == 0)  // traverse rows (ie. process each column)
    {
    out.zeros(1, n_cols);
    
    uword* out_mem = out.memptr();
    
    if(is_Mat<typename Proxy<T1>::stored_type>::value == true)
      {
      const unwrap<typename Proxy<T1>::stored_type> U(P.Q);
      
      for(uword col=0; col < n_cols; ++col)
        {
        const eT* colmem = U.M.colptr(col);
        
        for(uword row=0; row < n_rows; ++row)
          {
          if(colmem[row] != eT(0))  { out_mem[col] = uword(1); break; }
          }
        }
      }
    else
      {
      for(uword col=0; col < n_cols; ++col)
        {
        for(uword row=0; row < n_rows; ++row)
          {
          if(P.at(row,col) != eT(0))  { out_mem[col] = uword(1); break; }
          }
        }
      }
    }
  else
    {
    out.zeros(n_rows, 1);
    
    uword* out_mem = out.memptr();
    
    if(is_Mat<typename Proxy<T1>::stored_type>::value == true)
      {
      const unwrap<typename Proxy<T1>::stored_type> U(P.Q);
      
      for(uword col=0; col < n_cols; ++col)
        {
        const eT* colmem = U.M.colptr(col);
        
        for(uword row=0; row < n_rows; ++row)
          {
          if(colmem[row] != eT(0))  { out_mem[row] = uword(1); }
          }
        }
      }
    else
      {
      for(uword col=0; col < n_cols; ++col)
        {
        for(uword row=0; row < n_rows; ++row)
          {
          if(P.at(row,col) != eT(0))  { out_mem[row] = uword(1); }
          }
        }
      }
    }
  }
Exemplo n.º 12
0
inline
void
glue_histc::apply_noalias(Mat<uword>& C, const Mat<eT>& A, const Mat<eT>& B, const uword dim)
  {
  arma_extra_debug_sigprint();
  
  arma_debug_check( ((B.is_vec() == false) && (B.is_empty() == false)), "histc(): parameter 'edges' is not a vector" );
  
  const uword A_n_rows = A.n_rows;
  const uword A_n_cols = A.n_cols;
  
  const uword B_n_elem = B.n_elem;
  
  if( B_n_elem == uword(0) )  { C.reset(); return; }
  
  const eT*   B_mem       = B.memptr();
  const uword B_n_elem_m1 = B_n_elem - 1;
  
  if(dim == uword(0))
    {
    C.zeros(B_n_elem, A_n_cols);
    
    for(uword col=0; col < A_n_cols; ++col)
      {
      const eT*    A_coldata = A.colptr(col);
            uword* C_coldata = C.colptr(col);
      
      for(uword row=0; row < A_n_rows; ++row)
        {
        const eT x = A_coldata[row];
        
        for(uword i=0; i < B_n_elem_m1; ++i)
          {
               if( (B_mem[i]           <= x) && (x < B_mem[i+1]) )  { C_coldata[i]++;           break; }
          else if(  B_mem[B_n_elem_m1] == x                      )  { C_coldata[B_n_elem_m1]++; break; }    // for compatibility with Matlab
          }
        }
      }
    }
  else
  if(dim == uword(1))
    {
    C.zeros(A_n_rows, B_n_elem);
    
    if(A.n_rows == 1)
      {
      const uword  A_n_elem = A.n_elem;
      const eT*    A_mem    = A.memptr();
            uword* C_mem    = C.memptr();
      
      for(uword j=0; j < A_n_elem; ++j)
        {
        const eT x = A_mem[j];
        
        for(uword i=0; i < B_n_elem_m1; ++i)
          {
               if( (B_mem[i]           <= x) && (x < B_mem[i+1]) )  { C_mem[i]++;           break; }
          else if(  B_mem[B_n_elem_m1] == x                      )  { C_mem[B_n_elem_m1]++; break; }    // for compatibility with Matlab
          }
        }
      }
    else
      {
      for(uword row=0; row < A_n_rows; ++row)
      for(uword col=0; col < A_n_cols; ++col)
        {
        const eT x = A.at(row,col);
        
        for(uword i=0; i < B_n_elem_m1; ++i)
          {
               if( (B_mem[i]            <= x) && (x < B_mem[i+1]) )  { C.at(row,i)++;           break; }
          else if(  B_mem[B_n_elem_m1]  == x                      )  { C.at(row,B_n_elem_m1)++; break; }   // for compatibility with Matlab
          }
        }
      }
    }
  }
Exemplo n.º 13
0
inline
void
op_all::apply_helper(Mat<uword>& out, const Proxy<T1>& P, const uword dim)
  {
  arma_extra_debug_sigprint();
  
  const uword n_rows = P.get_n_rows();
  const uword n_cols = P.get_n_cols();
  
  typedef typename Proxy<T1>::elem_type eT;
  
  if(dim == 0)  // traverse rows (ie. process each column)
    {
    out.zeros(1, n_cols);
    
    if(out.n_elem == 0)  { return; }
    
    uword* out_mem = out.memptr();
    
    if(is_Mat<typename Proxy<T1>::stored_type>::value == true)
      {
      const unwrap<typename Proxy<T1>::stored_type> U(P.Q);
      
      for(uword col=0; col < n_cols; ++col)
        {
        const eT* colmem = U.M.colptr(col);
        
        uword count = 0;
        
        for(uword row=0; row < n_rows; ++row)
          {
          if(colmem[row] != eT(0))  { ++count; }
          }
        
        out_mem[col] = (n_rows == count) ? uword(1) : uword(0); 
        }
      }
    else
      {
      for(uword col=0; col < n_cols; ++col)
        {
        uword count = 0;
        
        for(uword row=0; row < n_rows; ++row)
          {
          if(P.at(row,col) != eT(0))  { ++count; }
          }
        
        out_mem[col] = (n_rows == count) ? uword(1) : uword(0); 
        }
      }
    }
  else
    {
    out.zeros(n_rows, 1);
    
    uword* out_mem = out.memptr();
    
    // internal dual use of 'out': keep the counts for each row
    
    if(is_Mat<typename Proxy<T1>::stored_type>::value == true)
      {
      const unwrap<typename Proxy<T1>::stored_type> U(P.Q);
      
      for(uword col=0; col < n_cols; ++col)
        {
        const eT* colmem = U.M.colptr(col);
        
        for(uword row=0; row < n_rows; ++row)
          {
          if(colmem[row] != eT(0))  { ++out_mem[row]; }
          }
        }
      }
    else
      {
      for(uword col=0; col < n_cols; ++col)
        {
        for(uword row=0; row < n_rows; ++row)
          {
          if(P.at(row,col) != eT(0))  { ++out_mem[row]; }
          }
        }
      }
    
    
    // see what the counts tell us
    
    for(uword row=0; row < n_rows; ++row)
      {
      out_mem[row] = (n_cols == out_mem[row]) ? uword(1) : uword(0);
      }
    
    }
  }
Exemplo n.º 14
0
inline
arma_warn_unused
typename enable_if2< is_arma_type<T1>::value, typename T1::pod_type >::result
norm
  (
  const T1&   X,
  const uword k = uword(2),
  const typename arma_real_or_cx_only<typename T1::elem_type>::result* junk = 0
  )
  {
  arma_extra_debug_sigprint();
  arma_ignore(junk);
  
  typedef typename T1::pod_type T;
  
  const Proxy<T1> P(X);
  
  if(P.get_n_elem() == 0)
    {
    return T(0);
    }
  
  const bool is_vec = (T1::is_row) || (T1::is_col) || (P.get_n_rows() == 1) || (P.get_n_cols() == 1);
  
  if(is_vec)
    {
    switch(k)
      {
      case 1:
        return op_norm::vec_norm_1(P);
        break;
      
      case 2:
        return op_norm::vec_norm_2(P);
        break;
      
      default:
        {
        arma_debug_check( (k == 0), "norm(): k must be greater than zero" );
        return op_norm::vec_norm_k(P, int(k));
        }
      }
    }
  else
    {
    switch(k)
      {
      case 1:
        return op_norm::mat_norm_1(P);
        break;
      
      case 2:
        return op_norm::mat_norm_2(P);
        break;
      
      default:
        arma_stop_logic_error("norm(): unsupported matrix norm type");
        return T(0);
      }
    }
  
  return T(0);  // prevent erroneous compiler warnings
  }
Exemplo n.º 15
0
inline
bool
op_chol::apply_direct(Mat<typename T1::elem_type>& out, const Base<typename T1::elem_type,T1>& A_expr, const uword layout)
  {
  arma_extra_debug_sigprint();
  
  out = A_expr.get_ref();
  
  arma_debug_check( (out.is_square() == false), "chol(): given matrix must be square sized" );
  
  if(out.is_empty())  { return true; }
  
  uword KD = 0;
  
  const bool is_band = (auxlib::crippled_lapack(out)) ? false : ((layout == 0) ? band_helper::is_band_upper(KD, out, uword(32)) : band_helper::is_band_lower(KD, out, uword(32)));
  
  const bool status = (is_band) ? auxlib::chol_band(out, KD, layout) : auxlib::chol(out, layout);
  
  return status;
  }
Exemplo n.º 16
0
inline
const SpSubview<eT>&
SpSubview<eT>::operator=(const Base<eT, T1>& in)
  {
  arma_extra_debug_sigprint();
  
  // this is a modified version of SpSubview::operator_equ_common(const SpBase)
  
  const SpProxy< SpMat<eT> > pa((*this).m);
  
  const unwrap<T1>     b_tmp(in.get_ref());
  const Mat<eT>&   b = b_tmp.M;
  
  arma_debug_assert_same_size(n_rows, n_cols, b.n_rows, b.n_cols, "insertion into sparse submatrix");
  
  const uword pa_start_row = (*this).aux_row1;
  const uword pa_start_col = (*this).aux_col1;
  
  const uword pa_end_row = pa_start_row + (*this).n_rows - 1;
  const uword pa_end_col = pa_start_col + (*this).n_cols - 1;
  
  const uword pa_n_rows = pa.get_n_rows();
  
  const uword b_n_elem = b.n_elem;
  const eT*   b_mem    = b.memptr();
  
  uword box_count = 0;
  
  for(uword i=0; i<b_n_elem; ++i)
    {
    box_count += (b_mem[i] != eT(0)) ? uword(1) : uword(0);
    }
  
  SpMat<eT> out(pa.get_n_rows(), pa.get_n_cols());
  
  const uword alt_count = pa.get_n_nonzero() - (*this).n_nonzero + box_count;
  
  // Resize memory to correct size.
  out.mem_resize(alt_count);
  
  typename SpProxy< SpMat<eT> >::const_iterator_type x_it  = pa.begin();
  typename SpProxy< SpMat<eT> >::const_iterator_type x_end = pa.end();
  
  uword b_row = 0;
  uword b_col = 0;
    
  bool x_it_ok = (x_it != x_end);
  bool y_it_ok = ( (b_row < b.n_rows) && (b_col < b.n_cols) );
  
  uword x_it_row = (x_it_ok) ? x_it.row() : 0;
  uword x_it_col = (x_it_ok) ? x_it.col() : 0;
  
  uword y_it_row = (y_it_ok) ? b_row + pa_start_row : 0;
  uword y_it_col = (y_it_ok) ? b_col + pa_start_col : 0;
    
  uword cur_val = 0;
  while(x_it_ok || y_it_ok)
    {
    const bool x_inside_box = (x_it_row >= pa_start_row) && (x_it_row <= pa_end_row) && (x_it_col >= pa_start_col) && (x_it_col <= pa_end_col);
    const bool y_inside_box = (y_it_row >= pa_start_row) && (y_it_row <= pa_end_row) && (y_it_col >= pa_start_col) && (y_it_col <= pa_end_col);
    
    const eT x_val = x_inside_box ? eT(0) : ( x_it_ok ? (*x_it) : eT(0) );
    
    const eT y_val = y_inside_box ? ( y_it_ok ? b.at(b_row,b_col) : eT(0) ) : eT(0);
    
    if( (x_it_row == y_it_row) && (x_it_col == y_it_col) )
      {
      if( (x_val != eT(0)) || (y_val != eT(0)) )  
        {
        access::rw(out.values[cur_val]) = (x_val != eT(0)) ? x_val : y_val;
        access::rw(out.row_indices[cur_val]) = x_it_row;
        ++access::rw(out.col_ptrs[x_it_col + 1]);
        ++cur_val;
        }
      
      if(x_it_ok)
        {
        ++x_it;
        
        if(x_it == x_end)  { x_it_ok = false; }
        }
      
      if(x_it_ok)
        {
        x_it_row = x_it.row();
        x_it_col = x_it.col();
        }
      else
        {
        x_it_row++;
        
        if(x_it_row >= pa_n_rows)  { x_it_row = 0; x_it_col++; }
        }
      
      if(y_it_ok)
        {
        b_row++;
        
        if(b_row >= b.n_rows)  { b_row = 0; b_col++; }
        
        if( (b_row > b.n_rows) || (b_col > b.n_cols) )  { y_it_ok = false; }
        }
      
      if(y_it_ok)
        {
        y_it_row = b_row + pa_start_row;
        y_it_col = b_col + pa_start_col;
        }
      else
        {
        y_it_row++;
        
        if(y_it_row >= pa_n_rows)  { y_it_row = 0; y_it_col++; }
        }
      }
    else
      {
      if((x_it_col < y_it_col) || ((x_it_col == y_it_col) && (x_it_row < y_it_row))) // if y is closer to the end
        {
        if(x_val != eT(0))
          {
          access::rw(out.values[cur_val]) = x_val;
          access::rw(out.row_indices[cur_val]) = x_it_row;
          ++access::rw(out.col_ptrs[x_it_col + 1]);
          ++cur_val;
          }
        
        if(x_it_ok)
          {
          ++x_it;
          
          if(x_it == x_end)  { x_it_ok = false; }
          }
        
        if(x_it_ok)
          {
          x_it_row = x_it.row();
          x_it_col = x_it.col();
          }
        else
          {
          x_it_row++;
          
          if(x_it_row >= pa_n_rows)  { x_it_row = 0; x_it_col++; }
          }
        }
      else
        {
        if(y_val != eT(0))
          {
          access::rw(out.values[cur_val]) = y_val;
          access::rw(out.row_indices[cur_val]) = y_it_row;
          ++access::rw(out.col_ptrs[y_it_col + 1]);
          ++cur_val;
          }
        
        if(y_it_ok)
          {
          b_row++;
          
          if(b_row >= b.n_rows)  { b_row = 0; b_col++; }
          
          if( (b_row > b.n_rows) || (b_col > b.n_cols) )  { y_it_ok = false; }
          }
        
        if(y_it_ok)
          {
          y_it_row = b_row + pa_start_row;
          y_it_col = b_col + pa_start_col;
          }
        else
          {
          y_it_row++;
          
          if(y_it_row >= pa_n_rows)  { y_it_row = 0; y_it_col++; }
          }
        }
      }
    }
  
  const uword out_n_cols = out.n_cols;
  
  uword* col_ptrs = access::rwp(out.col_ptrs);
  
  // Fix column pointers to be cumulative.
  for(uword c = 1; c <= out_n_cols; ++c)
    {
    col_ptrs[c] += col_ptrs[c - 1];
    }
  
  access::rw((*this).m).steal_mem(out);
  
  access::rw(n_nonzero) = box_count;
  
  return *this;
  }
Exemplo n.º 17
0
 arma_inline
 static
 bool
 eval(const uword n_elem)
   {
   #if defined(ARMA_USE_OPENMP)
     {
     const bool length_ok = (is_cx<eT>::yes || use_smaller_thresh) ? (n_elem >= (arma_config::mp_threshold/uword(2))) : (n_elem >= arma_config::mp_threshold);
     
     if(length_ok)
       {
       if(omp_in_parallel())  { return false; }
       }
     
     return length_ok;
     }
   #else
     {
     arma_ignore(n_elem);
     
     return false;
     }
   #endif
   }
Exemplo n.º 18
0
void DHT::humidity (uint16_t * ipart, uint16_t * dpart) const {
	decode_dht22 (uword (data_[kHl], data_[kHh]), ipart, dpart);
}
Exemplo n.º 19
0
inline
bool
op_all::all_vec_helper
  (
  const mtGlue<uword, T1, T2, glue_type>& X,
  const typename arma_glue_rel_only<glue_type>::result       junk1,
  const typename arma_not_cx<typename T1::elem_type>::result junk2,
  const typename arma_not_cx<typename T2::elem_type>::result junk3
  )
  {
  arma_extra_debug_sigprint();
  arma_ignore(junk1);
  arma_ignore(junk2);
  arma_ignore(junk3);
  
  typedef typename T1::elem_type eT1;
  typedef typename T2::elem_type eT2;
  
  typedef typename Proxy<T1>::ea_type ea_type1;
  typedef typename Proxy<T2>::ea_type ea_type2;
  
  const Proxy<T1> A(X.A);
  const Proxy<T2> B(X.B);
  
  arma_debug_assert_same_size(A, B, "relational operator");
  
  const uword n_elem = A.get_n_elem();
  
  uword count = 0;
  
  const bool prefer_at_accessor = Proxy<T1>::prefer_at_accessor || Proxy<T2>::prefer_at_accessor;
  
  if(prefer_at_accessor == false)
    {
    ea_type1 PA = A.get_ea();
    ea_type2 PB = B.get_ea();
    
    for(uword i=0; i<n_elem; ++i)
      {
      const eT1 tmp1 = PA[i];
      const eT2 tmp2 = PB[i];
      
           if(is_same_type<glue_type, glue_rel_lt    >::yes)  { count += (tmp1 <  tmp2) ? uword(1) : uword(0); }
      else if(is_same_type<glue_type, glue_rel_gt    >::yes)  { count += (tmp1 >  tmp2) ? uword(1) : uword(0); }
      else if(is_same_type<glue_type, glue_rel_lteq  >::yes)  { count += (tmp1 <= tmp2) ? uword(1) : uword(0); }
      else if(is_same_type<glue_type, glue_rel_gteq  >::yes)  { count += (tmp1 >= tmp2) ? uword(1) : uword(0); }
      else if(is_same_type<glue_type, glue_rel_eq    >::yes)  { count += (tmp1 == tmp2) ? uword(1) : uword(0); }
      else if(is_same_type<glue_type, glue_rel_noteq >::yes)  { count += (tmp1 != tmp2) ? uword(1) : uword(0); }
      else if(is_same_type<glue_type, glue_rel_and   >::yes)  { count += (tmp1 && tmp2) ? uword(1) : uword(0); }
      else if(is_same_type<glue_type, glue_rel_or    >::yes)  { count += (tmp1 || tmp2) ? uword(1) : uword(0); }
      }
    }
  else
    {
    const uword n_rows = A.get_n_rows();
    const uword n_cols = A.get_n_cols();
    
    for(uword col=0; col < n_cols; ++col)
    for(uword row=0; row < n_rows; ++row)
      {
      const eT1 tmp1 = A.at(row,col);
      const eT2 tmp2 = B.at(row,col);
      
           if(is_same_type<glue_type, glue_rel_lt    >::yes)  { if(tmp1 <  tmp2) { ++count; } }
      else if(is_same_type<glue_type, glue_rel_gt    >::yes)  { if(tmp1 >  tmp2) { ++count; } }
      else if(is_same_type<glue_type, glue_rel_lteq  >::yes)  { if(tmp1 <= tmp2) { ++count; } }
      else if(is_same_type<glue_type, glue_rel_gteq  >::yes)  { if(tmp1 >= tmp2) { ++count; } }
      else if(is_same_type<glue_type, glue_rel_eq    >::yes)  { if(tmp1 == tmp2) { ++count; } }
      else if(is_same_type<glue_type, glue_rel_noteq >::yes)  { if(tmp1 != tmp2) { ++count; } }
      else if(is_same_type<glue_type, glue_rel_and   >::yes)  { if(tmp1 && tmp2) { ++count; } }
      else if(is_same_type<glue_type, glue_rel_or    >::yes)  { if(tmp1 || tmp2) { ++count; } }
      }
    }
  
  return (n_elem == count);
  }