StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFtpcCalibMaker.cxx
1 // $Id: StFtpcCalibMaker.cxx,v 1.14 2017/03/02 18:28:38 jeromel Exp $
2 //
3 // $Log: StFtpcCalibMaker.cxx,v $
4 // Revision 1.14 2017/03/02 18:28:38 jeromel
5 // No changes, just empty lines cleanup
6 //
7 // Revision 1.13 2009/12/09 14:41:30 jcs
8 // delta_t0 and delta_gas can now both = 0
9 // new space point calculation always necessary since reconstruction done with data t0
10 //
11 // Revision 1.12 2009/11/18 12:09:50 jcs
12 // add USE_LOCAL_DRIFTMAP instructions
13 //
14 // Revision 1.11 2009/10/14 15:59:55 jcs
15 // changes to be able to vary the gas temperature in addition to varying t0 and
16 // gas composition
17 //
18 // Revision 1.10 2009/08/04 08:42:09 jcs
19 // The 'perfect' gain table and adjustAverageWest = adjustAverageEast = 0.0
20 // are used for laser run calibration
21 //
22 // Revision 1.9 2008/05/15 22:39:47 jcs
23 // re-activate helix fit
24 //
25 // Revision 1.8 2008/05/13 19:14:58 jcs
26 // get Laser t0 from Calibrations_ftpc/ftpcElectronics offline database table
27 // clean up comments
28 //
29 // Revision 1.7 2008/04/11 17:00:55 nav
30 // *** empty log message ***
31 //
32 // Revision 1.6 2007/04/28 17:56:08 perev
33 // Redundant StChain.h removed
34 //
35 // Revision 1.5 2007/01/22 13:08:15 jcs
36 // replace cout with LOG
37 //
38 // Revision 1.4 2007/01/19 08:53:54 jcs
39 // replace gMessMgr with LOG
40 //
41 // Revision 1.3 2006/04/04 14:34:39 jcs
42 // replace assert with a warning message and return kStWarn
43 //
44 // Revision 1.2 2006/04/04 10:57:03 jcs
45 // Fix memory leak
46 //
47 // Revision 1.1 2006/03/13 19:59:56 jcs
48 // commit initial version of the FTPC calibration maker
49 //
50 
52 // //
53 // StFtpcCalibMaker class for Makers //
54 // //
56 
57 #include <Stiostream.h>
58 #include <fstream>
59 #include <cmath>
60 #include "PhysicalConstants.h"
61 
62 #include "StMessMgr.h"
63 #include "StFtpcCalibMaker.h"
64 #include "StFtpcLaserCalib.hh"
65 #include "StFtpcLaserTrafo.hh"
66 
67 #include "St_DataSetIter.h"
68 #include "St_DataSet.h"
69 #include "StMessMgr.h"
70 
71 #include "TH1.h"
72 #include "TH2.h"
73 #include "TH3.h"
74 #include "TObjArray.h"
75 
76 ClassImp(StFtpcCalibMaker)
77 
78 //_____________________________________________________________________________
79 
80 
81 StFtpcCalibMaker::StFtpcCalibMaker(const char *name):
82  StMaker(name)
83 {
84 
85 }
86 
87 //_____________________________________________________________________________
88 
89 StFtpcCalibMaker::~StFtpcCalibMaker(){
90 
91 }
92 
93 //_____________________________________________________________________________
94 
95 void StFtpcCalibMaker::GetRunInfo(TString filename){
96 
97  StFtpcLaser *j = new StFtpcLaser();
98 
99  j->Init(filename);
100  j->GetTreeEntry(0);
101  LOG_DEBUG<<"StFtpcCalibMaker::GetRunInfo j->Run.run = "<<j->Run.run<<" j->Run.date = "<<j->Run.date<<" j->Run.time = "<<j->Run.time<<" j->Run.micropertimebin = "<<j->Run.micropertimebin<<" j->Run.normalizedNowPressure = "<<j->Run.normalizedNowPressure <<" j->Run.standardPressure = "<<j->Run.standardPressure<<" j->Run.baseTemperature = "<<j->Run.baseTemperature<<" j->Run.gasTemperatureWest = "<<j->Run.gasTemperatureWest<<" j->Run.gasTemperatureEast = "<<j->Run.gasTemperatureEast<<endm;
102  run = j->Run.run;
103  date = j->Run.date;
104  time = j->Run.time;
105  micropertime = j->Run.micropertimebin;
106  normalizedNowPressure = j->Run.normalizedNowPressure;
107  standardPressure = j->Run.standardPressure;
108  baseTemperature = j->Run.baseTemperature;
109  gasTemperatureWest = j->Run.gasTemperatureWest;
110  gasTemperatureEast = j->Run.gasTemperatureEast;
111  delete j; // destructor matches second constructor only
112  return;
113 
114 }
115 
116 //_____________________________________________________________________________
117 
124 Int_t StFtpcCalibMaker::DbInit(float mbfield)
125 {
126  St_DataSet *ftpc_db = NULL;
127 
128  if ( mbfield > 0.8 ) {
129  SetFlavor("ffp10kv","ftpcVDrift");
130  SetFlavor("ffp10kv","ftpcdVDriftdP");
131  SetFlavor("ffp10kv","ftpcDeflection");
132  SetFlavor("ffp10kv","ftpcdDeflectiondP");
133  LOG_INFO << "StFtpcCalibMaker::DbInit - flavor set to ffp10kv"<<endm;
134  }
135  else if ( mbfield > 0.2 ) {
136  SetFlavor("hfp10kv","ftpcVDrift");
137  SetFlavor("hfp10kv","ftpcdVDriftdP");
138  SetFlavor("hfp10kv","ftpcDeflection");
139  SetFlavor("hfp10kv","ftpcdDeflectiondP");
140  LOG_INFO << "StFtpcCalibMaker::DbInit - flavor set to hfp10kv"<<endm;
141  }
142  else if ( mbfield > -0.2 ) {
143  SetFlavor("zf10kv","ftpcVDrift");
144  SetFlavor("zf10kv","ftpcdVDriftdP");
145  SetFlavor("zf10kv","ftpcDeflection");
146  SetFlavor("zf10kv","ftpcdDeflectiondP");
147  LOG_INFO << "StFtpcCalibMaker::DbInit - flavor set to zf10kv"<<endm;
148  }
149  else if ( mbfield > -0.8 ) {
150  SetFlavor("hfn10kv","ftpcVDrift");
151  SetFlavor("hfn10kv","ftpcdVDriftdP");
152  SetFlavor("hfn10kv","ftpcDeflection");
153  SetFlavor("hfn10kv","ftpcdDeflectiondP");
154  LOG_INFO << "StFtpcCalibMaker::DbInit - flavor set to hfn10kv"<<endm;
155  }
156  else {
157  SetFlavor("ffn10kv","ftpcVDrift");
158  SetFlavor("ffn10kv","ftpcdVDriftdP");
159  SetFlavor("ffn10kv","ftpcDeflection");
160  SetFlavor("ffn10kv","ftpcdDeflectiondP");
161  LOG_INFO << "StFtpcCalibMaker::DbInit - flavor set to ffn10kv"<<endm;
162  }
163 
164  ftpc_db = GetDataBase("ftpc");
165  if (!ftpc_db) {
166  LOG_WARN << "StFtpcCalibMaker::DbInit - run parameter database StarDb/ftpc not found"<<endm;
167  return kStWarn;
168  }
169  St_DataSetIter local(ftpc_db);
170 
171  m_clusterpars = (St_ftpcClusterPars *)local("ftpcClusterPars");
172 
173  // USE_LOCAL_DRIFTMAP:
174  // To use the FTPC drift map tables in $PWD/StarDb instead of those
175  // in the MySQL offline database, uncomment the following 4 lines of code
176  //m_vdrift = (St_ftpcVDrift *)local("ftpcVDrift");
177  //m_deflection = (St_ftpcDeflection *)local("ftpcDeflection");
178  //m_dvdriftdp = (St_ftpcdVDriftdP *)local("ftpcdVDriftdP");
179  //m_ddeflectiondp = (St_ftpcdDeflectiondP *)local("ftpcdDeflectiondP");
180 
181  St_DataSet *ftpc_geometry_db = GetDataBase("Geometry/ftpc");
182  St_DataSetIter dblocal_geometry(ftpc_geometry_db);
183 
184  m_dimensions = (St_ftpcDimensions *)dblocal_geometry("ftpcDimensions");
185  m_padrow_z = (St_ftpcPadrowZ *)dblocal_geometry("ftpcPadrowZ");
186 
187  St_DataSet *ftpc_calibrations_db = GetDataBase("Calibrations/ftpc");
188  St_DataSetIter dblocal_calibrations(ftpc_calibrations_db);
189 
190  m_gas= (St_ftpcGas *)dblocal_calibrations("ftpcGas");
191  m_efield = (St_ftpcEField *)dblocal_calibrations("ftpcEField");
192 
193  // USE_LOCAL_DRIFTMAP:
194  // To use the FTPC drift map tables in $PWD/StarDb instead of those
195  // in the MySQL offline database, comment out the following 4 lines of code
196  m_vdrift = (St_ftpcVDrift *)dblocal_calibrations("ftpcVDrift");
197  m_deflection = (St_ftpcDeflection *)dblocal_calibrations("ftpcDeflection");
198  m_dvdriftdp = (St_ftpcdVDriftdP *)dblocal_calibrations("ftpcdVDriftdP");
199  m_ddeflectiondp = (St_ftpcdDeflectiondP *)dblocal_calibrations("ftpcdDeflectiondP");
200 
201  m_driftfield = (St_ftpcDriftField *)dblocal_calibrations("ftpcDriftField");
202  m_electronics = (St_ftpcElectronics *)dblocal_calibrations("ftpcElectronics");
203 
204  // create parameter reader
205  paramReader = new StFtpcParamReader(m_clusterpars);
206 
207  //LOG_DEBUG<<"StFtpcCalibMaker::DbInit - parameter reader created"<<endm;
208 
209  // create db reader
210  dbReader = new StFtpcDbReader(m_dimensions,
211  m_padrow_z,
212  m_efield,
213  m_vdrift,
214  m_deflection,
215  m_dvdriftdp,
216  m_ddeflectiondp,
217  m_electronics,
218  m_gas,
219  m_driftfield);
220 
221  //LOG_DEBUG<<"StFtpcCalibMaker::DbInit - Ftpc db reader created"<<endm;
222 
223  return kStOK;
224 
225 }
226 //_____________________________________________________________________________
227 
234 void StFtpcCalibMaker::DoLaserCalib(TString filename,int ftpc, int lsec, int straight, int gfit, int minz, int maxz, int minrad, int maxrad, char* t0, char* gas, float gastemp, float mbfield)
235 {
236 if (ftpc == 1) LOG_INFO<<"StFtpcCalibMaker::DoLaserCalib - entered for FTPC West"<<endm;
237 if (ftpc == 2) LOG_INFO<<"StFtpcCalibMaker::DoLaserCalib - entered for FTPC East"<<endm;
238  Int_t i=0;
239  Int_t ntracksold=0;
240  Int_t neventold=2;
241 
242  if (ftpc == 1){
243  Float_t adjustedAirPressureWest = normalizedNowPressure*((baseTemperature+STP_Temperature)/(gasTemperatureWest+gastemp+STP_Temperature));
244  deltap = adjustedAirPressureWest - standardPressure;
245  LOG_INFO << "d_TempWest = " << gastemp << " adjustedAirPressureWest = " << adjustedAirPressureWest << " deltap West = " << deltap << endm;
246  }
247  if (ftpc == 2) {
248  Float_t adjustedAirPressureEast = normalizedNowPressure*((baseTemperature+STP_Temperature)/(gasTemperatureEast+gastemp+STP_Temperature));
249  deltap = adjustedAirPressureEast - standardPressure;
250  LOG_INFO << "d_TempEast = " << gastemp << " adjustedAirPressureEast = " << adjustedAirPressureEast << " deltap East = " << deltap << endm;
251  }
252 
253  Bool_t laserRun = kTRUE;
254  dbReader->setLaserRun(laserRun);
255 
256  // get Laser t0 from Calibrations_ftpc/ftpcElectronics offline database table
257  tZero = dbReader->laserTZero();
258  LOG_INFO<<"StFtpcCalibMaker::DoLaserCalib() - laserTZero = "<<tZero<<endm;
259 
260  LOG_INFO<<"StFtpcCalibMaker::DoLaserCalib() ..."<<endm;
261 
262  // for (int step=0;step<10;step++)
263 
264  LOG_INFO<<" "<<endm;
265  LOG_INFO<<"StFtpcCalibMaker::DoLaserCalib() - Reading Magnetic-Field maps..."<<endm;
266  LOG_INFO<<" "<<endm;
267 // StMagUtilities *m_magf=new StMagUtilities(kMapped,mbfield,0);
268  StarMagField *m_magf=new StarMagField(StarMagField::kMapped, mbfield, kTRUE);
269  // analoges Problem bei Magnetfeld 0 shift nicht machen in Track !!!
270 
271  // new space point calculation only necessary if atof(t0)!=0 || atof(gas)!=0
272  //if (atof(t0)!=0 || atof(gas)!=0) {
273  // new space point calculation always necessary since reconstruction done with data t0
274  trafo=new StFtpcLaserTrafo(dbReader,paramReader,atof(t0),atof(gas),micropertime,deltap,mbfield,tZero);
275  if (trafo->calcpadtrans()) {
276  LOG_INFO<<"StFtpcCalibMaker::DoLaserCalib - calcpadtrans done !"<<endm;
277  }
278  else{
279  LOG_FATAL<<"StFtpcCalibMaker::DoLaserCalib - fatal error in calcpadtrans !"<<endm;
280  delete trafo;
281  return;
282  }
283  //}
284 
285  StFtpcLaserCalib *l=new StFtpcLaserCalib(ftpc,lsec,straight,gfit,minz,maxz,minrad,maxrad,atof(t0),atof(gas),gastemp,trafo,m_magf);
286 
287  l->Init(filename);
288  l->MakeOutput(filename,t0,gas,gastemp);
289 
290  Int_t maxentries=l->btcluster->GetEntries();
291  LOG_INFO<<" "<<endm;
292  LOG_INFO<<"StFtpcCalibMaker::DoLaserCalib() - processing Cluster-on-Track-Tree with "<<maxentries<<" clusters... please be patient"<<endm;
293  LOG_INFO<<" "<<endm;
294 
295  for (int k=0;k<=maxentries;k++)
296  {
297  if (k%(maxentries/10)==0 && k>0) {
298  //LOG_DEBUG<<"StFtpcCalibMaker::DoLaserCalib() - "<<k<<" cluster on tracks processed"<<endm;
299  }
300  l->GetTreeEntry(k);
301 
302  //calculate hardware sector
303  int hardsec = 6*(int)((l->tcluster.row-1)/2) + l->tcluster.sec;
304 
305  //LOG_DEBUG<<"StFtpcCalibMaker::DoLaserCalib() - hardsec "<< hardsec<<endm;
306 
307  // activate following code line to debug with 2 events
308  //if (l->tevent.nevent>2) break;
309 
310  //LOG_DEBUG<<"StFtpcCalibMaker::DoLaserCalib() - l->tevent.run = "<<l->tevent.run<<" l->tevent.nevent = "<<l->tevent.nevent<<endm;
311  //LOG_DEBUG<<"StFtpcCalibMaker::DoLaserCalib() - l->Run.run = "<<l->Run.run<<" l->Run.micropertimebin = "<<l->Run.micropertimebin<<" l->Run.deltapW = "<<l->Run.deltapW<<" l->Run.deltapE = "<<l->Run.deltapE<<endm;
312 
313  if (l->tevent.nevent==neventold)
314  {
315  if (l->tcluster.ntracks==ntracksold)
316  {
317  l->fillarray(l->thit.x,l->thit.y,l->thit.z,l->thit.ex,l->thit.ey,i,hardsec,l->tcluster.padpos,l->tcluster.padpossigma,l->tcluster.sec,l->tcluster.row,l->tcluster.timepos,l->tcluster.padlength,l->tcluster.timelength,l->tcluster.peakheight,l->tcluster.charge);
318  i++;
319  }
320  else
321  {
322  if (l->laser_straight(l->radius,i)==l->STRAIGHT || l->STRAIGHT==3)
323  if (l->laser_fit(i)==0) {}
324 
325  i=0;
326  l->fillarray(l->thit.x,l->thit.y,l->thit.z,l->thit.ex,l->thit.ey,i,hardsec,l->tcluster.padpos,l->tcluster.padpossigma,l->tcluster.sec,l->tcluster.row,l->tcluster.timepos,l->tcluster.padlength,l->tcluster.timelength,l->tcluster.peakheight,l->tcluster.charge);
327  i++;
328  }
329  }
330  else
331  {
332  if (l->laser_straight(l->radius,i)==l->STRAIGHT || l->STRAIGHT==3)
333  if (l->laser_fit(i)==0) {}
334 
335  i=0;
336  l->fillarray(l->thit.x,l->thit.y,l->thit.z,l->thit.ex,l->thit.ey,i,hardsec,l->tcluster.padpos,l->tcluster.padpossigma,l->tcluster.sec,l->tcluster.row,l->tcluster.timepos,l->tcluster.padlength,l->tcluster.timelength,l->tcluster.peakheight,l->tcluster.charge);
337  i++;
338  }
339 
340  neventold=l->tevent.nevent;
341  ntracksold=l->tcluster.ntracks;
342  }
343 
344  if (l->laser_straight(l->radius,i)==l->STRAIGHT || l->STRAIGHT==3)
345  if (l->laser_fit(i)==0) {}
346 
347  l->MakePs();
348  l->PositionLog();
349 
350  l->analyse_defl();
351 
352  delete l;
353 
354  if (t0!=0)
355  delete trafo;
356 
357  delete m_magf;
358 
359  LOG_INFO<<" "<<endm;
360  LOG_INFO<<"Laser calibration done :-) !"<<endm;
361  LOG_INFO<<" "<<endm;
362 }
363 
364 //_____________________________________________________________________________
365 
372 void StFtpcCalibMaker::DoT0Calib(TString filename, char* t0, char* gas, float mbfield)
373 {
374 
375 
376  Bool_t laserRun = kFALSE;
377  dbReader->setLaserRun(laserRun);
378 
379  tZero = dbReader->tZero();
380  LOG_INFO<<"StFtpcCalibMaker::DoT0Calib entered with filename "<<filename<<" t0 "<<t0<<" gas "<<gas<<" mbfield "<<mbfield<<" and tZero = "<<tZero<<endm;
381 
382  // new space point calculation only necessary if atof(t0)!=0 || atof(gas)!=0
383  // new space point calculation always necessary since reconstruction done with data t0
384  //if (atof(t0)!=0 || atof(gas)!=0) {
385  // FTPC West
386 
387  deltap = deltapW;
388  LOG_INFO<<"StFtpcCalibMaker::DoT0Calib deltap = deltapW = "<<deltap<<endm;
389 
390 
391  trafo = new StFtpcLaserTrafo(dbReader,paramReader,atof(t0),atof(gas),micropertime,deltap,mbfield,tZero);
392 
393  if (trafo->calcpadtrans()) {
394  LOG_INFO<<"StFtpcCalibMaker::DoT0Calib - calcpadtrans (west) done !"<<endm;
395  }
396  else {
397  LOG_FATAL<<"StFtpcCalibMaker::DoT0Calib - fatal error in calcpadtrans west !"<<endm;
398  delete trafo;
399  return;
400  }
401 
402  // FTPC East
403 
404  deltap = deltapE;
405  LOG_INFO<<"StFtpcCalibMaker::DoT0Calib deltap = deltapE = "<<deltap<<endm;
406 
407  trafo2 = new StFtpcLaserTrafo(dbReader,paramReader,atof(t0),atof(gas),micropertime,deltap,mbfield,tZero);
408 
409  if (trafo2->calcpadtrans()) {
410  LOG_INFO<<"StFtpcCalibMaker::DoT0Calib - calcpadtrans (east) done !"<<endm;
411  }
412  else {
413  LOG_FATAL<<"StFtpcCalibMaker::DoT0Calib - fatal error in calcpadtrans east !"<<endm;
414  delete trafo2;
415  return;
416  }
417 
418  //}
419 
420 
421  LOG_INFO<<"StFtpcCalibMaker::DoT0Calib() ..."<<endm;
422 
423  HistInit(4,filename,t0,gas);
424 
425  StFtpcLaser *l=new StFtpcLaser();
426 
427  l->Init(filename);
428 
429  Int_t maxentries=(int) l->bcluster->GetEntries();
430 
431  LOG_INFO<<" "<<endm;
432  LOG_INFO<<"StFtpcCalibMaker::DoT0Calib() - processing Cluster-Tree with "<<maxentries<<" clusters... please be patient"<<endm;
433  LOG_INFO<<" "<<endm;
434 
435  Float_t x,y,rad;//,phi;
436 
437  for (int k=0;k<=maxentries;k++) {
438 
439  if (k%(maxentries/10)==0 && k>0) {
440  //LOG_DEBUG<<k<<" cluster on tracks processed"<<endm;
441  }
442 
443  l->GetClusterTreeEntry(k);
444 
445  // new space point calculation always necessary since reconstruction done with data t0
446  //if (atof(t0)!=0 || atof(gas)!=0) {
447  if (l->cluster.sec<31)
448  trafo->padtrans(l->cluster.row,l->cluster.sec,l->cluster.timepos,l->cluster.padpos,&x,&y);
449  else
450  trafo2->padtrans(l->cluster.row,l->cluster.sec,l->cluster.timepos,l->cluster.padpos,&x,&y);
451  //}
452  //else {
453  // x=l->hit.x;
454  // y=l->hit.y;
455  //}
456 
457  rad=sqrt(x*x+y*y);
458 
459  if (l->cluster.sec<31) {
460  hradwall->Fill(rad);
461  hradw->Fill(rad);
462  htimew->Fill(l->cluster.timepos);
463  }
464  else {
465  hradeall->Fill(rad);
466  hrade->Fill(rad);
467  htimee->Fill(l->cluster.timepos);
468  }
469 
470  }
471 
472  delete l;
473 
474  if (t0!=0)
475  {
476  delete trafo;delete trafo2;
477  }
478 
479  MakeT0Ps(4,filename,t0,gas);
480 
481  // save histogram file
482  anaf->Write();
483  anaf->Close();
484 
485  LOG_INFO<<" "<<endm;
486  LOG_INFO<<"T0 calibration done :-) !"<<endm;
487  LOG_INFO<<" "<<endm;
488 }
489 
490 //_____________________________________________________________________________
491 
498 void StFtpcCalibMaker::HistInit(int nradbins,TString fname, char* t0, char* gas)
499 {
500 
501  TString outname=fname;
502  outname +="_";
503  outname +=t0;
504  outname +="_";
505  outname +=gas;
506  outname +="_t0.root";
507 
508  // DEBUG :
509  //LOG_DEBUG<<" "<<endm;LOG_DEBUG<<"HistInit() ..."<<endm;LOG_DEBUG<<" "<<endm;
510 
511  LOG_INFO<<"StFtpcCalibMaker::HistInit - Store histograms in ROOT-file : "<<outname<<endm;
512 
513  anaf=new TFile(outname,"RECREATE");
514 
515  hradeall=new TH1F("rad_east_all","radius FTPC East",nradbins*31,0.5,31.5);
516  hradwall=new TH1F("rad_west_all","radius FTPC West",nradbins*31,0.5,31.5);
517 
518  hrade=new TH1F("rad_east","radius FTPC East",nradbins*7,5,12);
519  hradw=new TH1F("rad_west","radius FTPC West",nradbins*7,5,12);
520 
521  htimee=new TH1F("time_east","Timepos. FTPC East",45,140,185);
522  htimew=new TH1F("time_west","Timepos. FTPC West",45,140,185);
523 
524 }
525 
526 //_____________________________________________________________________________
527 
534 void StFtpcCalibMaker::MakeT0Ps(int nradbins,TString psname, char* t0, char* gas)
535 {
536 
537  TString outname=psname;
538  outname +="_";
539  outname +=t0;
540  outname +="_";
541  outname +=gas;
542  outname +="_t0.ps";
543 
544  LOG_INFO<<" "<<endl;
545  LOG_INFO<<"StFtpcCalibMaker::MakeT0Ps - make ps file "<<outname<<endm;
546  LOG_INFO<<" "<<endm;
547 
548  TCanvas *c1 = new TCanvas("c1","ps",200,10,700,500);
549 
550  TPostScript *fps=new TPostScript(outname,112);
551 
552  fps->NewPage();
553 
554  c1->Divide(2,2);
555 
556  TLine *inner=new TLine();
557 
558  c1->cd(1);
559  hradwall->Draw();inner=new TLine(7.8,0,7.8,hradwall->GetMaximum());inner->SetLineColor(2);inner->Draw();
560  c1->cd(2);
561  hradeall->Draw();inner=new TLine(7.8,0,7.8,hradeall->GetMaximum());inner->SetLineColor(2);inner->Draw();
562  c1->cd(3);
563  hradw->Draw();inner=new TLine(7.8,0,7.8,hradw->GetMaximum());inner->SetLineColor(2);inner->Draw();
564  c1->cd(4);
565  hrade->Draw();inner=new TLine(7.8,0,7.8,hrade->GetMaximum());inner->SetLineColor(2);inner->Draw();
566 
567  c1->Update();
568 
569  //fps->NewPage();
570 
571  fps->NewPage();
572 
573  c1->cd(1);
574 
575  hradw->Scale(1/hradw->Integral(0,nradbins*7));
576  hrade->Scale(1/hrade->Integral(0,nradbins*7));
577 
578  hrade->Draw();hradw->SetLineColor(3);hradw->Draw("same");
579  inner=new TLine(7.8,0,7.8,hrade->GetMaximum());inner->SetLineColor(2);inner->Draw();
580 
581  c1->cd(2);
582  htimee->Scale(1/htimee->GetEntries());htimew->Scale(1/htimew->GetEntries());
583  htimee->DrawCopy();htimew->SetLineColor(3);htimew->DrawCopy("same");
584 
585  //c1->cd(4);
586  //htimee->DrawCopy();
587 
588  c1->Update();
589 
590  fps->Close();
591 
592  delete c1;
593 }
virtual Int_t DbInit(float mbfield)
The FTPC calibration maker.
void HistInit(int nradbins, TString fname, char *t0, char *gas)
void DoLaserCalib(TString filename, int ftpc, int lsec, int straight, int gfit, int minz, int maxz, int minrad, int maxrad, char *t0, char *gas, float gastemp, float mbfield)
void DoT0Calib(TString filename, char *t0, char *gas, float mbfield)
void MakeT0Ps(int nradbins, TString psname, char *t0, char *gas)
Definition: Stypes.h:42
Definition: Stypes.h:40