]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSGeoPlot.C
DP:Misalignment of CPV added
[u/mrichter/AliRoot.git] / ITS / AliITSGeoPlot.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include<Riostream.h>
3 #include<TROOT.h>
4 #include<TArrayI.h>
5 #include<TBranch.h>
6 #include<TCanvas.h>
7 #include<TClassTable.h>
8 #include<TClonesArray.h>
9 #include<TFile.h>
10 #include<TGeoManager.h>
11 #include<TH1.h>
12 #include<TH2.h>
13 #include <TInterpreter.h>
14 #include<TObject.h>
15 #include<TObjArray.h>
16 #include<TTree.h>
17 #include "AliRun.h"
18 #include "AliITS.h"
19 #include "AliITSgeom.h"
20 #include "AliITSDetTypeRec.h"
21 #include "AliITSRecPoint.h"
22 #include "AliITSRecPoint.h"
23 #include "AliITSdigit.h"
24 #include "AliITSdigitSSD.h"
25 #include "AliITShit.h"
26 #include "AliITSmodule.h" 
27 #include "AliITSsegmentation.h"
28 #include "AliITSsegmentationSPD.h" 
29 #include "AliITSsegmentationSDD.h"
30 #include "AliITSsegmentationSSD.h"
31 #include "AliRunLoader.h"
32 #include "AliITSLoader.h"
33 #include "AliHeader.h"
34 #endif
35 void GetHitsCoor(TObject *its, Int_t mod, TObjArray & histos, Int_t subd,Bool_t verb);
36 Int_t GetRecCoor(TObject *ge, TClonesArray *ITSrec, Int_t mod, TH2F *h2, TH1F *h1, Bool_t verb);
37 Int_t GetClusCoor(TObject *ge, TClonesArray *ITSrec, Int_t mod, TH2F *h2, TH1F *h1, Bool_t verb);
38 void GetDigits(TObject *tmps,TObject *ge,TClonesArray *ITSdigits, Int_t subd, Int_t mod, Bool_t verbose, TObjArray & histos);
39
40 Int_t AliITSGeoPlot (Int_t evesel=0, char *opt="All+ClustersV2", TString filename="galice.root", Int_t isfastpoints = 0) {
41   /*******************************************************************
42    *  This macro displays geometrical information related to the
43    *  hits, digits and rec points (or V2 clusters) in ITS.
44    *  There are histograms that are not displayed (i.e. energy
45    *  deposition) but that are saved onto a file (see below)
46    *
47    *  INPUT arguments:
48    *
49    *  Int_t evesel:  the selected event number. It's =0 by default
50    *
51    *  Options: Any combination of:
52    *    1) subdetector name:  "SPD", "SDD", "SSD", "All" (default)
53    *    2) Printouts:        "Verbose" Almost everything is printed out: 
54    *                          it is wise to redirect the output onto a file
55    *                    e.g.: .x AliITSGeoPlot.C("All+Verbose") > out.log 
56    *    3) Rec Points option: "Rec"   ---> Uses Rec. Points (default)
57    * 
58    *    4) ClustersV2 option: "ClustersV2" ---->Uses ClustersV2
59    *                           otherwise ---->uses hits and digits only
60    *    Examples:
61    *       .x AliITSGeoPlot();  (All subdetectors; no-verbose; no-recpoints)
62    *       .x AliITSGeoPlot("SPD+SSD+Verbose+Rec"); 
63    *   
64    *    filename:   It's "galice.root" by default. 
65    *    isfastpoints: integer. It is set to 0 by defaults. This means that
66    *                slow reconstruction is assumed. If fast recpoint are in
67    *                in use, isfastpoints must be set =1.
68    *
69    *  OUTPUT: It produces a root file with a list of histograms
70    *    
71    *  WARNING: spatial information for SSD/DIGITS is obtained by pairing
72    *          digits on p and n side originating from the same track, when
73    *          possible. This (mis)use of DIGITS is tolerated for debugging 
74    *          purposes only !!!!  The pairing in real life should be done
75    *          starting from info really available... 
76    * 
77    *  COMPILATION: this macro can be compiled. 
78    *      1)       You need to set your include path with
79    * gSystem->SetIncludePath("-I- -I$ALICE_ROOT/include -I$ALICE_ROOT/ITS -g");
80    *      3)       If you are using root instead of aliroot you need to
81    *               execute the macro loadlibs.C in advance
82    *      4)       To compile this macro from root (or aliroot):
83    *                 ---  .L AliITSGeoPlot.C++
84    *                 ---  AliITSGeoPlot();
85    *     
86    *  M.Masera  14/05/2001 18:30
87    *  Last rev. 31/05/2004 14:00 (Clusters V2 added)  E.C.          
88    ********************************************************************/
89
90   //Options
91   TString choice(opt);
92   Bool_t All = choice.Contains("All");
93   Bool_t verbose=choice.Contains("Verbose");
94   Bool_t userec=choice.Contains("Rec");
95   Bool_t useclustersv2=choice.Contains("ClustersV2");
96   Int_t retcode=1; //return code
97  
98   if (gClassTable->GetID("AliRun") < 0) {
99     gInterpreter->ExecuteMacro("loadlibs.C");
100   }
101   else { 
102     if(gAlice){
103       delete gAlice->GetRunLoader();
104       delete gAlice;
105       gAlice=0;
106     }
107   }
108   // retrives geometry 
109   TString geof(gSystem->DirName(filename));
110   geof += "/geometry.root";
111   TGeoManager::Import(geof.Data());
112   if (!gGeoManager) {
113     cout<<"geometry not found\n";
114     return -1;
115   }
116  
117   AliRunLoader* rl = AliRunLoader::Open(filename.Data());
118   if (rl == 0x0){
119     cerr<<"AliITSGeoPlot.C : Can not open session RL=NULL"<< endl;
120     return -1;
121   }
122   Int_t retval = rl->LoadgAlice();
123   if (retval){
124     cerr<<"AliITSGeoPlot.C : LoadgAlice returned error"<<endl;
125     return -1;
126   }
127   gAlice=rl->GetAliRun();
128
129   retval = rl->LoadHeader();
130   if (retval){
131     cerr<<"AliITSGeoPlot.C : LoadHeader returned error"<<endl;
132     return -1;
133   }
134
135   AliITSLoader* ITSloader =  (AliITSLoader*) rl->GetLoader("ITSLoader");
136
137   if(!ITSloader){
138     cerr<<"AliITSGeoPlot.C :  ITS loader not found"<<endl;
139     return -1;
140   }
141
142   ITSloader->LoadHits("read");
143   ITSloader->LoadDigits("read");
144   if(isfastpoints==1)ITSloader->SetRecPointsFileName("ITS.FastRecPoints.root");
145   ITSloader->LoadRecPoints("read");
146   rl->GetEvent(evesel);
147   Int_t nparticles = rl->GetHeader()->GetNtrack();
148   AliITS *ITS  = (AliITS*)gAlice->GetModule("ITS");
149   ITS->SetTreeAddress();
150   if(verbose) {
151     cout<<" "<<endl<<" "<<endl;
152     cout<<"******* Event processing started   *******"<<endl;
153     cout<<"In the following, hits with 'StatusEntering' flag active"<<endl;
154     cout<<"will not be processed"<<endl; 
155     cout << "Number of particles=  " << nparticles <<endl;
156   }
157
158   // HITS
159   TTree *TH = ITSloader->TreeH();
160   Stat_t ntracks = TH->GetEntries();
161   if(verbose)cout<<"Number of primary tracks= "<<ntracks<<endl;
162
163   // ITS
164   Int_t nmodules;
165   ITS->InitModules(-1,nmodules);
166   cout<<"Number of ITS modules= "<<nmodules<<endl;
167   cout<<"Filling modules... It takes a while, now. Please be patient"<<endl;
168   ITS->FillModules(0,0,nmodules," "," ");
169   cout<<"ITS modules .... DONE!"<<endl;
170   
171   AliITSDetTypeRec* detTypeRec = new AliITSDetTypeRec(ITSloader);
172   detTypeRec->SetDefaults();
173   // DIGITS
174   TTree *TD = ITSloader->TreeD();
175
176   //RECPOINTS (V2 clusters)
177   TTree *TR = ITSloader->TreeR();
178   TClonesArray *ITSrec  = detTypeRec->RecPoints();
179   TBranch *branch = 0;
180   if(userec && TR && ITSrec){
181     if(isfastpoints==1){
182       branch = ITSloader->TreeR()->GetBranch("ITSRecPointsF");
183       cout<<"using fast points\n";
184     }
185     else {
186       branch = ITSloader->TreeR()->GetBranch("ITSRecPoints");
187     }
188     if(branch)branch->SetAddress(&ITSrec);
189   }
190
191   if(userec && (!TR || !ITSrec || !branch)){
192     userec = kFALSE;
193     cout<<"\n ======================================================= \n";
194     cout<<"WARNING: there are no RECPOINTS on this file ! \n";
195     cout<<"======================================================= \n \n";
196   }
197   if(useclustersv2 && TR && ITSrec){
198     branch = ITSloader->TreeR()->GetBranch("ITSRecPoints");
199     if(branch)branch->SetAddress(&ITSrec);
200   }
201
202   if(useclustersv2 && (!TR || !ITSrec || !branch)){
203     useclustersv2 = kFALSE;
204     cout<<"\n ======================================================= \n";
205     cout<<"WARNING: there are no CLUSTERSV2 on this file ! \n";
206     cout<<"======================================================= \n \n";
207   }
208
209
210   //local variables
211   Int_t mod;   //module number
212   Int_t nbytes = 0; 
213   Double_t pos[3];  // Global position of the current module
214   Float_t ragdet; // Radius of detector (x y plane)
215   Int_t first; // first module
216   Int_t last; // last module
217   Int_t nrecp; //number of RecPoints for a given module
218
219   //List of histograms
220   TObjArray histos(26,0);  // contains the pointers to the histograms
221   // Book histograms SPD
222   TH2F *rspd = new TH2F("rspd","Radii of digits - SPD",50,-10.,10.,50,-10.,10.);
223   TH2F *rhspd = new TH2F("rhspd","Radii of hits - SPD",50,-10.,10.,50,-10.,10.);
224   TH2F *rmspd = new TH2F("rmspd","Radii of SPD modules",50,-10.,10.,50,-10.,10.);
225   TH1F *zspd = new TH1F("zspd","Z of digits - SPD",100,-30.,30.);
226   TH1F *zhspd = new TH1F("zhspd","Z of hits - SPD",100,-30.,30.);
227   TH1F *zmspd = new TH1F("zmspd","Z of SPD modules",100,-30,30.);
228
229   Char_t title1[50]="";
230   Char_t title2[50]="";
231   if(userec){ 
232     sprintf(title1,"Radii of recpoints - %s","SPD");
233     sprintf(title2,"Z of recpoints - %s","SPD");
234   }
235   if(useclustersv2){
236     sprintf(title1,"Radii of clustersV2 - %s","SPD");
237     sprintf(title2,"Z of clustersV2 - %s","SPD");
238   }
239   TH2F *rrspd = new TH2F("rrspd",title1,50,-10.,10.,50,-10.,10.);
240   TH1F *zrspd = new TH1F("zrspd",title2,100,-30.,30.);
241   TH1F *enespd = new TH1F("enespd","Energy deposition SPD (KeV)",100,0.,1000.);
242   histos.AddLast(rspd);  // 0
243   histos.AddLast(rhspd); // 1
244   histos.AddLast(rmspd); // 2
245   histos.AddLast(zspd);  // 3
246   histos.AddLast(zhspd); // 4
247   histos.AddLast(zmspd); // 5
248   histos.AddLast(rrspd); // 6
249   histos.AddLast(zrspd); // 7
250   histos.AddLast(enespd); // 8
251   // Book histograms SDD
252   TH2F *rsdd = new TH2F("rsdd","Radii of digits - SDD",50,-40.,40.,50,-40.,40.);
253   TH2F *rhsdd = new TH2F("rhsdd","Radii of hits - SDD",50,-40.,40.,50,-40.,40.);
254   TH2F *rmsdd = new TH2F("rmsdd","Radii of SDD modules",50,-40.,40.,50,-40.,40.);
255   TH1F *zsdd = new TH1F("zsdd","Z of digits - SDD",100,-40.,40.);
256   TH1F *zhsdd = new TH1F("zhsdd","Z of hits - SDD",100,-40.,40.);
257   TH1F *zmsdd = new TH1F("zmsdd","Z of SDD modules",100,-40,40.);
258   Char_t title3[50];
259   Char_t title4[50];
260   if(userec){ 
261     sprintf(title3,"Radii of recpoints - %s","SDD");
262     sprintf(title4,"Z of recpoints - %s","SDD");
263   }
264   if(useclustersv2){
265     sprintf(title3,"Radii of clustersV2 - %s","SDD");
266     sprintf(title4,"Z of clustersV2 - %s","SDD");
267   }
268   TH2F *rrsdd = new TH2F("rrsdd",title3,50,-40.,40.,50,-40.,40.);   
269   TH1F *zrsdd = new TH1F("zrsdd",title4,100,-40.,40.);
270   TH1F *enesdd = new TH1F("enesdd","Energy deposition SDD (KeV)",100,0.,1000.);
271   histos.AddLast(rsdd);  // 9
272   histos.AddLast(rhsdd); // 10
273   histos.AddLast(rmsdd); // 11
274   histos.AddLast(zsdd);  // 12
275   histos.AddLast(zhsdd); // 13
276   histos.AddLast(zmsdd); // 14
277   histos.AddLast(rrsdd); // 15
278   histos.AddLast(zrsdd); // 16
279   histos.AddLast(enesdd); // 17
280   // Book histogram SSD
281   TH2F *rssd = new TH2F("rssd","Radii of digits - SSD",50,-50.,50.,50,-50.,50.);
282   TH2F *rhssd = new TH2F("rhssd","Radii of hits - SSD",50,-50.,50.,50,-50.,50.);
283   TH2F *rmssd = new TH2F("rmssd","Radii of SSD modules",50,-50.,50.,50,-50.,50.);
284   TH1F *zssd = new TH1F("zssd","Z of digits - SSD",100,-70.,70.);
285   TH1F *zhssd = new TH1F("zhssd","Z of hits - SSD",100,-70.,70.);
286   TH1F *zmssd = new TH1F("zmssd","Z of SSD modules",100,-70,70.);
287   Char_t title5[50];
288   Char_t title6[50];
289   if(userec){ 
290     sprintf(title5,"Radii of recpoints - %s","SSD");
291     sprintf(title6,"Z of recpoints - %s","SSD");
292   }
293   if(useclustersv2){
294     sprintf(title5,"Radii of clustersV2 - %s","SSD");
295     sprintf(title6,"Z of clustersV2 - %s","SSD");
296   }
297
298   TH2F *rrssd = new TH2F("rrssd",title5,50,-50.,50.,50,-50.,50.);
299   TH1F *zrssd = new TH1F("zrssd",title6,100,-70.,70.);
300   TH1F *enessd = new TH1F("enessd","Energy deposition SSD (KeV)",100,0.,1000.);
301   histos.AddLast(rssd);  // 18
302   histos.AddLast(rhssd); // 19
303   histos.AddLast(rmssd); // 20
304   histos.AddLast(zssd);  // 21
305   histos.AddLast(zhssd); // 22
306   histos.AddLast(zmssd); // 23
307   histos.AddLast(rrssd); // 24
308   histos.AddLast(zrssd); // 25
309   histos.AddLast(enessd); // 26
310   //
311   // Loop on subdetectors
312   // 
313   cout<<"CALL GETITSGEOM \n \n \n";
314   AliITSgeom *geom = ITS->GetITSgeom();
315   TString detna; // subdetector name
316   for(Int_t subd=0;subd<3;subd++){
317     if(All || (choice.Contains("SPD") && subd==0) || (choice.Contains("SDD") && subd==1) || (choice.Contains("SSD") && subd==2)){
318       // Prepare array for the digits
319       TClonesArray *ITSdigits  = ITS->DigitsAddress(subd);
320       Bool_t usedigits = kTRUE;
321       if(!ITSdigits){
322         usedigits = kFALSE;
323         cout<<"\n ======================================================= \n";
324         cout<<"WARNING: there are no DIGITS on this file ! \n";
325         cout<<"======================================================= \n \n";
326       }
327       // Get segmentation model
328       if(subd==0)detna="SPD";
329       if(subd==1)detna="SDD";
330       if(subd==2)detna="SSD";
331       AliITSsegmentation *seg=(AliITSsegmentation*)detTypeRec->GetSegmentationModel(subd);
332       // Loop on modules
333       first = geom->GetStartDet(subd);
334       last = geom->GetLastDet(subd);
335       if(verbose){
336         cout<<"     "<<endl<<"-------------------------------------"<<endl;
337         cout<<"Start processing subdetector "<<detna<<endl;
338         cout<<detna<<" first module "<<first<<endl;
339         cout<<detna<<" last module "<<last<<endl;
340         cout<<" "<<endl<<" "<<endl;
341       }
342       for (mod=first; mod<=last; mod++){
343         geom->GetTrans(mod,pos);  // position of the module in the MRS
344         ragdet=sqrt(pos[0]*pos[0]+pos[1]*pos[1]);
345         // The following 2 histos are a check of the geometry
346         TH2F *bidi = (TH2F*)histos.At(2+subd*9);
347         TH1F *uni = (TH1F*)histos.At(5+subd*9);
348         bidi->Fill(pos[0],pos[1]);
349         uni->Fill(pos[2]);
350         if(verbose){
351           cout<<"=========================================================\n";
352           cout<<detna<<" module="<<mod<<endl;
353           cout<<"Mod. coordinates: "<<pos[0]<<", "<<pos[1]<<", ";
354           cout<<pos[2]<<" Radius "<<ragdet<<endl;
355         }
356
357         // Hits
358         GetHitsCoor(ITS,mod,histos,subd,verbose);
359
360         //RecPoints     
361         if(userec){
362           detTypeRec->ResetRecPoints();
363           branch->GetEvent(mod);
364           TH2F *bidi=(TH2F*)histos.At(6+subd*9);
365           TH1F *uni=(TH1F*)histos.At(7+subd*9);
366           nrecp=GetRecCoor(geom,ITSrec,mod,bidi,uni,verbose);
367         }
368         if(useclustersv2){
369           detTypeRec->ResetRecPoints();
370           branch->GetEvent(mod);
371           TH2F *bidi=(TH2F*)histos.At(6+subd*9);
372           TH1F *uni=(TH1F*)histos.At(7+subd*9);
373           nrecp=GetClusCoor(geom,ITSrec,mod,bidi,uni,verbose);
374           
375         }
376      
377         // Digits
378         if(usedigits){
379           detTypeRec->ResetDigits();
380           nbytes += TD->GetEvent(mod);
381           GetDigits(seg,geom,ITSdigits,subd,mod,verbose,histos);
382         }
383
384       } // End of loop on the modules
385       TH1F *h1tmp;
386       TH2F *h2tmp;
387       //  Plot the histograms
388       TCanvas *current=0; // current Canvas (1--> SPD, 2---> SDD, 3---> SSD)
389       if(subd==0){
390         // Prepare canvas 1
391         TCanvas *c1 = new TCanvas("c1","SPD",10,10,600,900);
392         c1->Divide(2,3);
393         current=c1;
394       }
395       if(subd==1){
396         // Prepare canvas 2
397         TCanvas *c2 = new TCanvas("c2","SDD",40,40,600,900);
398         c2->Divide(2,3);
399         current=c2;
400       }
401       if(subd==2){
402         // Prepare canvas 3
403         TCanvas *c3 = new TCanvas("c3","SSD",70,70,600,900);
404         c3->Divide(2,3);
405         current=c3;
406       }
407       current->cd(1);
408       h2tmp = (TH2F*)histos.At(9*subd);
409       h2tmp->Draw();
410       current->cd(2);
411       h1tmp=(TH1F*)histos.At(3+subd*9);
412       h1tmp->Draw();
413       current->cd(3);
414       h2tmp=(TH2F*)histos.At(1+9*subd);
415       h2tmp->Draw();
416       current->cd(4);
417       h1tmp=(TH1F*)histos.At(4+subd*9);
418       h1tmp->Draw();
419    
420       if(userec || useclustersv2){
421         current->cd(5);
422         h2tmp=(TH2F*)histos.At(6+9*subd);
423         h2tmp->Draw();
424         current->cd(6);
425         h1tmp=(TH1F*)histos.At(7+subd*9);
426         h1tmp->Draw();
427       }
428   
429       else {
430         current->cd(5);
431         h2tmp=(TH2F*)histos.At(2+9*subd);
432         h2tmp->Draw();
433         current->cd(6);
434         h2tmp=(TH2F*)histos.At(5+9*subd);
435         h2tmp->Draw();
436       }
437     } // if(All.....
438   } // end of loop on subdetectors
439   // Save the histograms
440   TFile *fh = new TFile("AliITSGeoPlot.root","recreate");
441   // The list is written to file as a single entry
442   TList *lihi = new TList();
443   // copy the pointers to the histograms to a TList object.
444   // The histograms concerning recpoints are not copied if
445   // 'userec' is false.
446   for(Int_t i=0;i<histos.GetEntriesFast();i++){
447     if(choice.Contains("All") || (choice.Contains("SPD") && i<8) || (choice.Contains("SDD") && i>7 && i<16) || (choice.Contains("SSD") && i>15)){
448       if(!(!userec && ((i+2)%9==0 || (i+1)%9==0)))lihi->Add(histos.At(i));
449     }
450   }
451   lihi->Write("Histograms ITS hits+digits+recpoints",TObject::kSingleKey);
452   fh->Close();
453
454   return retcode;
455 }
456
457
458 void GetHitsCoor(TObject *its, Int_t mod, TObjArray & histos, Int_t subd,Bool_t verb){
459   TH2F *h2=(TH2F*)histos.At(1+subd*9);
460   TH1F *h1=(TH1F*)histos.At(4+subd*9);
461   TH1F *ener=(TH1F*)histos.At(8+subd*9);
462   AliITS *ITS= (AliITS*)its;
463   AliITSmodule *modu = ITS->GetModule(mod);
464   TObjArray *fHits = modu->GetHits();
465   Int_t nhits = fHits->GetEntriesFast();
466   if(nhits>0){
467     if(verb){
468       cout<<"-------------------------------------------------------"<<endl;
469       cout<<"Number of HITS for module "<<mod<<": "<<nhits<<endl;
470     }
471     for (Int_t hit=0;hit<nhits;hit++) {
472       AliITShit *iHit = (AliITShit*) fHits->At(hit);
473       if(!iHit->StatusEntering()){
474         Float_t x=iHit->GetXG();
475         Float_t y=iHit->GetYG();
476         Float_t z=iHit->GetZG();
477         Float_t edep=iHit->GetIonization()*1000000;
478         h2->Fill(x,y);
479         h1->Fill(z);
480         ener->Fill(edep);
481         if(verb){
482           cout<<"hit # "<<hit<<" Global coordinates "<<x<<" "<<y<<" "<<z<<endl;
483           cout<<"track # "<<iHit->GetTrack()<<" energy deposition (KeV)= ";
484           cout<<edep<<endl;
485         }
486       }
487     }
488   }
489 }
490
491
492 Int_t GetClusCoor(TObject *ge, TClonesArray *ITSrec, Int_t mod, TH2F *h2, TH1F *h1, Bool_t verb){
493
494   AliITSgeom *geom = (AliITSgeom*)ge;
495   Int_t nrecp = ITSrec->GetEntries();
496   if(nrecp>0){
497     Float_t lc[3]; for(Int_t i=0; i<3; i++) lc[i]=0.;
498     Float_t gc[3]; for(Int_t i=0; i<3; i++) gc[i]=0.;
499     if(verb){
500       cout<<"-------------------------------------------------------"<<endl;
501       cout<<"Number of CLUSTERS for module "<<mod<<": "<<nrecp<<endl;
502     }
503     for(Int_t irec=0;irec<nrecp;irec++) {
504       AliITSRecPoint *recp = (AliITSRecPoint*)ITSrec->At(irec);
505       Double_t rot[9];     
506       geom->GetRotMatrix(mod,rot);
507       Int_t lay,lad,det;   
508       geom->GetModuleId(mod,lay,lad,det);
509       Float_t tx,ty,tz;    
510       geom->GetTrans(lay,lad,det,tx,ty,tz);     
511
512       Double_t alpha=TMath::ATan2(rot[1],rot[0])+TMath::Pi();
513       Double_t phi1=TMath::Pi()/2+alpha;
514       if(lay==1) phi1+=TMath::Pi();
515
516       Float_t cp=TMath::Cos(phi1), sp=TMath::Sin(phi1);
517       Float_t  r=tx*cp+ty*sp;
518       gc[0]= r*cp - recp->GetY()*sp;
519       gc[1]= r*sp + recp->GetY()*cp;
520       gc[2]=recp->GetZ();
521   
522       if(verb){
523         Float_t r=sqrt(gc[0]*gc[0]+gc[1]*gc[1]);
524         cout<<"Global coor. + radius "<<gc[0]<<" "<<gc[1]<<" "<<gc[2]<<" ";
525         cout<<r<<endl;
526         cout<<"Associated track "<<recp->GetLabel(0)<<endl;
527       }
528       h2->Fill(gc[0],gc[1]);
529       h1->Fill(gc[2]);
530
531     }
532   }
533   return nrecp;
534 }
535 Int_t GetRecCoor(TObject *ge, TClonesArray *ITSrec, Int_t mod, TH2F *h2, TH1F *h1, Bool_t verb){
536   AliITSgeom *geom = (AliITSgeom*)ge;
537   Int_t nrecp = ITSrec->GetEntries();
538   if(nrecp>0){
539     Float_t lc[3]; for(Int_t i=0; i<3; i++) lc[i]=0.;
540     Float_t gc[3]; for(Int_t i=0; i<3; i++) gc[i]=0.;
541     if(verb){
542       cout<<"-------------------------------------------------------"<<endl;
543       cout<<"Number of REC POINTS for module "<<mod<<": "<<nrecp<<endl;
544     }
545     for(Int_t irec=0;irec<nrecp;irec++) {
546       AliITSRecPoint *recp = (AliITSRecPoint*)ITSrec->At(irec);
547       lc[0]=recp->GetDetLocalX();
548       lc[2]=recp->GetDetLocalZ();
549       geom->LtoG(mod,lc,gc);
550       if(verb){
551         cout<<"recp # "<<irec<<" local coordinates. lx= "<<lc[0]<<" lz= ";
552         cout<<lc[2]<<endl;
553         Float_t r=sqrt(gc[0]*gc[0]+gc[1]*gc[1]);
554         cout<<"Global coor. + radius "<<gc[0]<<" "<<gc[1]<<" "<<gc[2]<<" ";
555         cout<<r<<endl;
556         cout<<"Associated track "<<recp->GetLabel(0)<<endl;
557       }
558       h2->Fill(gc[0],gc[1]);
559       h1->Fill(gc[2]);
560     }
561   }
562   return nrecp;
563 }
564
565 void GetDigits(TObject *tmps,TObject *ge,TClonesArray *ITSdigits, Int_t subd, Int_t mod, Bool_t verbose, TObjArray & histos){
566   AliITSsegmentation *seg = (AliITSsegmentation*)tmps;
567   AliITSgeom *geom = (AliITSgeom*)ge;
568   Int_t layer;
569   Int_t ladder;
570   Int_t detec;
571   if(subd==2){
572     geom->GetModuleId(mod,layer,ladder,detec);
573     seg->SetLayer(layer);
574   }
575   Float_t lcoor[3]; for(Int_t j=0; j<3; j++) lcoor[j]=0.;  //local coord dig.
576   Float_t gcoor[3]; for(Int_t j=0; j<3; j++) gcoor[j]=0.; // global coo. dig.
577   Float_t ragdig; // Radius digit
578   TArrayI ssdone(5000);  // used to match p and n side digits of strips
579   TArrayI pair(5000);    // as above 
580   Int_t ndigits = ITSdigits->GetEntries();
581   AliITSdigit *digs;
582   if(ndigits){ 
583     if(verbose){
584       cout<<"-------------------------------------------------------"<<endl;
585       cout<<"Number of DIGITS for module "<<mod<<": "<<ndigits<<endl;
586     }
587     // Get the coordinates of the module
588     if(subd==2){
589       for (Int_t digit=0;digit<ndigits;digit++){
590         ssdone[digit]=0;
591         pair[digit]=0;
592       }
593     }
594     for (Int_t digit=0;digit<ndigits;digit++) {
595       digs = (AliITSdigit*)ITSdigits->UncheckedAt(digit);
596       Int_t iz=digs->GetCoord1();  // cell number z
597       Int_t ix=digs->GetCoord2();  // cell number x
598       // Get local coordinates of the element 
599       if(subd<2){
600         seg->DetToLocal(ix,iz,lcoor[0],lcoor[2]);
601       }
602       else{
603         // SSD: if iz==0 ---> N side; if iz==1 P side
604         if(ssdone[digit]==0){
605           ssdone[digit]=1;
606           pair[digit]=-1;
607           Bool_t pside=(iz==1);
608           Bool_t impaired=kTRUE;
609           Int_t pstrip=0;
610           Int_t nstrip=0;
611           if(pside)pstrip=ix;
612           if(!pside)nstrip=ix;
613           for(Int_t digi2=0;digi2<ndigits;digi2++){
614             if(ssdone[digi2]==0 && impaired){
615               AliITSdigitSSD *dig2=(AliITSdigitSSD*)ITSdigits->UncheckedAt(digi2);
616               if(dig2->GetCoord1() != iz && dig2->GetTrack(0)==digs->GetTrack(0) && dig2->GetTrack(0)>0){
617                 ssdone[digi2]=2;
618                 pair[digit]=digi2;
619                 if(pside)nstrip=dig2->GetCoord2();
620                 if(!pside)pstrip=dig2->GetCoord2();
621                 impaired=kFALSE;
622               }
623             }
624           }
625           if(!impaired)seg->GetPadCxz(pstrip,nstrip,lcoor[0],lcoor[2]);
626         }
627       }
628       if(subd<2 || (subd==2 && ssdone[digit]==1)){
629         Int_t coor1=digs->GetCoord1();
630         Int_t coor2=digs->GetCoord2();
631         Int_t tra0=digs->GetTrack(0);
632         if(verbose){
633           cout<<"digit # "<<digit<<" fCoord1= "<<coor1<<" fCoord2= ";
634           cout<<coor2<<" track "<<tra0<<" "<<endl;
635           if(subd<2)cout<<"local coordinates -- x="<<lcoor[0]<<", z=";
636           cout<<lcoor[2]<<endl;
637           if(subd==2){
638             if(pair[digit]==-1){
639               cout<<"(digit  not paired)"<<endl;
640             }
641             else {
642               cout<<"local coordinates -- x="<<lcoor[0]<<", z="<<lcoor[2];
643               cout<<endl;
644               Int_t dtmp=pair[digit];
645               AliITSdigitSSD *dig2=(AliITSdigitSSD*)ITSdigits->UncheckedAt(dtmp);
646               Int_t coor1b=dig2->GetCoord1();
647               Int_t coor2b=dig2->GetCoord2();
648               Int_t tra0b=dig2->GetTrack(0);
649               cout<<"(digit paired with digit #"<<dtmp<<endl;
650               cout<<"with fCoord1= "<<coor1b<<" fCoord2= "<<coor2b;
651               cout<<" track "<<tra0b<<")"<<endl;
652             }
653           }
654         }
655         if(subd<2 || (subd==2 && pair[digit]!=-1)){
656           // Global coordinates of the element
657           //SDD and SPD use cm, SSD microns (GetPadCxz)
658           if(subd==2)for(Int_t j=0;j<3;j++)lcoor[j]=lcoor[j]/10000.;
659           lcoor[1]=0.;
660           geom->LtoG(mod,lcoor,gcoor);  // global coord. in cm
661           ragdig=sqrt(gcoor[0]*gcoor[0]+gcoor[1]*gcoor[1]);
662           if(verbose){
663             cout<<"global coordinates "<<gcoor[0]<<" "<<gcoor[1];
664             cout<<" "<<gcoor[2]<<" Radius "<<ragdig<<endl;
665           }   
666           //Fill histograms
667           TH2F *bidi = (TH2F*)histos.At(subd*9);
668           TH1F *uni = (TH1F*)histos.At(3+subd*9);
669           bidi->Fill(gcoor[0],gcoor[1]);
670           uni->Fill(gcoor[2]);
671         }
672       }
673     } // loop on digits for this module
674   } // if(ndigits>0....
675 }
676