hpstr
The Heavy Photon Search Toolkit for Reconstruction (hpstr) provides an interface to physics data from the HPS experiment saved in the LCIO format and converts it into an ROOT based format.
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Svt2DBlHistos.cxx
Go to the documentation of this file.
1#include "Svt2DBlHistos.h"
2#include <math.h>
3#include "TCanvas.h"
4
5Svt2DBlHistos::Svt2DBlHistos(const std::string& inputName, ModuleMapper* mmapper) {
6 m_name = inputName;
7 mmapper_ = mmapper;
8}
9
12
14 //Define vector of hybrid names using ModuleMapper
15 //Use this list to define multiple copies of histograms, one for each hybrid, from json file
16 std::vector<std::string> hybridNames;
17 mmapper_->getStrings(hybridNames);
18 if(debug_ > 0){
19 for(int i = 0; i< hybridNames.size(); i++)
20 std::cout << "Hybrid: " << hybridNames.at(i) << std::endl;
21 }
22 //Define histos
23 //All histogram keys in the JSON file that contain special tag will have multiple copies of that histogram template made, one for each string
24 std::string makeMultiplesTag = "SvtHybrids";
25 HistoManager::DefineHistos(hybridNames, makeMultiplesTag );
26
27}
28
29void Svt2DBlHistos::FillHistograms(std::vector<RawSvtHit*> *rawSvtHits_,float weight) {
30
31 int nhits = rawSvtHits_->size();
32 std::vector<std::string> hybridStrings={};
33 std::string histokey;
34 if(Event_number%10000 == 0) std::cout << "Event: " << Event_number
35 << " Number of RawSvtHits: " << nhits << std::endl;
36
37 //Following Block counts the total number of hits each hybrid records per event
38 int svtHybMulti[4][15] = {0};
39 for (int i = 0; i < nhits; i++)
40 {
41 RawSvtHit* rawSvtHit = rawSvtHits_->at(i);
42 int mod = rawSvtHit->getModule();
43 int lay = rawSvtHit->getLayer();
44 svtHybMulti[mod][lay]++;
45
46 }
47 for (int i =0; i < 4; i++)
48 {
49 for (int j = 1; j < 15; j++)
50 {
51 if (!(j<9 && i>1))
52 {
53 std::string swTag = mmapper_->getStringFromSw("ly"+std::to_string(j)+"_m"+std::to_string(i));
54 hybridStrings.push_back(swTag);
55 Fill1DHisto(swTag+ "_SvtHybridsHitN_h", svtHybMulti[i][j],weight);
56 }
57 }
58 }
59
60 Fill1DHisto("SvtHitMulti_h", nhits,weight);
61 //End of counting block
62
63 //Populates histograms for each hybrid
64 for (int i = 0; i < nhits; i++)
65 {
66 RawSvtHit* rawSvtHit = rawSvtHits_->at(i);
67 auto mod = std::to_string(rawSvtHit->getModule());
68 auto lay = std::to_string(rawSvtHit->getLayer());
69 std::string swTag= mmapper_->getStringFromSw("ly"+lay+"_m"+mod);
70
71 //Manually select which baselines (0 - 6) are included. THIS MUST MATCH THE JSON FILE!
72 int ss = 0;
73 histokey = swTag + "_SvtHybrids_s"+std::to_string(ss)+"_hh";
74 Fill2DHisto(histokey,
75 (float)rawSvtHit->getStrip(),
76 (float)rawSvtHit->getADCs()[ss],
77 weight);
78
79 ss = 3;
80 histokey = swTag + "_SvtHybrids_s"+std::to_string(ss)+"_hh";
81 Fill2DHisto(histokey,
82 (float)rawSvtHit->getStrip(),
83 (float)rawSvtHit->getADCs()[ss],
84 weight);
85
86 }
87
89}
void Fill2DHisto(const std::string &histoName, float valuex, float valuey, float weight=1.)
description
virtual void DefineHistos()
Definition of histograms from json config.
void Fill1DHisto(const std::string &histoName, float value, float weight=1.)
description
std::string m_name
description
std::string getStringFromSw(const std::string &key)
Get the String From Sw.
void getStrings(std::vector< std::string > &strings)
Get list of string modules.
int getStrip()
int getLayer()
Definition RawSvtHit.cxx:88
int * getADCs()
Definition RawSvtHit.cxx:76
int getModule()
Definition RawSvtHit.cxx:92
void DefineHistos()
description
int debug_
description
void FillHistograms(std::vector< RawSvtHit * > *rawSvtHits_, float weight=1.)
description
Svt2DBlHistos(const std::string &inputName, ModuleMapper *mmapper_)
Constructor.
int Event_number
description
ModuleMapper * mmapper_
description