Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-16 09:19:53

0001 /**
0002  * Using 16-bit short integers to store 32-bit floating-point values.
0003  * Joint work with Massachusetts Institute of Technology (MIT) and Institute for Information Industry (III)
0004  * Questions please contact fishyu@iii.org.tw
0005  * Date: May 22, 2021
0006  */
0007 
0008 #include "compressor.h"
0009 
0010 #include <TFile.h>
0011 #include <TTree.h>
0012 
0013 #include <iostream>
0014 #include <fstream>
0015 #include <string>
0016 #include <utility>
0017 
0018 void doCompression(std::ofstream &myfile, const TString& filename);
0019 Float_t computeSd(Float_t* avg, Int_t n_entries, TTree* t, const Float_t* gen_);
0020 
0021 void output_before_vs_after(
0022  const TString& filename,
0023  std::vector<UShort_t>* order, 
0024  std::vector<Float_t>* dict, 
0025  Int_t n_entries, 
0026  TTree* t, 
0027  const Float_t* gen_
0028 );
0029 
0030 uint64_t timeSinceEpochMillisec()
0031 {
0032   return std::chrono::duration_cast<std::chrono::milliseconds>(std::chrono::system_clock::now().time_since_epoch()).count();
0033 }
0034 
0035 void compress_clu_res_float32(std::string &filename)
0036 {
0037   std::ofstream logFile;
0038  logFile.open ("log.csv");
0039  logFile << "file,points,variable,dltSD,bits,dltdltSD/dltSD,avg,LB,UB,msec" << std::endl;
0040  doCompression(logFile, filename);
0041  logFile.close();
0042 }
0043 
0044 void doCompression(std::ofstream &myfile, const TString& filename)
0045 {
0046  TFile *f = new TFile(filename,"read");
0047  TTree *t = (TTree*)f->Get("ntp_trkres");
0048  Float_t genv1_;
0049  Float_t genv2_;
0050  //get branches
0051  t->SetBranchAddress("v1",&genv1_);
0052  t->SetBranchAddress("v2",&genv2_);
0053  Int_t n_entries = (Int_t)t->GetEntries();
0054 
0055  // start compress v1
0056  Float_t avg;
0057  Float_t dltSD = computeSd(&avg, n_entries, t, &genv1_);
0058  for (Int_t numBits = 3; numBits < 12; ++numBits) {
0059    std::vector<UShort_t> v1_order; // order, to replace v1 and save to the new ROOT file
0060    std::vector<Float_t> v1_dict; // dictionay, to save to the new ROOT file
0061    std::vector<size_t> v1_cnt; // how many data points in the dictionary
0062    uint64_t start = timeSinceEpochMillisec();
0063    Float_t dltdltSD = approx(&v1_order, &v1_dict, &v1_cnt, n_entries, t, &genv1_, (size_t) pow(2, numBits));
0064    myfile 
0065    << filename << "," 
0066    << v1_order.size() << ",dphi," 
0067    << dltSD << "," 
0068    << numBits << "," 
0069    << dltdltSD / dltSD << "," 
0070    << avg << "," 
0071    << avg - 4 * dltSD << "," 
0072    << avg + 4 * dltSD << "," 
0073    << timeSinceEpochMillisec() - start << std::endl;
0074    output_before_vs_after(filename + "_dltPHI_" + numBits + "bits.csv", &v1_order, &v1_dict, n_entries, t, &genv1_);
0075  }
0076  // end compress v1
0077  
0078  // start compress v2
0079  dltSD = computeSd(&avg, n_entries, t, &genv2_);
0080  for (Int_t numBits = 3; numBits < 12; ++numBits) {
0081    std::vector<UShort_t> v2_order; // order, to replace v2 and save to the new ROOT file
0082    std::vector<Float_t> v2_dict; // dictionay, to save to the new ROOT file
0083    std::vector<size_t> v2_cnt; // how many data points in the dictionary
0084    uint64_t start = timeSinceEpochMillisec();
0085    Float_t dltdltSD = approx(&v2_order, &v2_dict, &v2_cnt, n_entries, t, &genv2_, (size_t) pow(2, numBits));
0086    myfile 
0087      << filename << "," 
0088      << v2_order.size() << ",dz," 
0089      << dltSD << "," 
0090      << numBits << "," 
0091      << dltdltSD / dltSD << "," 
0092      << avg << "," 
0093      << avg - 4 * dltSD << "," 
0094      << avg + 4 * dltSD << "," 
0095      << timeSinceEpochMillisec() - start << std::endl;
0096      output_before_vs_after(filename + "_dltZ_" + numBits + "bits.csv", &v2_order, &v2_dict, n_entries, t, &genv2_);
0097   }
0098 }
0099 
0100 Float_t computeSd(Float_t* avg, Int_t n_entries, TTree* t, const Float_t* gen_)
0101 {
0102   Double_t sumSquared = 0;
0103   Double_t sum = 0;
0104   for (Int_t j = 0 ; j < n_entries; j++){
0105     t->GetEntry(j);
0106     
0107     sumSquared += (Double_t) *gen_ * (Double_t) *gen_;
0108     sum += (Double_t) *gen_;
0109   }
0110   
0111   Double_t _avg = sum / (Double_t) n_entries;
0112   *avg = _avg;
0113 
0114   return sqrt((sumSquared / (Double_t) n_entries) - _avg * _avg);
0115 }
0116 
0117 void output_before_vs_after(const TString& filename, std::vector<UShort_t>* order, std::vector<Float_t>* dict, Int_t n_entries, TTree* t, const Float_t* gen_)
0118 {
0119   std::ofstream myfile;
0120   myfile.open(filename);
0121   myfile << "before,after" << std::endl;
0122   for (Int_t j = 0 ; j < n_entries; j++){
0123     t->GetEntry(j);
0124     myfile << *gen_ << "," << dict->at(order->at(j)) << std::endl;
0125   }
0126   myfile.close();
0127 }
0128