CMS 3D CMS Logo

TB06TreeH2.cc
Go to the documentation of this file.
3 #include "TFile.h"
4 #include "TTree.h"
5 #include "TClonesArray.h"
6 
7 #include <iostream>
8 
10  : m_file(nullptr), m_tree(nullptr), m_data(nullptr), m_dataSize(0) {
11  TDirectory *dir = gDirectory;
12  m_file = new TFile(fileName.c_str(), "RECREATE");
13  m_file->cd();
14  m_tree = new TTree(treeName.c_str(), "Analysis tree");
15  m_tree->SetAutoSave(10000000);
16  dir->cd();
17 
18  // m_tree->cd () ;
19  m_data = new TClonesArray(TB06RecoH2::Class(), 1);
20  m_data->ExpandCreateFast(1);
21 
22  // m_tree->Branch ("EGCO", &m_data, 64000, 2) ;
23  m_tree->Branch("TB06O", &m_data, 64000, 2);
24  m_tree->Print();
25 }
26 
27 // -------------------------------------------------------------------
28 
30  std::cout << "[TB06TreeH2][dtor] saving TTree " << m_tree->GetName() << " with " << m_tree->GetEntries() << " entries"
31  << " on file: " << m_file->GetName() << std::endl;
32 
33  m_file->Write();
34  delete m_tree;
35  m_file->Close();
36  delete m_file;
37  delete m_data;
38 }
39 
40 // -------------------------------------------------------------------
41 
43 void TB06TreeH2::store(const int &tableIsMoving,
44  const int &run,
45  const int &event,
46  const int &S6adc,
47  const double &xhodo,
48  const double &yhodo,
49  const double &xslope,
50  const double &yslope,
51  const double &xquality,
52  const double &yquality,
53  const int &icMax,
54  const int &ietaMax,
55  const int &iphiMax,
56  const double &beamEnergy,
57  const double ampl[49],
58  const int &wcAXo,
59  const int &wcAYo,
60  const int &wcBXo,
61  const int &wcBYo,
62  const int &wcCXo,
63  const int &wcCYo,
64  const double &xwA,
65  const double &ywA,
66  const double &xwB,
67  const double &ywB,
68  const double &xwC,
69  const double &ywC,
70  const float &S1adc,
71  const float &S2adc,
72  const float &S3adc,
73  const float &S4adc,
74  const float &VM1,
75  const float &VM2,
76  const float &VM3,
77  const float &VM4,
78  const float &VM5,
79  const float &VM6,
80  const float &VM7,
81  const float &VM8,
82  const float &VMF,
83  const float &VMB,
84  const float &CK1,
85  const float &CK2,
86  const float &CK3,
87  const float &BH1,
88  const float &BH2,
89  const float &BH3,
90  const float &BH4,
91  const float &TOF1S,
92  const float &TOF2S,
93  const float &TOF1J,
94  const float &TOF2J) {
95  m_data->Clear();
96  TB06RecoH2 *entry = static_cast<TB06RecoH2 *>(m_data->AddrAt(0));
97 
98  entry->reset();
99  // reset (entry->myCalibrationMap) ;
100 
101  entry->tableIsMoving = tableIsMoving;
102  entry->run = run;
103  entry->event = event;
104  entry->S6ADC = S6adc;
105 
106  entry->MEXTLindex = icMax;
107  entry->MEXTLeta = ietaMax;
108  entry->MEXTLphi = iphiMax;
109  entry->MEXTLenergy = ampl[24];
110  entry->beamEnergy = beamEnergy;
111 
112  for (int eta = 0; eta < 7; ++eta)
113  for (int phi = 0; phi < 7; ++phi) {
114  // FIXME capire l'orientamento di phi!
115  // FIXME capire se eta, phi iniziano da 1 o da 0
116  entry->localMap[eta][phi] = ampl[eta * 7 + phi];
117  }
118 
119  //[Edgar] S1 uncleaned, uncalibrated energy
120  entry->S1uncalib_ = ampl[24];
121 
122  //[Edgar] S25 uncleaned, uncalibrated energy
123  for (int eta = 1; eta < 6; ++eta)
124  for (int phi = 1; phi < 6; ++phi) {
125  entry->S25uncalib_ += entry->localMap[eta][phi];
126  }
127 
128  //[Edgar] S49 uncleaned, uncalibrated energy
129  for (int eta = 0; eta < 7; ++eta)
130  for (int phi = 0; phi < 7; ++phi) {
131  entry->S49uncalib_ += entry->localMap[eta][phi];
132  }
133 
134  //[Edgar] S9 uncleaned, uncalibrated energy
135  for (int eta = 2; eta < 5; ++eta)
136  for (int phi = 2; phi < 5; ++phi) {
137  entry->S9uncalib_ += entry->localMap[eta][phi];
138  }
139 
140  entry->xHodo = xhodo;
141  entry->yHodo = yhodo;
142  entry->xSlopeHodo = xslope;
143  entry->ySlopeHodo = yslope;
144  entry->xQualityHodo = xquality;
145  entry->yQualityHodo = yquality;
146  entry->wcAXo_ = wcAXo;
147  entry->wcAYo_ = wcAYo;
148  entry->wcBXo_ = wcBXo;
149  entry->wcBYo_ = wcBYo;
150  entry->wcCXo_ = wcCXo;
151  entry->wcCYo_ = wcCYo;
152  entry->xwA_ = xwA;
153  entry->ywA_ = ywA;
154  entry->xwB_ = xwB;
155  entry->ywB_ = ywB;
156  entry->xwC_ = xwC;
157  entry->ywC_ = ywC;
158  entry->S1adc_ = S1adc;
159  entry->S2adc_ = S2adc;
160  entry->S3adc_ = S3adc;
161  entry->S4adc_ = S4adc;
162  entry->VM1_ = VM1;
163  entry->VM2_ = VM2;
164  entry->VM3_ = VM3;
165  entry->VM4_ = VM4;
166  entry->VM5_ = VM5;
167  entry->VM6_ = VM6;
168  entry->VM7_ = VM7;
169  entry->VM8_ = VM8;
170  entry->VMF_ = VMF;
171  entry->VMB_ = VMB;
172  entry->CK1_ = CK1;
173  entry->CK2_ = CK2;
174  entry->CK3_ = CK3;
175  entry->BH1_ = BH1;
176  entry->BH2_ = BH2;
177  entry->BH3_ = BH3;
178  entry->BH4_ = BH4;
179  entry->TOF1S_ = TOF1S;
180  entry->TOF2S_ = TOF2S;
181  entry->TOF1J_ = TOF1J;
182  entry->TOF2J_ = TOF2J;
183 
184  entry->convFactor = 0.;
185 
186  /*
187  // loop over the 5x5 see (1)
188  for (int xtal=0 ; xtal<25 ; ++xtal)
189  {
190  int ieta = xtal/5 + 3 ;
191  int iphi = xtal%5 + 8 ;
192  entry->myCalibrationMap[ieta][iphi] = ampl[xtal] ;
193  } // loop over the 5x5
194 
195  entry->electron_Tr_Pmag_ = beamEnergy ;
196 
197  entry->centralCrystalEta_ = ietaMax ;
198  entry->centralCrystalPhi_ = iphiMax ;
199  entry->centralCrystalEnergy_ = ampl[12] ;
200 
201  // this is a trick
202  entry->electron_Tr_Peta_ = xhodo ;
203  entry->electron_Tr_Pphi_ = yhodo ;
204  */
205  m_tree->Fill();
206 }
207 
208 // -------------------------------------------------------------------
209 
210 void TB06TreeH2::reset(float crystal[11][21]) {
211  for (int eta = 0; eta < 11; ++eta) {
212  for (int phi = 0; phi < 21; ++phi) {
213  crystal[eta][phi] = -999.;
214  }
215  }
216 }
217 
218 // -------------------------------------------------------------------
219 
221  TB06RecoH2 *entry = static_cast<TB06RecoH2 *>(m_data->AddrAt(0));
222 
223  std::cout << "[TB06TreeH2][check]reading . . . \n";
224  std::cout << "[TB06TreeH2][check] entry->run: " << entry->run << "\n";
225  std::cout << "[TB06TreeH2][check] entry->event: " << entry->event << "\n";
226  std::cout << "[TB06TreeH2][check] entry->tableIsMoving: " << entry->tableIsMoving << "\n";
227  std::cout << "[TB06TreeH2][check] entry->MEXTLeta: " << entry->MEXTLeta << "\n";
228  std::cout << "[TB06TreeH2][check] entry->MEXTLphi: " << entry->MEXTLphi << "\n";
229  std::cout << "[TB06TreeH2][check] entry->MEXTLenergy: " << entry->MEXTLenergy << "\n";
230 
231  for (int eta = 0; eta < 7; ++eta)
232  for (int phi = 0; phi < 7; ++phi)
233  std::cout << "[TB06TreeH2][check] entry->localMap[" << eta << "][" << phi << "]: " << entry->localMap[eta][phi]
234  << "\n";
235 
236  std::cout << "[TB06TreeH2][check] entry->xHodo: " << entry->xHodo << "\n";
237  std::cout << "[TB06TreeH2][check] entry->yHodo: " << entry->yHodo << "\n";
238  std::cout << "[TB06TreeH2][check] entry->xSlopeHodo: " << entry->xSlopeHodo << "\n";
239  std::cout << "[TB06TreeH2][check] entry->ySlopeHodo: " << entry->ySlopeHodo << "\n";
240  std::cout << "[TB06TreeH2][check] entry->xQualityHodo: " << entry->xQualityHodo << "\n";
241  std::cout << "[TB06TreeH2][check] entry->yQualityHodo: " << entry->yQualityHodo << "\n";
242  std::cout << "[TB06TreeH2][check] entry->convFactor: " << entry->convFactor << "\n";
243 
244  /* to be implemented with the right variables
245  std::cout << "[TB06TreeH2][check] ------------------------" << std::endl ;
246  std::cout << "[TB06TreeH2][check] " << entry->variable_name << std::endl ;
247  */
248 }
249 
250 /* (1) to fill the 25 crystals vector
251 
252  for (UInt_t icry=0 ; icry<25 ; ++icry)
253  {
254  UInt_t row = icry / 5 ;
255  Int_t column = icry % 5 ;
256  try
257  {
258  EBDetId tempo (maxHitId.ieta()+column-2,
259  maxHitId.iphi()+row-2,
260  EBDetId::ETAPHIMODE) ;
261 
262  Xtals5x5.push_back (tempo) ;
263  amplitude [icry] = hits->find (Xtals5x5[icry])->energy () ;
264 
265  }
266  catch ( std::runtime_error &e )
267  {
268  std::cout << "Cannot construct 5x5 matrix around EBDetId "
269  << maxHitId << std::endl ;
270  return ;
271  }
272  } // loop over the 5x5 matrix
273 
274 
275 */
TClonesArray * m_data
Definition: TB06TreeH2.h:81
void store(const int &tableIsMoving, const int &run, const int &event, const int &S6adc, const double &xhodo, const double &yhodo, const double &xslope, const double &yslope, const double &xquality, const double &yquality, const int &icMax, const int &ietaMax, const int &iphiMax, const double &beamEnergy, const double ampl[49], const int &wcAXo, const int &wcAYo, const int &wcBXo, const int &wcBYo, const int &wcCXo, const int &wcCYo, const double &xwA, const double &ywA, const double &xwB, const double &ywB, const double &xwC, const double &ywC, const float &S1adc, const float &S2adc, const float &S3adc, const float &S4adc, const float &VM1, const float &VM2, const float &VM3, const float &VM4, const float &VM5, const float &VM6, const float &VM7, const float &VM8, const float &VMF, const float &VMB, const float &CK1, const float &CK2, const float &CK3, const float &BH1, const float &BH2, const float &BH3, const float &BH4, const float &TOF1S, const float &TOF2S, const float &TOF1J, const float &TOF2J)
to be called at each loop
Definition: TB06TreeH2.cc:43
~TB06TreeH2()
dtor
Definition: TB06TreeH2.cc:29
TB06TreeH2(const std::string &fileName="TB06Tree.root", const std::string &treeName="Analysis")
ctor
Definition: TB06TreeH2.cc:9
void reset(float crystal[11][21])
Definition: TB06TreeH2.cc:210
TTree * m_tree
Definition: TB06TreeH2.h:79
void check()
Definition: TB06TreeH2.cc:220
TFile * m_file
Definition: TB06TreeH2.h:78
Definition: event.py:1