| 1 |
void ToFLevel2Ex(){ |
| 2 |
|
| 3 |
/* |
| 4 |
* |
| 5 |
* This script shows how use ToFLevel2 output without using PamLevel2 software |
| 6 |
* Notice that here no selection on image tracks is applied |
| 7 |
* All methods and all variables usage are shown |
| 8 |
* |
| 9 |
*/ |
| 10 |
|
| 11 |
|
| 12 |
/* Some plot just as example */ |
| 13 |
|
| 14 |
// Beta versus Rigidity |
| 15 |
TH2F* BetaRig= new TH2F("BetaRig","Beta versus Rigidity",1000,0.1,10.,520,-1.3,1.3); |
| 16 |
|
| 17 |
// dEdx == adc_c measurement |
| 18 |
TH1F* ADChisto= new TH1F("ADChisto","ADChisto",10,0.,20.); |
| 19 |
|
| 20 |
TH1F* TestHisto = new TH1F("TestHisto","TestHisto",500,0.,5000.); |
| 21 |
|
| 22 |
gStyle->SetOptStat(111111); |
| 23 |
|
| 24 |
TrkLevel2* trk_event = new TrkLevel2(); |
| 25 |
ToFLevel2* tof_event = new ToFLevel2(); |
| 26 |
|
| 27 |
TChain *T = new TChain("Tracker"); |
| 28 |
TChain *C = new TChain("ToF"); |
| 29 |
|
| 30 |
T->Add("./8212.Level2.root"); |
| 31 |
C->Add("./8212.Level2.root"); |
| 32 |
|
| 33 |
T->SetBranchAddress("TrkLevel2", &trk_event); // set tracker event adrss |
| 34 |
C->SetBranchAddress("ToFLevel2", &tof_event); // set ToF event adrss |
| 35 |
T->AddFriend(C); // make the trees friends |
| 36 |
|
| 37 |
Int_t nevent = T->GetEntries(); |
| 38 |
cout << " Events in Tracker Tree " << nevent << endl; |
| 39 |
Int_t ntofevent = C->GetEntries(); |
| 40 |
cout << " Events in ToF Tree " << ntofevent << endl; |
| 41 |
|
| 42 |
int j = 0; |
| 43 |
int k = 0; |
| 44 |
|
| 45 |
// for (Int_t ii=0; ii<nevent;ii++){ |
| 46 |
for (Int_t ii=0; ii<1000;ii++){ |
| 47 |
|
| 48 |
T->GetEntry(ii); |
| 49 |
|
| 50 |
/* |
| 51 |
* infos in ToFLevel2 class |
| 52 |
* |
| 53 |
*/ |
| 54 |
|
| 55 |
// j flag: it's a decoded value (array wioth 6 elements); |
| 56 |
// for each layer, it summed up if there |
| 57 |
// are more than one hit per layer |
| 58 |
// tof11_j = tof11_j + 2**(i-1) |
| 59 |
|
| 60 |
for(Int_t i=0; i<6; i++) |
| 61 |
Int_t tof11_j = tof_event->tof_j_flag[i]; |
| 62 |
|
| 63 |
/* |
| 64 |
* if you need informations from ToF and Tracker, |
| 65 |
* you have to read both the detector classes |
| 66 |
*/ |
| 67 |
|
| 68 |
/* ToFTrkVar is the track related class */ |
| 69 |
|
| 70 |
for (Int_t it=0; it<tof_event->ntrk() ; it++){ |
| 71 |
|
| 72 |
ToFTrkVar *tof = tof_event->GetToFTrkVar(it); |
| 73 |
|
| 74 |
/* trkseqno = index to get ToF alone or ToF+Tracker informations */ |
| 75 |
|
| 76 |
/* |
| 77 |
* trkseqno ==-1 -> ToF alone informations |
| 78 |
* trkseqno !=-1 correspond to the Tracker track index |
| 79 |
*/ |
| 80 |
|
| 81 |
// cout << " Index " << tof->trkseqno << endl; |
| 82 |
|
| 83 |
/* Get Tracker information if it exist */ |
| 84 |
|
| 85 |
Int_t trkgood = 0; |
| 86 |
|
| 87 |
if(tof->trkseqno != -1){ |
| 88 |
TrkTrack *trk = trk_event->GetStoredTrack(tof->trkseqno); |
| 89 |
|
| 90 |
/* Select "good" Tracker tracks" */ |
| 91 |
|
| 92 |
if(trk->chi2 > 0 && |
| 93 |
trk->chi2 < 100 && |
| 94 |
trk->GetNX() >= 4 && |
| 95 |
trk->GetNY() >= 3 && |
| 96 |
trk->xgood[0] && |
| 97 |
(trk->xgood[5] || trk->xgood[4]) && |
| 98 |
trk->al[4] != 0.) |
| 99 |
{ |
| 100 |
|
| 101 |
trkgood = 1; |
| 102 |
|
| 103 |
} |
| 104 |
|
| 105 |
} |
| 106 |
|
| 107 |
|
| 108 |
// Beta measurement based on ToF alone |
| 109 |
|
| 110 |
/* |
| 111 |
* beta is an array with 13 elements; |
| 112 |
* 12 correspond to the several planes combinations as described below; |
| 113 |
* the last one is a mean value |
| 114 |
*/ |
| 115 |
|
| 116 |
/* |
| 117 |
* |
| 118 |
* S11 - S31 = tof->beta[0] |
| 119 |
* S11 - S32 = tof->beta[1] |
| 120 |
* S12 - S31 = tof->beta[2] |
| 121 |
* S12 - S32 = tof->beta[3] |
| 122 |
* S21 - S31 = tof->beta[4] |
| 123 |
* S21 - S32 = tof->beta[5] |
| 124 |
* S22 - S31 = tof->beta[6] |
| 125 |
* S22 - S32 = tof->beta[7] |
| 126 |
* S11 - S21 = tof->beta[8] |
| 127 |
* S11 - S22 = tof->beta[9] |
| 128 |
* S12 - S21 = tof->beta[10] |
| 129 |
* S12 - S22 = tof->beta[11] |
| 130 |
* |
| 131 |
*/ |
| 132 |
|
| 133 |
|
| 134 |
Int_t ll = 0; |
| 135 |
Int_t mm = 0; |
| 136 |
|
| 137 |
/* |
| 138 |
* to get the same information for tracks from Tracker set: |
| 139 |
* if(tof->trkseqno==0) for the first tracker track; |
| 140 |
* if(tof->trkseqno==1) for the second tracker track; |
| 141 |
*/ |
| 142 |
|
| 143 |
if(tof->trkseqno == -1){ |
| 144 |
for (Int_t k=0; k<13 ; k++){ |
| 145 |
Float_t Beta = tof->beta[k]; |
| 146 |
// cout << "Beta " << Beta << endl; |
| 147 |
} |
| 148 |
} |
| 149 |
|
| 150 |
|
| 151 |
/* |
| 152 |
* read the pmtid for the tdc used in evaluate beta |
| 153 |
* for the first tracker track |
| 154 |
*/ |
| 155 |
|
| 156 |
Float_t tdcflagm[4][12]; |
| 157 |
|
| 158 |
for (Int_t j=0; j<4; j++){ |
| 159 |
for (Int_t k=0; k<12; k++){ |
| 160 |
tdcflagm[j][k]=0; |
| 161 |
} |
| 162 |
} |
| 163 |
|
| 164 |
|
| 165 |
|
| 166 |
if(tof->trkseqno == 0){ |
| 167 |
// npmttdc = number of pmt used for beta evaluation |
| 168 |
for (Int_t ipmtt=0; ipmtt<tof->npmttdc; ipmtt++){ |
| 169 |
// you can read the pmtid for the tdc used in evaluate beta |
| 170 |
Int_t pmttdc = tof->pmttdc[ipmtt]; |
| 171 |
// printf("pmttdc %d\n",pmttdc); |
| 172 |
Int_t tdcflag = tof->tdcflag[ipmtt]; |
| 173 |
if(tdcflag!=0) { |
| 174 |
// cout << " TDC " << tof->tdcflag[ipmtt] << " " << tdcflag << " " << pmttdc << endl; |
| 175 |
ll=0; |
| 176 |
mm=0; |
| 177 |
tof_event->GetPMTIndex(pmttdc,ll,mm); |
| 178 |
cout<<"pmttdc tdcflag "<<pmttdc<<" "<<tdcflag<<endl; |
| 179 |
tdcflagm[mm][ll]= tdcflag; |
| 180 |
} |
| 181 |
|
| 182 |
/* |
| 183 |
* Method to get Plane Index from a pmtid; |
| 184 |
* here is shown for a pmtid corresponding to |
| 185 |
* the tdc used in beta evaluation |
| 186 |
*/ |
| 187 |
Int_t pmt_plane=tof_event->GetPlaneIndex(pmttdc); |
| 188 |
|
| 189 |
} |
| 190 |
} |
| 191 |
|
| 192 |
|
| 193 |
ll = 0; |
| 194 |
mm = 0; |
| 195 |
|
| 196 |
/* |
| 197 |
* dEdx for each PMT with signal; |
| 198 |
* this value coprresponds to the signal from a single pmt |
| 199 |
*/ |
| 200 |
|
| 201 |
Float_t adc2[4][12]; |
| 202 |
|
| 203 |
if(tof->trkseqno == 0){ |
| 204 |
// npmtadc = number of pmt used for dEdx evaluation |
| 205 |
for (Int_t ipmta=0; ipmta<tof->npmtadc; ipmta++){ |
| 206 |
// dedx = adc_c |
| 207 |
Float_t dEdx = tof->dedx[ipmta]; |
| 208 |
// printf("dEdx %f %f\n",tof->dedx[ipmta],dEdx); |
| 209 |
|
| 210 |
// you can read the pmtid for the adc used in evaluate dedx=adc_c |
| 211 |
Int_t pmtadc = tof->pmtadc[ipmta]; |
| 212 |
// printf("pmtadc %d\n",tof->pmtadc[ipmta]); |
| 213 |
|
| 214 |
Int_t adcflag = tof->adcflag[ipmta]; |
| 215 |
// cout << " ADC " << adcflag << " " << pmtadc << " " << dEdx << endl; |
| 216 |
|
| 217 |
|
| 218 |
if(ipmta == 1) { |
| 219 |
ADChisto->Fill(tof->dedx[ipmta]); |
| 220 |
} |
| 221 |
|
| 222 |
|
| 223 |
/* |
| 224 |
to get pmt index in terms of matrix index |
| 225 |
from pmtid |
| 226 |
*/ |
| 227 |
tof_event->GetPMTIndex(pmtadc,ll,mm); |
| 228 |
adc2[mm][ll] = tof->dedx[ipmta]; |
| 229 |
|
| 230 |
|
| 231 |
/* |
| 232 |
to get pmt index from matrix index |
| 233 |
*/ |
| 234 |
tof_event->GetPMTid(ll,mm); |
| 235 |
Int_t ind = tof_event->GetPMTid(ll,mm); |
| 236 |
|
| 237 |
|
| 238 |
/* |
| 239 |
to get pmt name from pmtid |
| 240 |
*/ |
| 241 |
tof_event->GetPMTName(pmtadc); |
| 242 |
// cout << tof_event->GetPMTName(pmtadc) << endl; |
| 243 |
|
| 244 |
} |
| 245 |
} |
| 246 |
|
| 247 |
// Position measurement based on ToF hits |
| 248 |
|
| 249 |
/* |
| 250 |
* to get the same information for tracks from Tracker set: |
| 251 |
* if(tof->trkseqno==0) for the first tracker track; |
| 252 |
* if(tof->trkseqno==1) for the second tracker track; |
| 253 |
*/ |
| 254 |
|
| 255 |
if(tof->trkseqno==-1){ |
| 256 |
for (Int_t j=0; j<3; j++){ |
| 257 |
Float_t xtofpos = tof->xtofpos[j]; |
| 258 |
Float_t ytofpos = tof->ytofpos[j]; |
| 259 |
// printf("xtofpos %f %f %d\n",xtofpos,tof->xtofpos[j],ii); |
| 260 |
// printf("ytofpos %f %f %d\n",ytofpos,tof->ytofpos[j],ii); |
| 261 |
} |
| 262 |
} |
| 263 |
|
| 264 |
// Beta versus Rigidity |
| 265 |
|
| 266 |
/* |
| 267 |
* for the first tracker track tof->trkseqno == 0 |
| 268 |
* for the second tracker track tof->trkseqno == 1 |
| 269 |
* and so on... |
| 270 |
*/ |
| 271 |
|
| 272 |
if(tof->trkseqno == 0 && trkgood==1){ |
| 273 |
BetaRig->Fill( trk->GetRigidity(), tof->beta[4]); |
| 274 |
} |
| 275 |
|
| 276 |
|
| 277 |
/* |
| 278 |
* Method to get the dEdx on a ToF plane. |
| 279 |
* It is evaluated as sum of all adc signals on the plane |
| 280 |
* after corrections; its usage is not advised; |
| 281 |
* This method will be modified in the next ToFLevel2 release |
| 282 |
*/ |
| 283 |
|
| 284 |
if(tof->trkseqno == 0 && trkgood==1){ |
| 285 |
for (Int_t plane=0; plane<13; plane++){ |
| 286 |
Float_t MethoddEdx=tof_event->GetdEdx(it,plane); // it is the index for |
| 287 |
// tof_event->ntrk() |
| 288 |
} |
| 289 |
} |
| 290 |
|
| 291 |
|
| 292 |
// fill adc and tdc using getmatrix |
| 293 |
Float_t adc[4][12]; |
| 294 |
Float_t tdc[4][12]; |
| 295 |
tof_event->GetMatrix(it,adc,tdc); // it is the index for tof_event->ntrk$ |
| 296 |
|
| 297 |
for (Int_t j=0; j<4; j++){ |
| 298 |
for (Int_t k=0; k<12; k++){ |
| 299 |
Int_t ind = tof_event->GetPMTid(j,k); |
| 300 |
|
| 301 |
if (tdc[j][k] < 4095) cout<<ind<<" "<<tdc[j][k]<<" "<<tdcflagm[j][k]<<endl; |
| 302 |
} |
| 303 |
} |
| 304 |
|
| 305 |
cout<<"==================="<<endl; |
| 306 |
|
| 307 |
|
| 308 |
|
| 309 |
} // end loop over tracks |
| 310 |
|
| 311 |
Int_t nnpmt = tof_event->npmt(); |
| 312 |
|
| 313 |
Float_t adct[48]; |
| 314 |
Float_t tdct[48]; |
| 315 |
|
| 316 |
/* |
| 317 |
* ToFPMT is the PMT related class |
| 318 |
*/ |
| 319 |
|
| 320 |
|
| 321 |
for (Int_t ipmt=0; ipmt<tof_event->npmt() ; ipmt++){ |
| 322 |
/* |
| 323 |
tof_event->npmt() is the number of PMT with at least a TDC signal |
| 324 |
*/ |
| 325 |
|
| 326 |
ToFPMT *tofpmt = tof_event->GetToFPMT(ipmt); |
| 327 |
// for each PMT with a signal, tofpmt->adc gives the raw adc |
| 328 |
adct[ipmt]=tofpmt->adc; |
| 329 |
// for each PMT with a signal, tofpmt->tdc_tw gives the tdc after Time-walk correction |
| 330 |
tdct[ipmt]=tofpmt->tdc_tw; |
| 331 |
TestHisto->Fill(adct[2]); |
| 332 |
} |
| 333 |
|
| 334 |
|
| 335 |
/* |
| 336 |
* Some useful Methods (Examples): |
| 337 |
*/ |
| 338 |
|
| 339 |
Int_t nplane = 6; |
| 340 |
Int_t npaddle = 8; |
| 341 |
|
| 342 |
for (Int_t idpl=0; idpl< nplane; idpl++){ |
| 343 |
|
| 344 |
/* |
| 345 |
* Method to get the number of hit paddles on a ToF plane. |
| 346 |
* @param plane Plane ID (11, 12, 21, 22, 31, 32) or Plane index (0,1,2,3,4,5). |
| 347 |
* it is evaluated for each event in the run |
| 348 |
*/ |
| 349 |
Int_t NHit = tof_event->GetNHitPaddles(idpl); |
| 350 |
// printf("Number of Hit, idPlane %d %d %d %d\n",NHit,idpl,ii,it); |
| 351 |
|
| 352 |
/** |
| 353 |
* Method to get the plane ID (11 12 21 22 31 32) from the plane index (0 1 2 3 4 5) |
| 354 |
*/ |
| 355 |
Int_t ToFPlaneID=tof_event->GetToFPlaneID(idpl); |
| 356 |
// printf("ToFPlaneID, idPlane %d %d\n",ToFPlaneID,idpl); |
| 357 |
} |
| 358 |
|
| 359 |
|
| 360 |
/* |
| 361 |
* Method to know if a given ToF paddle was hit, that is there is a TDC signal |
| 362 |
* from both PMTs |
| 363 |
* @param plane Plane ID (11, 12, 21, 22, 31, 32) or Plane index (0,1,2,3,4,5). |
| 364 |
* @param paddle_id Paddle ID. |
| 365 |
* @return 1 if the paddle was hit. |
| 366 |
*/ |
| 367 |
|
| 368 |
for (Int_t idpad=0; idpad< npaddle; idpad++){ |
| 369 |
Int_t HitPad = tof_event->HitPaddle(idpl,idpad); |
| 370 |
// printf("Hit Pad %d %d %d\n",HitPad,idpl,idpad); |
| 371 |
} |
| 372 |
|
| 373 |
/** |
| 374 |
* Method to get the plane index (0 1 2 3 4 5) from the plane ID (11 12 21 22 31 32) |
| 375 |
*/ |
| 376 |
|
| 377 |
Int_t plane_id=11; |
| 378 |
Int_t ToFPlaneIndex = tof_event->GetToFPlaneIndex(plane_id); |
| 379 |
// printf("Plane Index %d %d\n",ToFPlaneIndex,plane_id); |
| 380 |
|
| 381 |
|
| 382 |
|
| 383 |
} // end loop over events |
| 384 |
|
| 385 |
TCanvas *BetaRigCanvas = new TCanvas("BetaRigCanvas","BetaRigCanvas", 700, 800); |
| 386 |
BetaRigCanvas->cd(); |
| 387 |
BetaRig->Draw(); |
| 388 |
|
| 389 |
TCanvas *SCanvas = new TCanvas("SCanvas","SCanvas", 700, 800); |
| 390 |
SCanvas->cd(); |
| 391 |
ADChisto->Draw(); |
| 392 |
|
| 393 |
TCanvas *TestCanvas = new TCanvas("TestCanvas","TestCanvas", 700, 800); |
| 394 |
TestCanvas->cd(); |
| 395 |
TestHisto->Draw(); |
| 396 |
} |