2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
19 //_________________________________________________________________________
20 // An analysis task to check the PMD data in simulated data
23 //////////////////////////////////////////////////////////////////////////////
34 #include "AliPMDQATask.h"
35 #include "AliPMDUtility.h"
39 //______________________________________________________________________________
40 AliPMDQATask::AliPMDQATask(const char *name) :
41 AliAnalysisTask(name,""),
84 // Input slot #0 works with an Ntuple
85 DefineInput(0, TChain::Class());
86 // Output slot #0 writes into a TH1 container
87 DefineOutput(0, TObjArray::Class()) ;
90 //______________________________________________________________________________
91 AliPMDQATask::~AliPMDQATask()
94 fOutputContainer->Clear() ;
95 delete fOutputContainer ;
138 //______________________________________________________________________________
139 void AliPMDQATask::ConnectInputData(const Option_t*)
141 // Initialisation of branch container and histograms
143 AliInfo(Form("*** Initialization of %s", GetName())) ;
146 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
148 AliError(Form("Input 0 for %s not found\n", GetName()));
152 // One should first check if the branch address was taken by some other task
153 char ** address = (char **)GetBranchAddress(0, "ESD");
155 fESD = (AliESD*)(*address);
158 SetBranchAddress(0, "ESD", &fESD);
162 //________________________________________________________________________
163 void AliPMDQATask::CreateOutputObjects()
167 fhPMDP1 = new TH2F("fhPMDP1","XY of Clusters",100,-100.,100.,100,-100.,100.);
168 fhPMDC2 = new TH1F("fhPMDC2","CPV PHI",200,-1,9);
169 fhPMDP2 = new TH1F("fhPMDP2","PRE PHI",200,-1,9);
170 fhPMDC3 = new TH1F("fhPMDC3","CPV Clus",30,0.,500.);
171 fhPMDP3 = new TH1F("fhPMDP3","PRE N-gammalike",20,0.,500.);
172 fhPMDP4 = new TH1F("fhPMDP4","PRE EDEP",30,0.,1000.);
173 fhPMDC5 = new TH1F("fhPMDC5","CPV n-cell",20,0.,100.);
174 fhPMDP5 = new TH1F("fhPMDP5","PMD n-cell",20,0.,100.);
175 fhPMDCP0 = new TH2F("fhPMDCP0","PRE CLUS Quad.1 vs 2",150,0.,300.,150,0.,300.);
176 fhPMDCP1 = new TH2F("fhPMDCP1","PRE CLUS Quad.3 vs 4",150,0.,300.,150,0.,300.);
177 fhPMDCP2 = new TH2F("fhPMDCP2","PRE EDEP Quad.3 vs 4",50,0.,300.,50,0.,300.);
178 fhPMDCP3 = new TH2F("fhPMDCP3","PRE EDEP vs Tot Clus ",10,0.,1000.,10,0.,300.);
179 fhPMDCP4 = new TH2F("fhPMDCP4","PRE Clus vs CPV Clus ",150,0.,200.,150,0.,200.);
181 fhPMDSM1 = new TH2F("fhPMDSM1","PRE Cluster XY",200,-100,100,200,-100,100);
182 fhPMDSM2 = new TH2F("fhPMDSM2","",999,-100.0,100.0,999,-100.0,100.0);
183 fhPMDSM3 = new TH2F("fhPMDSM3","",999,-100.0,100.0,999,-100.0,100.0);
184 fhPMDSM4 = new TH2F("fhPMDSM4","",999,-100.0,100.0,999,-100.0,100.0);
185 fhPMDSM5 = new TH2F("fhPMDSM5","",999,-100.0,100.0,999,-100.0,100.0);
186 fhPMDSM6 = new TH2F("fhPMDSM6","",999,-100.0,100.0,999,-100.0,100.0);
187 fhPMDSM7 = new TH2F("fhPMDSM7","",999,-100.0,100.0,999,-100.0,100.0);
188 fhPMDSM8 = new TH2F("fhPMDSM8","",999,-100.0,100.0,999,-100.0,100.0);
189 fhPMDSM9 = new TH2F("fhPMDSM9","",999,-100.0,100.0,999,-100.0,100.0);
190 fhPMDSM10 = new TH2F("fhPMDSM10","",999,-100.0,100.0,999,-100.0,100.0);
191 fhPMDSM11 = new TH2F("fhPMDSM11","",999,-100.0,100.0,999,-100.0,100.0);
192 fhPMDSM12 = new TH2F("fhPMDSM12","",999,-100.0,100.0,999,-100.0,100.0);
193 fhPMDSM13 = new TH2F("fhPMDSM13","",999,-100.0,100.0,999,-100.0,100.0);
194 fhPMDSM14 = new TH2F("fhPMDSM14","",999,-100.0,100.0,999,-100.0,100.0);
195 fhPMDSM15 = new TH2F("fhPMDSM15","",999,-100.0,100.0,999,-100.0,100.0);
196 fhPMDSM16 = new TH2F("fhPMDSM16","",999,-100.0,100.0,999,-100.0,100.0);
197 fhPMDSM17 = new TH2F("fhPMDSM17","",999,-100.0,100.0,999,-100.0,100.0);
198 fhPMDSM18 = new TH2F("fhPMDSM18","",999,-100.0,100.0,999,-100.0,100.0);
199 fhPMDSM19 = new TH2F("fhPMDSM19","",999,-100.0,100.0,999,-100.0,100.0);
200 fhPMDSM20 = new TH2F("fhPMDSM20","",999,-100.0,100.0,999,-100.0,100.0);
201 fhPMDSM21 = new TH2F("fhPMDSM21","",999,-100.0,100.0,999,-100.0,100.0);
202 fhPMDSM22 = new TH2F("fhPMDSM22","",999,-100.0,100.0,999,-100.0,100.0);
203 fhPMDSM23 = new TH2F("fhPMDSM23","",999,-100.0,100.0,999,-100.0,100.0);
204 fhPMDSM24 = new TH2F("fhPMDSM24","",999,-100.0,100.0,999,-100.0,100.0);
205 fhPMDSM = new TH1F("fhPMDSM","Plot of all 24 Super Modules",24,0,24);
207 // create output container
209 fOutputContainer = new TObjArray(38) ;
210 fOutputContainer->SetName("PMD") ;
212 fOutputContainer->AddAt(fhPMDP1, 0 );
213 fOutputContainer->AddAt(fhPMDC2, 1 );
214 fOutputContainer->AddAt(fhPMDP2, 2 );
215 fOutputContainer->AddAt(fhPMDC3, 3 );
216 fOutputContainer->AddAt(fhPMDP3, 4 );
217 fOutputContainer->AddAt(fhPMDP4, 5 );
218 fOutputContainer->AddAt(fhPMDC5, 6 );
219 fOutputContainer->AddAt(fhPMDP5, 7 );
220 fOutputContainer->AddAt(fhPMDCP0, 8 );
221 fOutputContainer->AddAt(fhPMDCP1, 9);
222 fOutputContainer->AddAt(fhPMDCP2, 10 );
223 fOutputContainer->AddAt(fhPMDCP3, 11 );
224 fOutputContainer->AddAt(fhPMDCP4, 12 );
226 fOutputContainer->AddAt(fhPMDSM1, 13 );
227 fOutputContainer->AddAt(fhPMDSM2, 14 );
228 fOutputContainer->AddAt(fhPMDSM3, 15 );
229 fOutputContainer->AddAt(fhPMDSM4, 16 );
230 fOutputContainer->AddAt(fhPMDSM5, 17 );
231 fOutputContainer->AddAt(fhPMDSM6, 18 );
232 fOutputContainer->AddAt(fhPMDSM7, 19 );
233 fOutputContainer->AddAt(fhPMDSM8, 20 );
234 fOutputContainer->AddAt(fhPMDSM9, 21 );
235 fOutputContainer->AddAt(fhPMDSM10, 22 );
236 fOutputContainer->AddAt(fhPMDSM11, 23 );
237 fOutputContainer->AddAt(fhPMDSM12, 24 );
238 fOutputContainer->AddAt(fhPMDSM13, 25 );
239 fOutputContainer->AddAt(fhPMDSM14, 26 );
240 fOutputContainer->AddAt(fhPMDSM15, 27 );
241 fOutputContainer->AddAt(fhPMDSM16, 28 );
242 fOutputContainer->AddAt(fhPMDSM17, 29 );
243 fOutputContainer->AddAt(fhPMDSM18, 30 );
244 fOutputContainer->AddAt(fhPMDSM19, 31 );
245 fOutputContainer->AddAt(fhPMDSM20, 32 );
246 fOutputContainer->AddAt(fhPMDSM21, 33 );
247 fOutputContainer->AddAt(fhPMDSM22, 34 );
248 fOutputContainer->AddAt(fhPMDSM23, 35 );
249 fOutputContainer->AddAt(fhPMDSM24, 36 );
250 fOutputContainer->AddAt(fhPMDSM, 37 );
253 //______________________________________________________________________________
254 void AliPMDQATask::Exec(Option_t *)
256 // Processing of one event
258 Long64_t entry = fChain->GetReadEntry() ;
261 AliError("fESD is not connected to the input!") ;
265 if ( !((entry-1)%100) )
266 AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ;
268 // ************************ PMD
270 AliPMDUtility *cc = new AliPMDUtility();
279 Float_t preCluQUAD[4] ;
280 Float_t cpvCluQUAD[4] ;
281 Float_t preADCQUAD[4] ;
282 Float_t cpvADCQUAD[4] ;
283 Float_t preCelQUAD[4] ;
284 Float_t cpvCelQUAD[4] ;
286 Int_t npmdCl = fESD->GetNumberOfPmdTracks();
288 // ****** The loop over PMD clusters
290 for (Int_t kk = 0; kk < 4 ; kk++) {
291 cpvCluQUAD[kk] = 0.0 ;
292 preCluQUAD[kk] = 0.0 ;
293 cpvCelQUAD[kk] = 0.0 ;
294 preCelQUAD[kk] = 0.0 ;
295 preADCQUAD[kk] = 0.0 ;
300 AliESDPmdTrack * pmdtr = fESD->GetPmdTrack(npmdCl);
301 Int_t det = pmdtr->GetDetector();
302 Float_t clsX = pmdtr->GetClusterX();
303 Float_t clsY = pmdtr->GetClusterY();
304 Float_t clsZ = pmdtr->GetClusterZ();
305 Float_t ncell = pmdtr->GetClusterCells();
306 Float_t adc = pmdtr->GetClusterADC();
308 cc->SetXYZ(clsX,clsY,clsZ);
311 Float_t eta = cc->GetEta();
312 Float_t phi = cc->GetPhi();
314 // Calculating S.Module Number from Cluster .
316 CalculateSMN(clsX, clsY, smn);
319 if(smn >= 0 && smn <= 5) {
321 cpvADCQUAD[0] =+ adc ;
322 cpvCelQUAD[0] =+ ncell ;
324 if(smn >= 6 && smn <=11) {
326 cpvADCQUAD[1] =+ adc ;
327 cpvCelQUAD[1] =+ ncell ;
329 if(smn >=12 && smn <=17) {
331 cpvADCQUAD[2] =+ adc ;
332 cpvCelQUAD[2] =+ ncell ;
334 if(smn >=18 && smn <=23) {
336 cpvADCQUAD[3] =+ adc ;
337 cpvCelQUAD[3] =+ ncell ;
340 if(eta >= 2.3 && eta <= 3.5)
347 if(smn >= 0 && smn <= 5) {
349 preADCQUAD[0] =+ adc ;
350 preCelQUAD[0] =+ ncell ;
352 if(smn >= 6 && smn <=11) {
354 preADCQUAD[1] =+ adc ;
355 preCelQUAD[1] =+ ncell ;
357 if(smn >=12 && smn <=17) {
359 preADCQUAD[2] =+ adc ;
360 preCelQUAD[2] =+ ncell ;
362 if(smn >=18 && smn <=23) {
364 preADCQUAD[3] =+ adc ;
365 preCelQUAD[3] =+ ncell ;
369 if(smn == 0) fhPMDSM1->Fill(-clsX,clsY);
370 if(smn == 0) fhPMDSM1->Fill(-clsX,clsY);
371 if(smn == 1) fhPMDSM2->Fill(-clsX,clsY);
372 if(smn == 2) fhPMDSM3->Fill(-clsX,clsY);
373 if(smn == 3) fhPMDSM4->Fill(-clsX,clsY);
374 if(smn == 4) fhPMDSM5->Fill(-clsX,clsY);
375 if(smn == 5) fhPMDSM6->Fill(-clsX,clsY);
376 if(smn == 6) fhPMDSM7->Fill(-clsX,clsY);
377 if(smn == 7) fhPMDSM8->Fill(-clsX,clsY);
378 if(smn == 8) fhPMDSM9->Fill(-clsX,clsY);
379 if(smn == 9) fhPMDSM10->Fill(-clsX,clsY);
380 if(smn ==10) fhPMDSM11->Fill(-clsX,clsY);
381 if(smn ==11) fhPMDSM12->Fill(-clsX,clsY);
382 if(smn ==12) fhPMDSM13->Fill(-clsX,clsY);
383 if(smn ==13) fhPMDSM14->Fill(-clsX,clsY);
384 if(smn ==14) fhPMDSM15->Fill(-clsX,clsY);
385 if(smn ==15) fhPMDSM16->Fill(-clsX,clsY);
386 if(smn ==16) fhPMDSM17->Fill(-clsX,clsY);
387 if(smn ==17) fhPMDSM18->Fill(-clsX,clsY);
388 if(smn ==18) fhPMDSM19->Fill(-clsX,clsY);
389 if(smn ==19) fhPMDSM20->Fill(-clsX,clsY);
390 if(smn ==20) fhPMDSM21->Fill(-clsX,clsY);
391 if(smn ==21) fhPMDSM22->Fill(-clsX,clsY);
392 if(smn ==22) fhPMDSM23->Fill(-clsX,clsY);
393 if(smn ==23) fhPMDSM24->Fill(-clsX,clsY);
395 if(eta >= 2.3 && eta <= 3.5)
399 fhPMDP1->Fill(clsX,clsY);
402 for (Int_t k = 0 ; k < 4 ; k++) {
403 totCPVClus =+ cpvCluQUAD [k] ;
404 totPREClus =+ preCluQUAD [k] ;
405 totCPVCell =+ cpvCelQUAD [k] ;
406 totPRECell =+ preCelQUAD [k] ;
407 totPREEdep =+ preADCQUAD [k] ;
409 Float_t totCPVpreClus = totPREClus + totCPVClus ;
411 // if(eta >= 2.3 && eta <= 3.5) {
412 fhPMDC3->Fill(totCPVClus);
413 fhPMDP3->Fill(totPREClus);
414 fhPMDP4->Fill(totPREEdep);
415 fhPMDP5->Fill(totPRECell);
416 fhPMDCP0->Fill(preCluQUAD[0],preCluQUAD[1]);
417 fhPMDCP1->Fill(preCluQUAD[2],preCluQUAD[3]);
418 fhPMDCP2->Fill(preADCQUAD[2],preADCQUAD[3]);
419 fhPMDCP3->Fill(totPREEdep,totCPVpreClus);
420 fhPMDCP4->Fill(totPREClus,totCPVClus);
428 PostData(0, fOutputContainer);
431 //______________________________________________________________________________
432 void AliPMDQATask::Terminate(Option_t *)
434 // Processing when the event loop is ended
435 fOutputContainer = (TObjArray*)GetOutputData(0);
437 fhPMDP1 = (TH2F*)fOutputContainer->At(0);
438 fhPMDC2 = (TH1F*)fOutputContainer->At(1);
439 fhPMDP2 = (TH1F*)fOutputContainer->At(2);
440 fhPMDC3 = (TH1F*)fOutputContainer->At(3);
441 fhPMDP3 = (TH1F*)fOutputContainer->At(4);
442 fhPMDP4 = (TH1F*)fOutputContainer->At(5);
443 fhPMDC5 = (TH1F*)fOutputContainer->At(6);
444 fhPMDP5 = (TH1F*)fOutputContainer->At(7);
445 fhPMDCP0 = (TH2F*)fOutputContainer->At(8);
446 fhPMDCP1 = (TH2F*)fOutputContainer->At(9);
447 fhPMDCP2 = (TH2F*)fOutputContainer->At(10);
448 fhPMDCP3 = (TH2F*)fOutputContainer->At(11);
449 fhPMDCP4 = (TH2F*)fOutputContainer->At(12);
451 fhPMDSM1 = (TH2F*)fOutputContainer->At(13);
452 fhPMDSM2 = (TH2F*)fOutputContainer->At(14);
453 fhPMDSM3 = (TH2F*)fOutputContainer->At(15);
454 fhPMDSM4 = (TH2F*)fOutputContainer->At(16);
455 fhPMDSM5 = (TH2F*)fOutputContainer->At(17);
456 fhPMDSM6 = (TH2F*)fOutputContainer->At(18);
457 fhPMDSM7 = (TH2F*)fOutputContainer->At(19);
458 fhPMDSM8 = (TH2F*)fOutputContainer->At(20);
459 fhPMDSM9 = (TH2F*)fOutputContainer->At(21);
460 fhPMDSM10 = (TH2F*)fOutputContainer->At(22);
461 fhPMDSM11 = (TH2F*)fOutputContainer->At(23);
462 fhPMDSM12 = (TH2F*)fOutputContainer->At(24);
463 fhPMDSM13 = (TH2F*)fOutputContainer->At(25);
464 fhPMDSM14 = (TH2F*)fOutputContainer->At(26);
465 fhPMDSM15 = (TH2F*)fOutputContainer->At(27);
466 fhPMDSM16 = (TH2F*)fOutputContainer->At(28);
467 fhPMDSM17 = (TH2F*)fOutputContainer->At(29);
468 fhPMDSM18 = (TH2F*)fOutputContainer->At(30);
469 fhPMDSM19 = (TH2F*)fOutputContainer->At(31);
470 fhPMDSM20 = (TH2F*)fOutputContainer->At(32);
471 fhPMDSM21 = (TH2F*)fOutputContainer->At(33);
472 fhPMDSM22 = (TH2F*)fOutputContainer->At(34);
473 fhPMDSM23 = (TH2F*)fOutputContainer->At(35);
474 fhPMDSM24 = (TH2F*)fOutputContainer->At(36);
475 fhPMDSM = (TH1F*)fOutputContainer->At(37);
477 Bool_t problem = kFALSE ;
478 AliInfo(Form(" *** %s Report:", GetName())) ;
480 gStyle->SetOptStat(110000);
481 gStyle->SetOptFit(1);
483 TCanvas *cPMD0 = new TCanvas("cPMD0","PMD ESD Test #1", 10,10, 600, 600);
484 cPMD0->Range(-100, -100,100 ,100 );
485 fhPMDSM1->SetMarkerColor(2);
487 fhPMDSM1->GetXaxis()->SetTitle("Cluster X");
488 fhPMDSM1->GetYaxis()->SetTitle("Cluster Y");
489 fhPMDSM2->SetMarkerColor(2);
490 fhPMDSM2->Draw("same");
491 fhPMDSM3->SetMarkerColor(2);
492 fhPMDSM3->Draw("same");
493 fhPMDSM4->SetMarkerColor(2);
494 fhPMDSM4->Draw("same");
495 fhPMDSM5->SetMarkerColor(2);
496 fhPMDSM5->Draw("same");
497 fhPMDSM6->SetMarkerColor(2);
498 fhPMDSM6->Draw("same");
499 fhPMDSM7->SetMarkerColor(4);
500 fhPMDSM7->Draw("same");
501 fhPMDSM8->SetMarkerColor(4);
502 fhPMDSM8->Draw("same");
503 fhPMDSM9->SetMarkerColor(4);
504 fhPMDSM9->Draw("same");
505 fhPMDSM10->SetMarkerColor(4);
506 fhPMDSM10->Draw("same");
507 fhPMDSM11->SetMarkerColor(4);
508 fhPMDSM11->Draw("same");
509 fhPMDSM12->SetMarkerColor(4);
510 fhPMDSM12->Draw("same");
511 fhPMDSM13->SetMarkerColor(6);
512 fhPMDSM13->Draw("same");
513 fhPMDSM14->SetMarkerColor(6);
514 fhPMDSM14->Draw("same");
515 fhPMDSM15->SetMarkerColor(6);
516 fhPMDSM15->Draw("same");
517 fhPMDSM16->SetMarkerColor(6);
518 fhPMDSM16->Draw("same");
519 fhPMDSM17->SetMarkerColor(6);
520 fhPMDSM17->Draw("same");
521 fhPMDSM18->SetMarkerColor(6);
522 fhPMDSM18->Draw("same");
523 fhPMDSM19->SetMarkerColor(8);
524 fhPMDSM19->Draw("same");
525 fhPMDSM20->SetMarkerColor(8);
526 fhPMDSM20->Draw("same");
527 fhPMDSM21->SetMarkerColor(8);
528 fhPMDSM21->Draw("same");
529 fhPMDSM22->SetMarkerColor(8);
530 fhPMDSM22->Draw("same");
531 fhPMDSM23->SetMarkerColor(8);
532 fhPMDSM23->Draw("same");
533 fhPMDSM24->SetMarkerColor(8);
534 fhPMDSM24->Draw("same");
537 DrawPMDBoundarySM1();
538 DrawPMDBoundarySM2();
539 DrawPMDBoundarySM3();
540 DrawPMDBoundarySM4();
541 cPMD0->Print("ClusterXY.eps");
543 TCanvas *cPMD1 = new TCanvas("cPMD1"," PMD ESD Test #2",10, 10, 600,600);
546 cPMD1->SetFillColor(0);
547 fhPMDC2->SetLineColor(4);
550 fhPMDP2->SetLineColor(2);
552 cPMD1->Print("CPVPREphi.eps");
554 TCanvas *cPMD2 = new TCanvas("cPMD2","PMD ESD test #3",10, 10, 600, 600);
556 fhPMDSM->SetFillColor(2);
558 cPMD2->Print("AllSMN.eps");
560 TCanvas *cPMD3 = new TCanvas("cPMD3", "PMD ESD test #4",10, 10, 600, 600);
563 fhPMDCP0->SetMarkerColor(9);
566 fhPMDCP1->SetMarkerColor(6);
569 fhPMDP3->SetLineColor(2);
572 fhPMDCP4->SetMarkerColor(3);
574 cPMD3->Print("CPVPREClus.eps");
576 TCanvas *cPMD4 = new TCanvas("cPMD4","PMD ESD test #5", 10, 10, 600, 600);
579 fhPMDC3->SetLineColor(4);
582 fhPMDP4->SetLineColor(2);
584 cPMD4->Print("CPVPREAdc.eps");
587 sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ;
588 gROOT->ProcessLine(line);
590 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
594 report="Problems found, please check!!!";
598 AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report)) ;
601 //______________________________________________________________________________
602 void AliPMDQATask::CalculateSMN( Float_t clsX, Float_t clsY, Int_t & smn)
604 Double_t xcon[96] = {75.133, 54.204, 53.254, 32.326, 31.376,10.447,
605 75.133, 54.204, 53.254, 32.326, 31.376,10.447,
606 75.133, 54.204, 53.254, 32.326, 31.376,10.447,
607 75.133, 54.204, 53.254, 32.326, 31.376,10.447,
608 -75.133, -54.204, -53.254, -32.326, -31.376,-10.447,
609 -75.133, -54.204, -53.254, -32.326, -31.376,-10.447,
610 -75.133, -54.204, -53.254, -32.326, -31.376,-10.447,
611 -75.133, -54.204, -53.254, -32.326, -31.376,-10.447,
612 9.167, -32.543, -33.493, -75.133,
613 9.167, -32.543, -33.493, -75.133,
614 9.167, -32.543, -33.493, -75.133,
615 9.167, -32.543, -33.493, -75.133,
616 9.167, -32.543, -33.493, -75.133,
617 9.167, -32.543, -33.493, -75.133,
618 -9.167, 32.543, 33.493, 75.133,
619 -9.167, 32.543, 33.493, 75.133,
620 -9.167, 32.543, 33.493, 75.133,
621 -9.167, 32.543, 33.493, 75.133,
622 -9.167, 32.543, 33.493, 75.133,
623 -9.167, 32.543, 33.493, 75.133};
625 Double_t ycon[96] = {86.475, 86.475, 86.475, 86.475, 86.475, 86.475,
626 38.225, 38.225, 38.225, 38.225, 38.225, 38.225,
627 37.325, 37.325, 37.325, 37.325, 37.325, 37.325,
628 -10.925, -10.925, -10.925, -10.925, -10.925, -10.925,
629 -86.475, -86.475, -86.475, -86.475, -86.475, -86.475,
630 -38.225, -38.225, -38.225, -38.225, -38.225, -38.225,
631 -37.325, -37.325, -37.325, -37.325, -37.325, -37.325,
632 10.925, 10.925, 10.925, 10.925, 10.925, 10.925,
633 86.475, 86.475, 86.475, 86.475,
634 62.225, 62.225, 62.225, 62.225,
635 61.325, 61.325, 61.325, 61.325,
636 37.075, 37.075, 37.075, 37.075,
637 36.175, 36.175, 36.175, 36.175,
638 11.925, 11.925, 11.925 , 11.925,
639 -86.475, -86.475, -86.475, -86.475,
640 -62.225, -62.225, -62.225, -62.225,
641 -61.325, -61.325, -61.325, -61.325,
642 -37.075, -37.075, -37.075, -37.075,
643 -36.175, -36.175, -36.175, -36.175,
644 -11.925, -11.925, -11.925 , -11.925 };
646 if((clsX <= xcon[0]) && (clsX >= xcon[1]) &&
647 (clsY <= ycon[0]) && (clsY >= ycon[6])) smn = 0 ;
649 else if((clsX <=xcon[2]) && (clsX >= xcon[3]) &&
650 (clsY <= ycon[1]) && (clsY >= ycon[7]))smn = 1 ;
652 else if((clsX <=xcon[4]) && (clsX >= xcon[5]) &&
653 (clsY <= ycon[3]) && (clsY >= ycon[8]))smn = 2 ;
655 else if((clsX <= xcon[0]) && (clsX >= xcon[1]) &&
656 (clsY <= ycon[12]) && (clsY >= ycon[18])) smn = 3 ;
658 else if((clsX <=xcon[2]) && (clsX >= xcon[3]) &&
659 (clsY <= ycon[12]) && (clsY >= ycon[18]))smn = 4 ;
661 else if((clsX <=xcon[4]) && (clsX >= xcon[5]) &&
662 (clsY <= ycon[12]) && (clsY >= ycon[18]))smn = 5 ;
663 //------------------------------------------------------------------
664 else if((clsX >= xcon[24]) && (clsX <= xcon[25]) &&
665 (clsY >= ycon[24]) && (clsY <= ycon[30])) smn = 6 ;
667 else if((clsX >=xcon[26]) && (clsX <= xcon[27]) &&
668 (clsY >= ycon[25]) && (clsY <= ycon[31]))smn = 7 ;
670 else if((clsX >=xcon[28]) && (clsX <= xcon[29]) &&
671 (clsY >= ycon[26]) && (clsY <= ycon[32]))smn = 8 ;
673 else if((clsX >= xcon[24]) && (clsX <= xcon[25]) &&
674 (clsY >= ycon[36]) && (clsY <= ycon[42])) smn = 9 ;
676 else if((clsX >=xcon[26]) && (clsX <= xcon[27]) &&
677 (clsY >= ycon[36]) && (clsY <= ycon[42]))smn = 10;
679 else if((clsX >=xcon[28]) && (clsX <= xcon[29]) &&
680 (clsY >= ycon[36]) && (clsY <= ycon[42]))smn = 11;
681 //------------------------------------------------------------------
682 else if((clsX <= xcon[48]) && (clsX >= xcon[49]) &&
683 (clsY <= ycon[48]) && (clsY >= ycon[52])) smn = 12 ;
685 else if((clsX <=xcon[50]) && (clsX >= xcon[51]) &&
686 (clsY <= ycon[48]) && (clsY >= ycon[52]))smn = 13 ;
688 else if((clsX <=xcon[48]) && (clsX >= xcon[49]) &&
689 (clsY <= ycon[56]) && (clsY >= ycon[60]))smn = 14 ;
691 else if((clsX <=xcon[50]) && (clsX >= xcon[51]) &&
692 (clsY <= ycon[56]) && (clsY >= ycon[60]))smn = 15 ;
694 else if((clsX <=xcon[48]) && (clsX >= xcon[49]) &&
695 (clsY <= ycon[64]) && (clsY >= ycon[68]))smn = 16 ;
697 else if((clsX <=xcon[50]) && (clsX >= xcon[51]) &&
698 (clsY <= ycon[64]) && (clsY >= ycon[68]))smn = 17 ;
699 //--------------------------------------------------------------
700 else if((clsX >= xcon[72]) && (clsX <= xcon[73]) &&
701 (clsY >= ycon[72]) && (clsY <= ycon[76])) smn = 18 ;
703 else if((clsX >=xcon[74]) && (clsX <= xcon[75]) &&
704 (clsY >= ycon[72]) && (clsY <= ycon[76]))smn = 19 ;
706 else if((clsX >=xcon[72]) && (clsX <= xcon[73]) &&
707 (clsY >= ycon[80]) && (clsY <= ycon[84]))smn = 20 ;
709 else if((clsX >=xcon[74]) && (clsX <= xcon[75]) &&
710 (clsY >= ycon[80]) && (clsY <= ycon[84]))smn = 21;
712 else if((clsX >= xcon[72]) && (clsX <= xcon[73]) &&
713 (clsY >= ycon[88]) && (clsY <= ycon[92])) smn = 22 ;
715 else if((clsX >=xcon[74]) && (clsX <= xcon[75]) &&
716 (clsY >= ycon[88]) && (clsY <= ycon[92]))smn = 23 ;
720 //______________________________________________________________________________
721 void AliPMDQATask::DrawPMDBoundary()
723 // Draw PMD boundaries
725 gStyle->SetLineWidth(2);
726 gStyle->SetLineColor(2);
728 l = new TLine(75.1333, 86.475, -75.1333, 86.475); l->Draw("same");
729 l = new TLine(-75.1333, 86.470,-75.1333, -86.475); l->Draw("same");
730 l = new TLine(-75.1333, -86.475,75.1333, -86.475); l->Draw("same");
731 l = new TLine(75.1333, -86.475,75.1333, 86.475); l->Draw("same");
734 //______________________________________________________________________________
735 void AliPMDQATask::DrawPMDBoundarySM1()
737 // Draw boundaries of Super Module 1
739 gStyle->SetLineWidth(1);
740 gStyle->SetLineColor(4);
742 l = new TLine(-75.1333, 86.475, -10.447, 86.475); l->Draw("same");
743 l = new TLine(-10.447, 86.475, -10.446, -10.925); l->Draw("same");
744 l = new TLine(-10.446, -10.925, -75.1333,-10.925); l->Draw("same");
745 l = new TLine(-75.1333,-10.925, -75.1333, 86.475); l->Draw("same");
748 //______________________________________________________________________________
749 void AliPMDQATask::DrawPMDBoundarySM2()
751 // Draw boundaries of Super Module 2
753 gStyle->SetLineWidth(1);
754 gStyle->SetLineColor(4);
756 l = new TLine(75.1333, -86.475, 10.446, -86.475); l->Draw("same");
757 l = new TLine(10.446, -86.475, 10.446, 10.925); l->Draw("same");
758 l = new TLine(10.446, 10.925, 75.1333, 10.925); l->Draw("same");
759 l = new TLine(75.1333, 10.925, 75.1333, -86.475); l->Draw("same");
763 //______________________________________________________________________________
764 void AliPMDQATask::DrawPMDBoundarySM3()
766 // Draw boundaries of Super Module 3
768 gStyle->SetLineWidth(1);
769 gStyle->SetLineColor(1);
771 l = new TLine( -9.167, 86.475, 75.1333, 86.475); l->Draw("same");
772 l = new TLine(75.1333,86.475, 75.1333, 11.925); l->Draw("same");
773 l = new TLine(75.1333,11.925, -9.167, 11.925); l->Draw("same");
774 l = new TLine( -9.167, 11.925, -9.167, 86.475); l->Draw("same");
777 //______________________________________________________________________________
778 void AliPMDQATask::DrawPMDBoundarySM4()
780 // Draw boundaries of Super Module 4
782 gStyle->SetLineWidth(1);
783 gStyle->SetLineColor(1);
785 l = new TLine(9.167, -86.475, -75.1333,-86.475); l->Draw("same");
786 l = new TLine(-75.1333,-86.475, -75.1333,-11.925); l->Draw("same");
787 l = new TLine(-75.1333,-11.925, 9.167, -11.925); l->Draw("same");
788 l = new TLine(9.167, -11.925, 9.167, -86.475); l->Draw("same");