CGEM BOSS 6.6.5.f
BESIII Offline Software System
Loading...
Searching...
No Matches
bak_ReadCosmicRayData-00-00-05/src/TestInputOutput.cxx
Go to the documentation of this file.
1// **************************************************************************/
2// authors: L. Lavezzi (univ. of Torino & INFN, Italy)
3// R. Farinelli (univ. of Ferrara & INFN, Italy)
4//
5
6
7//include system lib
8#include <iostream>
9#include <iomanip>
10#include <string>
11#include <cmath>
12
13// Include files
14#include "GaudiKernel/AlgFactory.h"
15#include "GaudiKernel/DataObject.h"
16#include "GaudiKernel/IEventProcessor.h"
17
18#include "GaudiKernel/Incident.h"
19#include "GaudiKernel/IIncidentSvc.h"
20#include "GaudiKernel/Memory.h"
21
22#include <csignal>
23
24//for save digit & cluster
25#include "GaudiKernel/ISvcLocator.h"
26#include "GaudiKernel/IDataProviderSvc.h"
27#include "GaudiKernel/Bootstrap.h"
28#include "GaudiKernel/RegistryEntry.h"
29#include "GaudiKernel/MsgStream.h"
30
31#include "CgemRawEvent/CgemDigi.h"
32#include "CgemRecEvent/RecCgemCluster.h"
33
34#include "Identifier/CgemID.h"
35
36#include "RawEvent/RawDataUtil.h"
37#include "RawEvent/DigiEvent.h"
38#include "ReconEvent/ReconEvent.h"
39#include "EventModel/EventHeader.h"
40#include "GaudiKernel/SmartDataPtr.h"
41
42#include "ReadCosmicRayData/TestInputOutput.h"
43
44//using namespace std;
45
46TestInputOutput::TestInputOutput(const std::string& name, ISvcLocator* pSvcLocator):
47 Algorithm(name,pSvcLocator){
48
49 declareProperty("Dir_file", Dir_file = "Cosmic_data_01.root");
50 declareProperty("TreeDigi", TreeDigi = "t1");
51 declareProperty("CosmicRayDataSetID", CosmicRayDataSetID = "CR201909");
52
53
54}
55
57
59 MsgStream log(msgSvc(), name());
60 log << MSG::INFO << "TestInputOutput initialize()" << endreq;
61
62 output = new TFile("IOtest_histo.root", "RECREATE");
63
64 // --------------- CgemBoss
65 // x
66 hnhit_L1_x = new TH1F("hnhit_L1_x", "number of hits on layer 1, x srtips", MAXNOFHITS, 0, MAXNOFHITS);
67 hnhit_L2_x = new TH1F("hnhit_L2_x", "number of hits on layer 2, x srtips", MAXNOFHITS, 0, MAXNOFHITS);
68 hcharge_L1_x = new TH1F("hcharge_L1_x", "hcharge on layer 1, x srtips", 50, 0, 50);
69 hcharge_L2_x = new TH1F("hcharge_L2_x", "hcharge on layer 2, x strips", 50, 0, 50);
70 htime_L1_x = new TH1F("htime_L1_x", "htime_L1, x strips", 500, -300, 300);
71 htime_L2_x = new TH1F("htime_L2_x", "htime_L2, x strips", 500, -300, 300);
72 hstripid_L1_x = new TH1F("hstripid_L1_x", "strip IDs of hits on layer 1, x srtips", CgemID::getXSTRIP_MAX(0), 0, CgemID::getXSTRIP_MAX(0));
73 hstripid_L2_x = new TH1F("hstripid_L2_x", "strip IDs of hits on layer 2, x srtips", CgemID::getXSTRIP_MAX(1), 0, CgemID::getXSTRIP_MAX(1));
74
75 // v
76 hnhit_L1_v = new TH1F("hnhit_L1_v", "number of hits on layer 1, v srtips", MAXNOFHITS, 0, MAXNOFHITS);
77 hnhit_L2_v = new TH1F("hnhit_L2_v", "number of hits on layer 2, v srtips", MAXNOFHITS, 0, MAXNOFHITS);
78 hcharge_L1_v = new TH1F("hcharge_L1_v", "hcharge on layer 1, v srtips", 50, 0, 50);
79 hcharge_L2_v = new TH1F("hcharge_L2_v", "hcharge on layer 2, v strips", 50, 0, 50);
80 htime_L1_v = new TH1F("htime_L1_v", "htime_L1, v strips", 500, -300, 300);
81 htime_L2_v = new TH1F("htime_L2_v", "htime_L2, v strips", 500, -300, 300);
82 hstripid_L1_v = new TH1F("hstripid_L1_v", "strip IDs of hits on layer 1, v srtips", CgemID::getVSTRIP_MAX(0), 0, CgemID::getVSTRIP_MAX(0));
83 hstripid_L2_v = new TH1F("hstripid_L2_v", "strip IDs of hits on layer 2, v srtips", CgemID::getVSTRIP_MAX(1), 0, CgemID::getVSTRIP_MAX(1));
84
85 // ----------------- GRAAL
86 // x
87 hnhit_L1_x_GRAAL = new TH1F("hnhit_L1_x_GRAAL", "GRAAL: number of hits on layer 1, x srtips", MAXNOFHITS, 0, MAXNOFHITS);
88 hnhit_L2_x_GRAAL = new TH1F("hnhit_L2_x_GRAAL", "GRAAL: number of hits on layer 2, x srtips", MAXNOFHITS, 0, MAXNOFHITS);
89 hcharge_L1_x_GRAAL = new TH1F("hcharge_L1_x_GRAAL", "GRAAL: hcharge on layer 1, x srtips", 50, 0, 50);
90 hcharge_L2_x_GRAAL = new TH1F("hcharge_L2_x_GRAAL", "GRAAL: hcharge on layer 2, x strips", 50, 0, 50);
91 htime_L1_x_GRAAL = new TH1F("htime_L1_x_GRAAL", "GRAAL: htime_L1, x strips", 500, -300, 300);
92 htime_L2_x_GRAAL = new TH1F("htime_L2_x_GRAAL", "GRAAL: htime_L2, x strips", 500, -300, 300);
93 hstripid_L1_x_GRAAL = new TH1F("hstripid_L1_x_GRAAL", "GRAAL: strip IDs of hits on layer 1, x srtips",
95 hstripid_L2_x_GRAAL = new TH1F("hstripid_L2_x_GRAAL", "GRAAL: strip IDs of hits on layer 2, x srtips",
97
98 // v
99 hnhit_L1_v_GRAAL = new TH1F("hnhit_L1_v_GRAAL", "GRAAL: number of hits on layer 1, v srtips", MAXNOFHITS, 0, MAXNOFHITS);
100 hnhit_L2_v_GRAAL = new TH1F("hnhit_L2_v_GRAAL", "GRAAL: number of hits on layer 2, v srtips", MAXNOFHITS, 0, MAXNOFHITS);
101 hcharge_L1_v_GRAAL = new TH1F("hcharge_L1_v_GRAAL", "GRAAL: hcharge on layer 1, v srtips", 50, 0, 50);
102 hcharge_L2_v_GRAAL = new TH1F("hcharge_L2_v_GRAAL", "GRAAL: hcharge on layer 2, v strips", 50, 0, 50);
103 htime_L1_v_GRAAL = new TH1F("htime_L1_v_GRAAL", "GRAAL: htime_L1, v strips", 500, -300, 300);
104 htime_L2_v_GRAAL = new TH1F("htime_L2_v_GRAAL", "GRAAL: htime_L2, v strips", 500, -300, 300);
105 hstripid_L1_v_GRAAL = new TH1F("hstripid_L1_v_GRAAL", "GRAAL: strip IDs of hits on layer 1, v srtips",
107 hstripid_L2_v_GRAAL = new TH1F("hstripid_L2_v_GRAAL", "GRAAL: strip IDs of hits on layer 2, v srtips",
109
110
111 TString TDir_file(Dir_file);
112 f = new TFile(TDir_file);
113 TString TTreeDigi(TreeDigi);
114 Tdigi = (TTree*)f->Get(TTreeDigi);
115
116 // Get Cgem digi tree
117 Tdigi->SetBranchAddress("nGemHit", &m_nGemHit); // nof GEM hits
118 // information on the IDs
119 Tdigi->SetBranchAddress("GemHit_plane", m_plane); // plane
120 Tdigi->SetBranchAddress("GemHit_view", m_view); // view (axial or stereo strips)
121 Tdigi->SetBranchAddress("GemHit_strip", m_strip); // strip no.
122 // physical information
123 Tdigi->SetBranchAddress("GemHit_q", m_charge); // charge (fC)
124 Tdigi->SetBranchAddress("GemHit_time", m_time); // time (ns)
125
126 Ind_Entry_D = 0;
127
128 return StatusCode::SUCCESS;
129}
130
131StatusCode TestInputOutput::execute(){
132 MsgStream log(msgSvc(), name());
133 cout << "->TestInputOutput::execute" << endl;
134
135
136 bool gotit = GetReferenceFromGRAAL();
137 if(gotit == false)
138 cout<<"Could not accesss Reference Histo!"<<endl;
139
140
141 //interface to event data service
142 ISvcLocator* svcLocator = Gaudi::svcLocator();
143 StatusCode sc=svcLocator->service("EventDataSvc", m_evtSvc);
144 if (sc.isFailure())
145 cout<<"Could not accesss EventDataSvc!"<<endl;
146
147 //retrieve CGEM digits from TDS
148 SmartDataPtr<CgemDigiCol> aDigiCol(m_evtSvc,"/Event/Digi/CgemDigiCol");
149 if(!aDigiCol)
150 cout<<"Could not retrieve CGEM digi collection"<<endl;
151
152 int nhit = aDigiCol->size();
153 int nhit_L1_x = 0;
154 int nhit_L2_x = 0;
155 int nhit_L1_v = 0;
156 int nhit_L2_v = 0;
157 CgemDigiCol::iterator iDigiCol;
158 for(iDigiCol=aDigiCol->begin(); iDigiCol!=aDigiCol->end(); iDigiCol++)
159 {
160 const Identifier ident = (*iDigiCol)->identify();
161
162 /**
163 cout<<"testing layer: "<<CgemID::layer(ident)
164 <<" sheet: "<<CgemID::sheet(ident)
165 <<" strip: "<<CgemID::strip(ident)
166 <<" charge: "<<(*iDigiCol)->getCharge_fc()
167 <<" time: "<<(*iDigiCol)->getTime_ns()<<endl;
168 **/
169 int layerid = CgemID::layer(ident);
170 int stripid = CgemID::strip(ident);
171 bool is_xstrip = CgemID::is_xstrip(ident);
172 double charge = (*iDigiCol)->getCharge_fc();
173 double time = (*iDigiCol)->getTime_ns();
174
175 if(is_xstrip == 1) {
176
177 if(layerid == 0) { // LAYER 1, x strips
178 hcharge_L1_x->Fill(charge);
179 htime_L1_x->Fill(time);
180 hstripid_L1_x->Fill(stripid);
181 nhit_L1_x++;
182 }
183 else if(layerid == 1) { // LAYER 2, x strips
184 hcharge_L2_x->Fill(charge);
185 htime_L2_x->Fill(time);
186 hstripid_L2_x->Fill(stripid);
187 nhit_L1_v++;
188 }
189 }
190 else {
191 if(layerid ==0) { // LAYER 1, v strips
192 hcharge_L1_v->Fill(charge);
193 htime_L1_v->Fill(time);
194 hstripid_L1_v->Fill(stripid);
195 nhit_L2_x++;
196 }
197 else if(layerid == 1) { // LAYER 2, v strips
198 hcharge_L2_v->Fill(charge);
199 htime_L2_v->Fill(time);
200 hstripid_L2_v->Fill(stripid);
201 nhit_L2_v++;
202 }
203 }
204 }
205
206 hnhit_L1_x->Fill(nhit_L1_x);
207 hnhit_L1_v->Fill(nhit_L1_v);
208 hnhit_L2_x->Fill(nhit_L2_x);
209 hnhit_L2_v->Fill(nhit_L2_v);
210
211 cout<<"end of retrieve CGEM digi collection " << nhit <<endl;
212
213 // if(Ind_Entry_D==No_Entries_D-1)
214 // {
215 // log << MSG::INFO << "scheduling a event processing stop...." << endreq;
216 // SmartIF<IEventProcessor> ep(serviceLocator());
217 // if (ep) ep->stopRun();
218 // }
219
220 return StatusCode::SUCCESS;
221}
222
223StatusCode TestInputOutput::finalize(){
224 MsgStream log(msgSvc(),name());
225 log << MSG::INFO << "TestInputOutput finalize()" << endreq;
226 output->Write();
227
228 CheckIO();
229 return StatusCode::SUCCESS;
230}
231
232
233bool TestInputOutput::GetReferenceFromGRAAL(){
234
235 if(CosmicRayDataSetID =="CR201909") {
236
237 Tdigi->GetEntry(Ind_Entry_D);
238 Ind_Entry_D++;
239
240
241 int nhit_L1_x = 0;
242 int nhit_L1_v =0;
243 int nhit_L2_x =0;
244 int nhit_L2_v =0;
245
246 for(int ihit = 0; ihit < m_nGemHit; ihit++) {
247
248 if(m_view[ihit] == 2) {
249 if(m_plane[ihit] == 0) { // LAYER 1, x strips
250 hcharge_L1_x_GRAAL->Fill(m_charge[ihit]);
251 htime_L1_x_GRAAL->Fill(m_time[ihit]);
252 hstripid_L1_x_GRAAL->Fill(m_strip[ihit]);
253 nhit_L1_x++;
254 }
255 else if(m_plane[ihit] == 1 || m_plane[ihit] == 2) { // LAYER 2, x strips
256 hcharge_L2_x_GRAAL->Fill(m_charge[ihit]);
257 htime_L2_x_GRAAL->Fill(m_time[ihit]);
258 if(m_strip[ihit] < CgemID::getXSTRIP_MAX(1)) hstripid_L2_x_GRAAL->Fill(m_strip[ihit]);
259 else hstripid_L2_x_GRAAL->Fill(m_strip[ihit] - CgemID::getXSTRIP_MAX(1));
260 nhit_L1_v++;
261 }
262 }
263 else if(m_view[ihit] == 3) {
264 if(m_plane[ihit] == 0) { // LAYER 1, v strips
265 hcharge_L1_v_GRAAL->Fill(m_charge[ihit]);
266 htime_L1_v_GRAAL->Fill(m_time[ihit]);
267 hstripid_L1_v_GRAAL->Fill(m_strip[ihit]);
268 nhit_L2_x++;
269 }
270 else if(m_plane[ihit] == 1 || m_plane[ihit] == 2) { // LAYER 2, v strips
271 hcharge_L2_v_GRAAL->Fill(m_charge[ihit]);
272 htime_L2_v_GRAAL->Fill(m_time[ihit]);
273 if(m_strip[ihit] < CgemID::getVSTRIP_MAX(1)) hstripid_L2_v_GRAAL->Fill(m_strip[ihit]);
274 else hstripid_L2_v_GRAAL->Fill(m_strip[ihit] - CgemID::getVSTRIP_MAX(1));
275 nhit_L2_v++;
276 }
277 }
278 }
279
280 hnhit_L1_x_GRAAL->Fill(nhit_L1_x);
281 hnhit_L2_x_GRAAL->Fill(nhit_L2_x);
282 hnhit_L1_v_GRAAL->Fill(nhit_L1_v);
283 hnhit_L2_v_GRAAL->Fill(nhit_L2_v);
284
285 }
286 else return false;
287
288
289 return true;
290
291}
292
293
294void TestInputOutput::CheckIO() {
295
296 cout << endl;
297 cout << "tests: 1 = OK; 0 = WRONG" << endl;
298 cout << "1 - check entries" << endl;
299 // nhit
300 cout << "nhit L1 x " << (hnhit_L1_x->GetEntries() == hnhit_L1_x_GRAAL->GetEntries() ) << endl;
301 cout << "nhit L2 x " << (hnhit_L2_x->GetEntries() == hnhit_L2_x_GRAAL->GetEntries() ) << endl;
302 cout << "nhit L1 v " << (hnhit_L1_v->GetEntries() == hnhit_L1_v_GRAAL->GetEntries() ) << endl;
303 cout << "nhit L2 v " << (hnhit_L2_v->GetEntries() == hnhit_L2_v_GRAAL->GetEntries() ) << endl;
304
305 //stripid
306 cout << "stripid L1 x " << (hstripid_L1_x->GetEntries() == hstripid_L1_x_GRAAL->GetEntries() ) << endl;
307 cout << "stripid L2 x " << (hstripid_L2_x->GetEntries() == hstripid_L2_x_GRAAL->GetEntries() ) << endl;
308 cout << "stripid L1 v " << (hstripid_L1_v->GetEntries() == hstripid_L1_v_GRAAL->GetEntries() ) << endl;
309 cout << "stripid L2 v " << (hstripid_L2_v->GetEntries() == hstripid_L2_v_GRAAL->GetEntries() ) << endl;
310
311 // time
312 cout << "time L1 x " << (htime_L1_x->GetEntries() == htime_L1_x_GRAAL->GetEntries() ) << endl;
313 cout << "time L2 x " << (htime_L2_x->GetEntries() == htime_L2_x_GRAAL->GetEntries() ) << endl;
314 cout << "time L1 v " << (htime_L1_v->GetEntries() == htime_L1_v_GRAAL->GetEntries() ) << endl;
315 cout << "time L2 v " << (htime_L2_v->GetEntries() == htime_L2_v_GRAAL->GetEntries() ) << endl;
316
317 // charge
318 cout << "charge L1 x " << (hcharge_L1_x->GetEntries() == hcharge_L1_x_GRAAL->GetEntries() ) << endl;
319 cout << "charge L2 x " << (hcharge_L2_x->GetEntries() == hcharge_L2_x_GRAAL->GetEntries() ) << endl;
320 cout << "charge L1 v " << (hcharge_L1_v->GetEntries() == hcharge_L1_v_GRAAL->GetEntries() ) << endl;
321 cout << "charge L2 v " << (hcharge_L2_v->GetEntries() == hcharge_L2_v_GRAAL->GetEntries() ) << endl;
322
323 cout << "2 - check mean values" << endl;
324 // nhit
325 cout << "nhit L1 x " << (hnhit_L1_x->GetMean() == hnhit_L1_x_GRAAL->GetMean() ) << endl;
326 cout << "nhit L2 x " << (hnhit_L2_x->GetMean() == hnhit_L2_x_GRAAL->GetMean() ) << endl;
327 cout << "nhit L1 v " << (hnhit_L1_v->GetMean() == hnhit_L1_v_GRAAL->GetMean() ) << endl;
328 cout << "nhit L2 v " << (hnhit_L2_v->GetMean() == hnhit_L2_v_GRAAL->GetMean() ) << endl;
329
330 //stripid
331 cout << "stripid L1 x " << (hstripid_L1_x->GetMean() == hstripid_L1_x_GRAAL->GetMean() ) << endl;
332 cout << "stripid L2 x " << (hstripid_L2_x->GetMean() == hstripid_L2_x_GRAAL->GetMean() ) << endl;
333 cout << "stripid L1 v " << (hstripid_L1_v->GetMean() == hstripid_L1_v_GRAAL->GetMean() ) << endl;
334 cout << "stripid L2 v " << (hstripid_L2_v->GetMean() == hstripid_L2_v_GRAAL->GetMean() ) << endl;
335
336 // time
337 cout << "time L1 x " << (htime_L1_x->GetMean() == htime_L1_x_GRAAL->GetMean() ) << endl;
338 cout << "time L2 x " << (htime_L2_x->GetMean() == htime_L2_x_GRAAL->GetMean() ) << endl;
339 cout << "time L1 v " << (htime_L1_v->GetMean() == htime_L1_v_GRAAL->GetMean() ) << endl;
340 cout << "time L2 v " << (htime_L2_v->GetMean() == htime_L2_v_GRAAL->GetMean() ) << endl;
341
342 // charge
343 cout << "charge L1 x " << (hcharge_L1_x->GetMean() == hcharge_L1_x_GRAAL->GetMean() ) << endl;
344 cout << "charge L2 x " << (hcharge_L2_x->GetMean() == hcharge_L2_x_GRAAL->GetMean() ) << endl;
345 cout << "charge L1 v " << (hcharge_L1_v->GetMean() == hcharge_L1_v_GRAAL->GetMean() ) << endl;
346 cout << "charge L2 v " << (hcharge_L2_v->GetMean() == hcharge_L2_v_GRAAL->GetMean() ) << endl;
347
348 cout << "3 - check RMS values" << endl;
349 // nhit
350 cout << "nhit L1 x " << (hnhit_L1_x->GetRMS() == hnhit_L1_x_GRAAL->GetRMS() ) << endl;
351 cout << "nhit L2 x " << (hnhit_L2_x->GetRMS() == hnhit_L2_x_GRAAL->GetRMS() ) << endl;
352 cout << "nhit L1 v " << (hnhit_L1_v->GetRMS() == hnhit_L1_v_GRAAL->GetRMS() ) << endl;
353 cout << "nhit L2 v " << (hnhit_L2_v->GetRMS() == hnhit_L2_v_GRAAL->GetRMS() ) << endl;
354
355 //stripid
356 cout << "stripid L1 x " << (hstripid_L1_x->GetRMS() == hstripid_L1_x_GRAAL->GetRMS() ) << endl;
357 cout << "stripid L2 x " << (hstripid_L2_x->GetRMS() == hstripid_L2_x_GRAAL->GetRMS() ) << endl;
358 cout << "stripid L1 v " << (hstripid_L1_v->GetRMS() == hstripid_L1_v_GRAAL->GetRMS() ) << endl;
359 cout << "stripid L2 v " << (hstripid_L2_v->GetRMS() == hstripid_L2_v_GRAAL->GetRMS() ) << endl;
360
361 // time
362 cout << "time L1 x " << (htime_L1_x->GetRMS() == htime_L1_x_GRAAL->GetRMS() ) << endl;
363 cout << "time L2 x " << (htime_L2_x->GetRMS() == htime_L2_x_GRAAL->GetRMS() ) << endl;
364 cout << "time L1 v " << (htime_L1_v->GetRMS() == htime_L1_v_GRAAL->GetRMS() ) << endl;
365 cout << "time L2 v " << (htime_L2_v->GetRMS() == htime_L2_v_GRAAL->GetRMS() ) << endl;
366
367 // charge
368 cout << "charge L1 x " << (hcharge_L1_x->GetRMS() == hcharge_L1_x_GRAAL->GetRMS() ) << endl;
369 cout << "charge L2 x " << (hcharge_L2_x->GetRMS() == hcharge_L2_x_GRAAL->GetRMS() ) << endl;
370 cout << "charge L1 v " << (hcharge_L1_v->GetRMS() == hcharge_L1_v_GRAAL->GetRMS() ) << endl;
371 cout << "charge L2 v " << (hcharge_L2_v->GetRMS() == hcharge_L2_v_GRAAL->GetRMS() ) << endl;
372
373
374 cout << "4 - check differences" << endl;
375 // nhit
376 hnhit_L1_x->Add(hnhit_L1_x_GRAAL, -1);
377 hnhit_L2_x->Add(hnhit_L2_x_GRAAL, -1);
378 hnhit_L1_v->Add(hnhit_L1_v_GRAAL, -1);
379 hnhit_L2_v->Add(hnhit_L2_v_GRAAL, -1);
380 cout << "diff nhit L1 x: min= " << hnhit_L1_x->GetMinimum() << " max= " << hnhit_L1_x->GetMaximum() << endl;
381 cout << "diff nhit L2 x: min= " << hnhit_L2_x->GetMinimum() << " max= " << hnhit_L2_x->GetMaximum() << endl;
382 cout << "diff nhit L1 v: min= " << hnhit_L1_v->GetMinimum() << " max= " << hnhit_L1_v->GetMaximum() << endl;
383 cout << "diff nhit L2 v: min= " << hnhit_L2_v->GetMinimum() << " max= " << hnhit_L2_v->GetMaximum() << endl;
384
385 // stripid
386 hstripid_L1_x->Add(hstripid_L1_x_GRAAL, -1);
387 hstripid_L2_x->Add(hstripid_L2_x_GRAAL, -1);
388 hstripid_L1_v->Add(hstripid_L1_v_GRAAL, -1);
389 hstripid_L2_v->Add(hstripid_L2_v_GRAAL, -1);
390 cout << "diff stripid L1 x: min= " << hstripid_L1_x->GetMinimum() << " max= " << hstripid_L1_x->GetMaximum() << endl;
391 cout << "diff stripid L2 x: min= " << hstripid_L2_x->GetMinimum() << " max= " << hstripid_L2_x->GetMaximum() << endl;
392 cout << "diff stripid L1 v: min= " << hstripid_L1_v->GetMinimum() << " max= " << hstripid_L1_v->GetMaximum() << endl;
393 cout << "diff stripid L2 v: min= " << hstripid_L2_v->GetMinimum() << " max= " << hstripid_L2_v->GetMaximum() << endl;
394
395 // time
396 htime_L1_x->Add(htime_L1_x_GRAAL, -1);
397 htime_L2_x->Add(htime_L2_x_GRAAL, -1);
398 htime_L1_v->Add(htime_L1_v_GRAAL, -1);
399 htime_L2_v->Add(htime_L2_v_GRAAL, -1);
400 cout << "diff time L1 x: min= " << htime_L1_x->GetMinimum() << " max= " << htime_L1_x->GetMaximum() << endl;
401 cout << "diff time L2 x: min= " << htime_L2_x->GetMinimum() << " max= " << htime_L2_x->GetMaximum() << endl;
402 cout << "diff time L1 v: min= " << htime_L1_v->GetMinimum() << " max= " << htime_L1_v->GetMaximum() << endl;
403 cout << "diff time L2 v: min= " << htime_L2_v->GetMinimum() << " max= " << htime_L2_v->GetMaximum() << endl;
404
405
406 // charge
407 hcharge_L1_x->Add(hcharge_L1_x_GRAAL, -1);
408 hcharge_L2_x->Add(hcharge_L2_x_GRAAL, -1);
409 hcharge_L1_v->Add(hcharge_L1_v_GRAAL, -1);
410 hcharge_L2_v->Add(hcharge_L2_v_GRAAL, -1);
411 cout << "diff charge L1 x: min= " << hcharge_L1_x->GetMinimum() << " max= " << hcharge_L1_x->GetMaximum() << endl;
412 cout << "diff charge L2 x: min= " << hcharge_L2_x->GetMinimum() << " max= " << hcharge_L2_x->GetMaximum() << endl;
413 cout << "diff charge L1 v: min= " << hcharge_L1_v->GetMinimum() << " max= " << hcharge_L1_v->GetMaximum() << endl;
414 cout << "diff charge L2 v: min= " << hcharge_L2_v->GetMinimum() << " max= " << hcharge_L2_v->GetMaximum() << endl;
415
416
417
418
419
420}
Double_t time
static int strip(const Identifier &id)
static value_type getXSTRIP_MAX(unsigned int f_layer)
static int layer(const Identifier &id)
static value_type getVSTRIP_MAX(unsigned int f_layer)
static bool is_xstrip(const Identifier &id)
TestInputOutput(const std::string &name, ISvcLocator *pSvcLocator)