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