1 #ifndef TMVA_RSTANDARDSCALER
2 #define TMVA_RSTANDARDSCALER
13 namespace Experimental {
16 class RStandardScaler {
18 std::vector<T> fMeans;
22 RStandardScaler() =
default;
23 RStandardScaler(std::string_view title, std::string_view filename);
24 void Fit(
const RTensor<T>& x);
25 std::vector<T> Compute(
const std::vector<T>& x);
26 RTensor<T> Compute(
const RTensor<T>& x);
27 std::vector<T> GetMeans()
const {
return fMeans; }
28 std::vector<T> GetStds()
const {
return fStds; }
29 void Save(std::string_view title, std::string_view filename);
33 inline RStandardScaler<T>::RStandardScaler(std::string_view title, std::string_view filename) {
34 auto file = TFile::Open(filename.data(),
"READ");
35 RStandardScaler<T>* obj;
36 file->GetObject(title.data(), obj);
37 fMeans = obj->GetMeans();
38 fStds = obj->GetStds();
44 inline void RStandardScaler<T>::Save(std::string_view title, std::string_view filename) {
45 auto file = TFile::Open(filename.data(),
"UPDATE");
46 file->WriteObject<RStandardScaler<T>>(
this, title.data());
52 inline void RStandardScaler<T>::Fit(
const RTensor<T>& x) {
53 const auto shape = x.GetShape();
54 if (shape.size() != 2)
55 throw std::runtime_error(
"Can only fit to input tensor of rank 2.");
57 fMeans.resize(shape[1]);
59 fStds.resize(shape[1]);
62 for (std::size_t i = 0; i < shape[0]; i++) {
63 for (std::size_t j = 0; j < shape[1]; j++) {
67 for (std::size_t i = 0; i < shape[1]; i++) {
68 fMeans[i] /= shape[0];
72 for (std::size_t i = 0; i < shape[0]; i++) {
73 for (std::size_t j = 0; j < shape[1]; j++) {
74 fStds[j] += (x(i, j) - fMeans[j]) * (x(i, j) - fMeans[j]);
77 for (std::size_t i = 0; i < shape[1]; i++) {
78 fStds[i] = std::sqrt(fStds[i] / (shape[0] - 1));
83 inline std::vector<T> RStandardScaler<T>::Compute(
const std::vector<T>& x) {
84 const auto size = x.size();
85 if (size != fMeans.size())
86 throw std::runtime_error(
"Size of input vector is not equal to number of fitted variables.");
88 std::vector<T> y(size);
89 for (std::size_t i = 0; i < size; i++) {
90 y[i] = (x[i] - fMeans[i]) / fStds[i];
97 inline RTensor<T> RStandardScaler<T>::Compute(
const RTensor<T>& x) {
98 const auto shape = x.GetShape();
99 if (shape.size() != 2)
100 throw std::runtime_error(
"Can only compute output for input tensor of rank 2.");
101 if (shape[1] != fMeans.size())
102 throw std::runtime_error(
"Second dimension of input tensor is not equal to number of fitted variables.");
105 for (std::size_t i = 0; i < shape[0]; i++) {
106 for (std::size_t j = 0; j < shape[1]; j++) {
107 y(i, j) = (x(i, j) - fMeans[j]) / fStds[j];
117 #endif // TMVA_RSTANDARDSCALER