1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141
|
/*
Copyright (C) 2014 Matteo Fasiolo matteo.fasiolo@gmail.com
This program is free software; you can redistribute it and/or
modify it under the terms of the GNU General Public License
as published by the Free Software Foundation; either version 2
of the License, or (at your option) any later version.
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
(www.gnu.org/copyleft/gpl.html)
You should have received a copy of the GNU General Public License
along with this program; if not, write to the Free Software
Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307,
USA. */
#include "mvnfast.h"
/*
* Fast computation of Mahalanobis distance
*
* See ?maha() for a description of the arguments and output.
*/
/*
* Interface to R
*/
RcppExport SEXP mahaCpp(SEXP X, SEXP mu, SEXP sigma, SEXP ncores, SEXP isChol)
{
using namespace Rcpp;
try{
arma::mat X_ = as<arma::mat>(X);
arma::vec mu_ = as<arma::vec>(mu);
arma::mat sigma_ = as<arma::mat>(sigma);
unsigned int ncores_ = as<unsigned int>(ncores);
bool isChol_ = as<bool>(isChol);
NumericVector dist = wrap( mahaInt(X_, mu_, sigma_, ncores_, isChol_) );
dist.attr( "dim" ) = R_NilValue;
return dist;
} catch( std::exception& __ex__){
forward_exception_to_r(__ex__);
} catch(...){
::Rf_error( "c++ exception (unknown reason)" );
}
return wrap(NA_REAL);
}
/*
* Internal C++ function for Mahalanobis distance
*/
arma::vec mahaInt(arma::mat & X,
arma::vec & mu,
arma::mat & sigma,
const unsigned int ncores,
const bool isChol = false)
{
using namespace arma;
// Some sanity checks
if(ncores == 0) Rcpp::stop("ncores has to be positive.");
if(mu.n_elem != sigma.n_cols) Rcpp::stop("The mean vector has a different dimensions from the covariance matrix.");
if(X.n_cols != sigma.n_cols) Rcpp::stop("The number of columns of X is different from the dimension of the covariance matrix.");
// Calculate transposed cholesky dec. unless sigma is alread a cholesky dec.
mat cholDec;
if( isChol == false ) {
cholDec = trimatl(chol(sigma).t());
}
else{
cholDec = trimatl(sigma.t());
if(any(cholDec.diag() <= 0.0)) Rcpp::stop("The supplied cholesky decomposition has values <= 0.0 on the main diagonal.");
}
vec D = cholDec.diag();
vec out(X.n_rows);
#ifdef _OPENMP
#pragma omp parallel num_threads(ncores) if(ncores > 1)
{
#endif
// Declaring some private variables
uint32_t d = X.n_cols;
uint32_t n = X.n_rows;
vec tmp(d);
double acc;
uint32_t icol, irow, ii;
// For each of the "n" random vectors, forwardsolve the corresponding linear system.
// Forwardsolve because I'm using the lower triangle Cholesky.
#ifdef _OPENMP
#pragma omp for schedule(static)
#endif
for(icol = 0; icol < n; icol++)
{
for(irow = 0; irow < d; irow++)
{
acc = 0.0;
for(ii = 0; ii < irow; ii++) acc += tmp.at(ii) * cholDec.at(irow, ii);
tmp.at(irow) = ( X.at(icol, irow) - mu.at(irow) - acc ) / D.at(irow);
}
out.at(icol) = sum(square(tmp));
}
#ifdef _OPENMP
}
#endif
return out;
}
/*
#Equivalent R function:
.fastMahalanobis <- function(X, mean, mcov)
{
dec <- chol(mcov)
tmp <- forwardsolve(t(dec), t(X) - mean )
colSums( tmp ^ 2 )
}
*/
|