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