]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITS.cxx
Release version of ITS code
[u/mrichter/AliRoot.git] / ITS / AliITS.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17 $Log$
18 Revision 1.9.2.3  2000/02/02 13:42:09  barbera
19 fixed AliITS.cxx for new AliRun structure. Added ITS hits list to list of hits which will have their track numbers updated
20
21 Revision 1.9.2.2  2000/01/23 03:03:13  nilsen
22 //fixed FillModule. Removed fi(fabs(xl)<dx....
23
24 Revision 1.9.2.1  2000/01/12 19:03:32  nilsen
25 This is the version of the files after the merging done in December 1999.
26 See the ReadMe110100.txt file for details
27
28 Revision 1.9  1999/11/14 14:33:25  fca
29 Correct problems with distructors and pointers, thanks to I.Hrivnacova
30
31 Revision 1.8  1999/09/29 09:24:19  fca
32 Introduction of the Copyright and cvs Log
33
34 */
35
36 ///////////////////////////////////////////////////////////////////////////////
37 //
38 //      An overview of the basic philosophy of the ITS code development
39 // and analysis is show in the figure below.
40 //Begin_Html
41 /*
42 <img src="picts/ITS/ITS_Analysis_schema.gif">
43 </pre>
44 <br clear=left>
45 <font size=+2 color=red>
46 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
47 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
48 </font>
49 <pre>
50 */
51 //End_Html
52 //
53 //  AliITS. Inner Traking System base class.
54 //  This class contains the base procedures for the Inner Tracking System
55 //
56 //Begin_Html
57 /*
58 <img src="picts/ITS/AliITS_Class_Diagram.gif">
59 </pre>
60 <br clear=left>
61 <font size=+2 color=red>
62 <p>This show the class diagram of the different elements that are part of
63 the AliITS class.
64 </font>
65 <pre>
66 */
67 //End_Html
68 //
69 // Version: 0
70 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
71 //
72 // Version: 1
73 // Modified and documented by Bjorn S. Nilsen
74 // July 11 1999
75 //
76 // Version: 2
77 // Modified and documented by A. Bologna
78 // October 18 1999
79 //
80 // AliITS is the general base class for the ITS. Also see AliDetector for
81 // futher information.
82 //
83 ///////////////////////////////////////////////////////////////////////////////
84  
85 #include <TMath.h>
86 #include <TRandom.h>
87 #include <TVector.h>
88 #include <TObjArray.h>
89 #include <TROOT.h>
90 #include <TObjectTable.h>
91
92
93
94 #include "AliRun.h"
95 #include "AliITS.h"
96 #include "AliITSMap.h"
97 #include "AliITSDetType.h"
98 #include "AliITSClusterFinder.h"
99 #include "AliITSsimulation.h"
100 #include "AliITSsegmentationSPD.h"
101 #include "AliITSresponseSPD.h"
102 #include "AliITSsegmentationSDD.h"
103 #include "AliITSresponseSDD.h"
104 #include "AliITSsegmentationSSD.h"
105 #include "AliITSresponseSSD.h"
106
107
108
109 ClassImp(AliITS)
110  
111 //_____________________________________________________________________________
112 AliITS::AliITS() : AliDetector() {
113   //
114   // Default initialiser for ITS
115   //     The default constructor of the AliITS class. In addition to
116   // creating the AliITS class it zeros the variables fIshunt (a member
117   // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
118   // fITSpoints, fIdSens, and fIdName. The AliDetector default constructor
119   // is also called.
120   //
121
122   fIshunt     = 0;
123   fEuclidOut  = 0;
124
125   //fNDetTypes = fgkNTYPES;
126   fIdN        = 0;
127   fIdName     = 0;
128   fIdSens     = 0;
129   fITSmodules = 0;
130   //
131   fDetTypes   = 0;
132   //
133   fDtype  = 0;
134   fNdtype = 0;
135   fCtype  = 0;
136   fNctype = 0;
137   fRecPoints = 0;
138   fNRecPoints = 0;
139   fTreeC = 0;
140   //
141   fITSgeom=0;
142 }
143
144 //_____________________________________________________________________________
145 AliITS::AliITS(const char *name, const char *title):AliDetector(name,title){
146   //
147   // Default initialiser for ITS
148   //     The constructor of the AliITS class. In addition to creating the
149   // AliITS class, it allocates memory for the TClonesArrays fHits and
150   // fDigits, and for the TObjArray fITSpoints. It also zeros the variables
151   // fIshunt (a member of AliDetector class), fEuclidOut, and fIdN, and zeros
152   // the pointers fIdSens and fIdName. To help in displaying hits via the ROOT
153   // macro display.C AliITS also sets the marker color to red. The variables
154   // passes with this constructor, const char *name and *title, are used by
155   // the constructor of AliDetector class. See AliDetector class for a
156   // description of these parameters and its constructor functions.
157   //
158
159   fHits       = new TClonesArray("AliITShit", 1560);
160   gAlice->AddHitList(fHits);
161
162   //fNDetTypes = fgkNTYPES;
163
164   fNdtype = new Int_t[fgkNTYPES];
165   fDtype = new TObjArray(fgkNTYPES);
166
167   fNctype = new Int_t[fgkNTYPES];
168   fCtype = new TObjArray(fgkNTYPES);
169
170
171   fRecPoints = 0;
172   fNRecPoints = 0;
173
174   fTreeC = 0;
175
176   fITSmodules = 0; 
177
178   fIshunt     = 0;
179   fEuclidOut  = 0;
180   fIdN        = 0;
181   fIdName     = 0;
182   fIdSens     = 0;
183  
184   fDetTypes = new TObjArray(fgkNTYPES);  
185
186   Int_t i;
187   for(i=0;i<fgkNTYPES;i++) {
188     (*fDetTypes)[i]=new AliITSDetType(); 
189     fNdtype[i]=0;
190     fNctype[i]=0;
191    }
192   //
193
194   SetMarkerColor(kRed);
195
196   fITSgeom=0;
197 }
198 //___________________________________________________________________________
199 AliITS::AliITS(AliITS &source){
200   // copy constructor
201   if(this==&source) return;
202   printf("Error: You are not allowed to make a copy of the AliITS\n");
203   exit(1);
204 }
205 //____________________________________________________________________________
206 AliITS& AliITS::operator=(AliITS &source){
207   // assignment operator
208   if(this==&source) return *this;
209   printf("Error: You are not allowed to make a copy of the AliITS\n");
210   exit(1);
211   return *this; //fake return
212 }
213 //____________________________________________________________________________
214 void AliITS::ClearModules(){
215
216   //clear the modules TObjArray
217   Int_t i;
218
219   if(fITSmodules!=0) {
220         Int_t indSPD = fITSgeom->GetModuleIndex(2,fITSgeom->GetNladders(2),
221                                                 fITSgeom->GetNdetectors(2));
222         Int_t indSDD = fITSgeom->GetModuleIndex(4,fITSgeom->GetNladders(4),
223                                                 fITSgeom->GetNdetectors(4));
224       for(i=0;i<fITSmodules->GetEntriesFast();i++){
225             if(!fITSmodules->At(i)) continue;
226             if(i<indSPD)
227               delete (AliITSmodule *) fITSmodules->At(i);
228             else if(i<indSDD)
229               delete (AliITSmodule *) fITSmodules->At(i);
230             else
231               delete (AliITSmodule *) fITSmodules->At(i);
232       } // end for i
233   }// end if fITSmodules!=0
234
235 }
236 //_____________________________________________________________________________
237 AliITS::~AliITS(){
238   //
239   // Default distructor for ITS
240   //     The default destructor of the AliITS class. In addition to deleting
241   // the AliITS class it deletes the memory pointed to by the fHits, fDigits,
242   // fIdSens, fIdName, and fITSpoints.
243   //
244
245
246   delete fHits;
247   delete fDigits;
248   delete fRecPoints;
249   if(fIdName!=0) delete[] fIdName;
250   if(fIdSens!=0) delete[] fIdSens;
251   if(fITSmodules!=0) {
252       this->ClearModules();
253       delete fITSmodules;
254   }// end if fITSmodules!=0
255
256   //
257   Int_t i;
258   if(fDtype) {
259     for (i=0;i<fgkNTYPES;i++) {
260       delete (*fDtype)[i];
261       fNdtype[i]=0;
262     }
263   }
264
265   for (i=0;i<fgkNTYPES;i++) {
266       delete (*fCtype)[i];
267       fNctype[i]=0;
268   }
269
270   //
271
272   if (fDetTypes) {
273     fDetTypes->Delete();
274     delete fDetTypes;
275   }
276
277   if (fTreeC) delete fTreeC;
278
279 }
280
281 //___________________________________________
282 AliITSDetType* AliITS::DetType(Int_t id)
283 {
284   //return pointer to id detector type
285     return ((AliITSDetType*) (*fDetTypes)[id]);
286
287 }
288 //___________________________________________
289 void AliITS::SetClasses(Int_t id, const char *digit, const char *cluster)
290 {
291   //set the digit and cluster classes to be used for the id detector type
292     ((AliITSDetType*) (*fDetTypes)[id])->ClassNames(digit,cluster);
293
294 }
295 //___________________________________________
296 void AliITS::SetResponseModel(Int_t id, AliITSresponse *response)
297 {
298   //set the response model for the id detector type
299
300     ((AliITSDetType*) (*fDetTypes)[id])->ResponseModel(response);
301
302 }
303
304 //___________________________________________
305 void AliITS::SetSegmentationModel(Int_t id, AliITSsegmentation *seg)
306 {
307   //set the segmentation model for the id detector type
308
309     ((AliITSDetType*) (*fDetTypes)[id])->SegmentationModel(seg);
310
311 }
312
313 //___________________________________________
314 void AliITS::SetSimulationModel(Int_t id, AliITSsimulation *sim)
315 {
316   //set the simulation model for the id detector type
317
318    ((AliITSDetType*) (*fDetTypes)[id])->SimulationModel(sim);
319
320 }
321 //___________________________________________
322 void AliITS::SetReconstructionModel(Int_t id, AliITSClusterFinder *reconst)
323 {
324   //set the cluster finder model for the id detector type
325
326    ((AliITSDetType*) (*fDetTypes)[id])->ReconstructionModel(reconst);
327
328 }
329
330 //_____________________________________________________________________________
331 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
332   //
333   // Add an ITS hit
334   //     The function to add information to the AliITShit class. See the
335   // AliITShit class for a full description. This function allocates the
336   // necessary new space for the hit information and passes the variable
337   // track, and the pointers *vol and *hits to the AliITShit constructor
338   // function.
339   //
340   TClonesArray &lhits = *fHits;
341   new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
342 }
343 //_____________________________________________________________________________
344 void AliITS::AddRealDigit(Int_t id, Int_t *digits) 
345 {
346   // add a real digit - as coming from data
347
348   TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
349   new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
350
351 }
352 //_____________________________________________________________________________
353 void AliITS::AddSimDigit(Int_t id, AliITSdigit *d) 
354 {
355
356   // add a simulated digit
357
358   TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
359
360   switch(id)
361   {
362   case 0:
363      new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
364      break;
365   case 1:
366      new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
367      break;
368   case 2:
369      new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
370      break;
371   }
372
373 }
374
375 //_____________________________________________________________________________
376 void AliITS::AddSimDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,Int_t *hits,Float_t *charges){
377
378   // add a simulated digit to the list
379
380   TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
381   switch(id)
382   {
383   case 0:
384      new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks,hits);
385      break;
386   case 1:
387      new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,hits,charges);
388      break;
389   case 2:
390      new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks,hits);
391      break;
392   }
393  
394 }
395
396 //_____________________________________________________________________________
397 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c) 
398 {
399
400   // add a cluster to the list
401
402   TClonesArray &lcl = *((TClonesArray*)(*fCtype)[id]);
403
404   switch(id)
405   {
406   case 0:
407      new(lcl[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
408      break;
409   case 1:
410      new(lcl[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
411      break;
412   case 2:
413      new(lcl[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
414      break;
415   }
416
417 }
418
419
420 //_____________________________________________________________________________
421 void AliITS::AddRecPoint(const AliITSRecPoint &r)
422 {
423   //
424   // Add a reconstructed space point to the list
425   //
426   TClonesArray &lrecp = *fRecPoints;
427   new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
428 }
429  
430
431 //____________________________________________
432 void AliITS::ResetDigits()
433 {
434     //
435     // Reset number of digits and the digits array for thE ITS detector
436     //
437
438     if (!fDtype) return;
439
440     Int_t i;
441     for (i=0;i<fgkNTYPES;i++ ) {
442         if ((*fDtype)[i])    ((TClonesArray*)(*fDtype)[i])->Clear();
443         if (fNdtype)  fNdtype[i]=0;
444     }
445 }
446
447 //____________________________________________
448 void AliITS::ResetDigits(Int_t i)
449 {
450     //
451     // Reset number of digits and the digits array for this branch
452     //
453   if ((*fDtype)[i])    ((TClonesArray*)(*fDtype)[i])->Clear();
454   if (fNdtype)  fNdtype[i]=0;
455 }
456
457
458 //____________________________________________
459 void AliITS::ResetClusters()
460 {
461     //
462     // Reset number of clusters and the clusters array for ITS
463     //
464
465     Int_t i;
466     for (i=0;i<fgkNTYPES;i++ ) {
467         if ((*fCtype)[i])    ((TClonesArray*)(*fCtype)[i])->Clear();
468         if (fNctype)  fNctype[i]=0;
469     }
470
471 }
472
473 //____________________________________________
474 void AliITS::ResetClusters(Int_t i)
475 {
476     //
477     // Reset number of clusters and the clusters array for this branch
478     //
479         if ((*fCtype)[i])    ((TClonesArray*)(*fCtype)[i])->Clear();
480         if (fNctype)  fNctype[i]=0;
481
482 }
483
484
485 //____________________________________________
486 void AliITS::ResetRecPoints()
487 {
488     //
489     // Reset number of rec points and the rec points array 
490     //
491     fNRecPoints = 0;
492     if (fRecPoints) fRecPoints->Clear();
493
494 }
495
496 //_____________________________________________________________________________
497 Int_t AliITS::DistancetoPrimitive(Int_t , Int_t ){
498   //
499   // Distance from mouse to ITS on the screen. Dummy routine
500   //     A dummy routine used by the ROOT macro display.C to allow for the
501   // use of the mouse (pointing device) in the macro. In general this should
502   // never be called. If it is it returns the number 9999 for any value of
503   // x and y.
504   //
505   return 9999;
506 }
507
508 //_____________________________________________________________________________
509 void AliITS::Init(){
510   //
511   // Initialise ITS after it has been built
512   //     This routine initializes the AliITS class. It is intended to be called
513   // from the Init function in AliITSv?. Besides displaying a banner
514   // indicating that it has been called it initializes the array fIdSens
515   // and sets the default segmentation, response, digit and raw cluster classes
516   // Therefore it should be called after a call to CreateGeometry.
517   //
518
519
520   SetDefaults();
521
522   Int_t i;
523   printf("\n");
524   for(i=0;i<35;i++) printf("*");
525   printf(" ITS_INIT ");
526   for(i=0;i<35;i++) printf("*");
527   printf("\n");
528   //
529   //
530   for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
531   //
532   for(i=0;i<80;i++) printf("*");
533   printf("\n");
534 }
535
536 //_____________________________________________________________________________
537 void AliITS::SetDefaults()
538 {
539   // sets the default segmentation, response, digit and raw cluster classes
540
541
542   AliITSDetType *iDetType;
543
544   //SPD 
545
546   AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(fITSgeom);
547   AliITSresponseSPD *resp0=new AliITSresponseSPD();
548   iDetType=DetType(0); 
549   if (!iDetType->GetSegmentationModel()) SetSegmentationModel(0,seg0); 
550   if (!iDetType->GetResponseModel()) SetResponseModel(0,resp0); 
551   // set digit and raw cluster classes to be used
552   const char *kData0=resp0->DataType();
553   if (strstr(kData0,"real")) {
554       iDetType->ClassNames("AliITSdigit","AliITSRawClusterSPD");
555   } else iDetType->ClassNames("AliITSdigitSPD","AliITSRawClusterSPD");
556
557   // SDD                                          //
558   AliITSresponseSDD *resp1=new AliITSresponseSDD();
559   AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(fITSgeom,resp1);
560   iDetType=DetType(1); 
561   if (!iDetType->GetSegmentationModel()) SetSegmentationModel(1,seg1); 
562   if (!iDetType->GetResponseModel()) SetResponseModel(1,resp1); 
563   const char *kData1=resp1->DataType();
564   const char *kopt=resp1->ZeroSuppOption();
565   if ((!strstr(kopt,"2D")) && (!strstr(kopt,"1D")) || strstr(kData1,"real") ) {
566       iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
567   } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
568
569   // SSD
570   AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(fITSgeom);
571   AliITSresponseSSD *resp2=new AliITSresponseSSD();
572   iDetType=DetType(2); 
573   if (!iDetType->GetSegmentationModel()) SetSegmentationModel(2,seg2); 
574   if (!iDetType->GetResponseModel()) SetResponseModel(2,resp2); 
575   const char *kData2=resp2->DataType();
576   if (strstr(kData2,"real")) {
577       iDetType->ClassNames("AliITSdigit","AliITSRawClusterSSD");
578   } else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
579
580   if (fgkNTYPES>3) {
581     Warning("SetDefaults","Only the three basic detector types are initialised!");
582   } 
583
584 }
585
586
587 //_____________________________________________________________________________
588 void AliITS::SetDefaultSimulation()
589 {
590   // to be written
591
592 }
593 //_____________________________________________________________________________
594 void AliITS::SetDefaultClusterFinders()
595 {
596   // to be written
597
598 }
599 //_____________________________________________________________________________
600
601 void AliITS::MakeTreeC(Option_t *option)
602 {
603   // create a separate tree to store the clusters
604
605      char *optC = strstr(option,"C");
606      if (optC && !fTreeC) fTreeC = new TTree("TC","Clusters in ITS");
607
608      Int_t buffersize = 4000;
609      char branchname[30];
610
611      char *det[3] = {"SPD","SDD","SSD"};
612
613      // one branch for Clusters per type of detector
614      Int_t i;
615      for (i=0; i<fgkNTYPES ;i++) {
616         if (fgkNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
617         else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
618         if (fCtype   && fTreeC) {
619            TreeC()->Branch(branchname,&((*fCtype)[i]), buffersize);
620            printf("Making Branch %s for Clusters of detector type %d\n",branchname,i+1);
621         }       
622      }
623
624 }
625
626 //_____________________________________________________________________________
627 void AliITS::GetTreeC(Int_t event)
628 {
629
630   // get the clusters tree for this event and set the branch address
631     char treeName[20];
632     char branchname[30];
633
634     char *det[3] = {"SPD","SDD","SSD"};
635
636     ResetClusters();
637     if (fTreeC) {
638           delete fTreeC;
639     }
640
641     sprintf(treeName,"TreeC%d",event);
642     fTreeC = (TTree*)gDirectory->Get(treeName);
643
644
645     TBranch *branch;
646     if (fTreeC) {
647         Int_t i;
648         for (i=0; i<fgkNTYPES; i++) {
649            if (fgkNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
650            else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
651            if (fCtype) {
652                 branch = fTreeC->GetBranch(branchname);
653                 if (branch) branch->SetAddress(&((*fCtype)[i]));
654             }
655         }
656     } else {
657         printf("ERROR: cannot find Clusters Tree for event:%d\n",event);
658     }
659
660 }
661 //_____________________________________________________________________________
662 void AliITS::MakeBranch(Option_t* option){
663   //
664   // Creates Tree branches for the ITS.
665   //
666
667
668   Int_t buffersize = 4000;
669   char branchname[30];
670   sprintf(branchname,"%s",GetName());
671
672   AliDetector::MakeBranch(option);
673
674
675 // one branch for digits per type of detector
676   
677    char *det[3] = {"SPD","SDD","SSD"};
678
679    char digclass[40];
680    char clclass[40];
681
682    Int_t i;
683    for (i=0; i<fgkNTYPES ;i++) {
684        AliITSDetType *iDetType=DetType(i); 
685        iDetType->GetClassNames(digclass,clclass);
686        //printf("i, digclass, recclass %d %s %s\n",i,digclass,clclass); 
687        // digits
688        (*fDtype)[i] = new TClonesArray(digclass,10000); 
689        // clusters
690        (*fCtype)[i] = new TClonesArray(clclass,10000); 
691    }
692
693
694   for (i=0; i<fgkNTYPES ;i++) {
695       if (fgkNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
696       else  sprintf(branchname,"%sDigits%d",GetName(),i+1);
697       
698       if (fDtype   && gAlice->TreeD()) {
699           gAlice->TreeD()->Branch(branchname,&((*fDtype)[i]), buffersize);
700           printf("Making Branch %s for digits of type %d\n",branchname,i+1);
701       } 
702   }
703
704   // only one branch for rec points for all detector types
705   sprintf(branchname,"%sRecPoints",GetName());
706
707   fRecPoints=new TClonesArray("AliITSRecPoint",10000);
708
709   if (fRecPoints && gAlice->TreeR()) {
710     gAlice->TreeR()->Branch(branchname,&fRecPoints, buffersize);
711     printf("Making Branch %s for reconstructed space points\n",branchname);
712   }     
713
714
715 }
716
717 //___________________________________________
718 void AliITS::SetTreeAddress()
719 {
720
721   // Set branch address for the Trees.
722
723   char branchname[30];
724   AliDetector::SetTreeAddress();
725
726   char *det[3] = {"SPD","SDD","SSD"};
727
728   TBranch *branch;
729   TTree *treeD = gAlice->TreeD();
730   TTree *treeR = gAlice->TreeR();
731
732   Int_t i;
733   if (treeD) {
734       for (i=0; i<fgkNTYPES; i++) {
735           if (fgkNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
736           else  sprintf(branchname,"%sDigits%d",GetName(),i+1);
737           if (fDtype) {
738               branch = treeD->GetBranch(branchname);
739               if (branch) branch->SetAddress(&((*fDtype)[i]));
740           }
741       }
742   }
743
744  
745   if (treeR) {
746     sprintf(branchname,"%sRecPoints",GetName());
747     if (fRecPoints) {
748       branch = treeR->GetBranch(branchname);
749       if (branch) branch->SetAddress(&fRecPoints);
750     }
751   }
752   
753
754 }
755
756 //____________________________________________________________________________
757 void AliITS::InitModules(Int_t size,Int_t &nmodules){
758
759   //initialize the modules array
760
761   if(fITSmodules){ 
762     //this->ClearModules();
763       delete fITSmodules;
764   }
765
766     Int_t nl,indexMAX,index;
767     Int_t indSPD,indSDD;
768
769     if(size<=0){ // default to using data stored in AliITSgeom
770         if(fITSgeom==0) {
771             printf("Error in AliITS::InitModule fITSgeom not defined\n");
772             return;
773         } // end if fITSgeom==0
774         nl = fITSgeom->GetNlayers();
775         indexMAX = fITSgeom->GetModuleIndex(nl,fITSgeom->GetNladders(nl),
776                                             fITSgeom->GetNdetectors(nl))+1;
777         nmodules = indexMAX;
778         fITSmodules = new TObjArray(indexMAX);
779         indSPD = fITSgeom->GetModuleIndex(2,fITSgeom->GetNladders(2),
780                                             fITSgeom->GetNdetectors(2));
781         indSDD = fITSgeom->GetModuleIndex(4,fITSgeom->GetNladders(4),
782                                             fITSgeom->GetNdetectors(4));
783         for(index=0;index<indexMAX;index++){
784             if(index<=indSPD)
785                 fITSmodules->AddAt( new AliITSmodule(index),index);
786             else if(index<=indSDD)
787                 fITSmodules->AddAt( new AliITSmodule(index),index);
788             else
789                 fITSmodules->AddAt( new AliITSmodule(index),index);
790         } // end for index
791     }else{
792         fITSmodules = new TObjArray(size);
793         nmodules = size;
794     } // end i size<=0
795 }
796
797 //____________________________________________________________________________
798 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t nmodules,Option_t *option,Text_t *filename){
799
800   // fill the modules with the sorted by module hits; add hits from background
801   // if option=Add
802
803     static TTree *trH1;                 //Tree with background hits
804     static TClonesArray *fHits2;        //List of hits for one track only
805
806     static Bool_t first=kTRUE;
807     static TFile *file;
808     char *addBgr = strstr(option,"Add");
809
810
811     if (addBgr ) {
812         if(first) {
813             cout<<"filename "<<filename<<endl;
814             file=new TFile(filename);
815             cout<<"I have opened "<<filename<<" file "<<endl;
816             fHits2     = new TClonesArray("AliITShit",1000  );
817         }           
818         first=kFALSE;
819         file->cd();
820         file->ls();
821         // Get Hits Tree header from file
822         if(fHits2) fHits2->Clear();
823         if(trH1) delete trH1;
824         trH1=0;
825         
826         char treeName[20];
827         sprintf(treeName,"TreeH%d",bgrev);
828         trH1 = (TTree*)gDirectory->Get(treeName);
829         //printf("TrH1 %p of treename %s for event %d \n",trH1,treeName,bgrev);
830         
831         if (!trH1) {
832             printf("ERROR: cannot find Hits Tree for event:%d\n",bgrev);
833         }
834         // Set branch addresses
835         TBranch *branch;
836         char branchname[20];
837         sprintf(branchname,"%s",GetName());
838         if (trH1 && fHits2) {
839             branch = trH1->GetBranch(branchname);
840             if (branch) branch->SetAddress(&fHits2);
841         }
842
843         // test
844         //Int_t ntracks1 =(Int_t)TrH1->GetEntries();
845         //printf("background - ntracks1 - %d\n",ntracks1);
846    }
847
848     Int_t npart = gAlice->GetEvent(evnt);
849     if(npart<=0) return;
850     TClonesArray *itsHits = this->Hits();
851     Int_t lay,lad,det,index;
852     AliITShit *itsHit=0;
853     AliITSmodule *mod=0;
854
855     TTree *iTH = gAlice->TreeH();
856     Int_t ntracks =(Int_t) iTH->GetEntries();
857
858     Int_t t,h;
859     for(t=0; t<ntracks; t++){
860         gAlice->ResetHits();
861         iTH->GetEvent(t);
862         Int_t nhits = itsHits->GetEntriesFast();
863         if (!nhits) continue;
864         for(h=0; h<nhits; h++){
865             itsHit = (AliITShit *)itsHits->UncheckedAt(h);
866             itsHit->GetDetectorID(lay,lad,det);
867             index = fITSgeom->GetModuleIndex(lay,lad,det);
868             mod = this->GetModule(index);
869             if(lay == 1 || lay == 2)
870                 mod->AddHit((AliITShit *) itsHit,t,h);
871             else if(lay == 3 || lay == 4)
872                     mod->AddHit((AliITShit *) itsHit,t,h);
873             else if(lay == 5 || lay ==6)
874                 mod->AddHit((AliITShit *)itsHit,t,h);
875             else
876                 mod->AddHit(itsHit,t,h);
877
878         } // end loop over hits 
879     } // end loop over tracks
880
881     // open the file with background
882     
883     if (addBgr ) {
884           Int_t track,i;
885           ntracks =(Int_t)trH1->GetEntries();
886             //printf("background - ntracks1 %d\n",ntracks);
887             //printf("background - Start loop over tracks \n");     
888             //   Loop over tracks
889
890             for (track=0; track<ntracks; track++) {
891
892                 if (fHits2)       fHits2->Clear();
893                 trH1->GetEvent(track);
894                 //   Loop over hits
895                 for(i=0;i<fHits2->GetEntriesFast();++i) {
896
897                     itsHit=(AliITShit*) (*fHits2)[i];
898                     itsHit->GetDetectorID(lay,lad,det);
899                     index = fITSgeom->GetModuleIndex(lay,lad,det);
900                     mod = this->GetModule(index);
901
902                     if(lay == 1 || lay == 2)
903                        mod->AddHit((AliITShit *) itsHit,track,i);
904                     else if(lay == 3 || lay == 4)
905                             mod->AddHit((AliITShit *) itsHit,track,i);
906                     else if(lay == 5 || lay ==6)
907                             mod->AddHit((AliITShit *)itsHit,track,i);
908                     else
909                        mod->AddHit(itsHit,track,i);
910
911
912                }  // end loop over hits
913             } // end loop over tracks
914
915             TTree *fAli=gAlice->TreeK();
916             TFile *fileAli=0;
917             
918             if (fAli) fileAli =fAli->GetCurrentFile();
919             file->cd();
920
921     } // end if add
922
923     //gObjectTable->Print();
924
925 }
926
927
928 //____________________________________________________________________________
929 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size, Option_t *option, Option_t *opt,Text_t *filename)
930 {
931     // keep galice.root for signal and name differently the file for 
932     // background when add! otherwise the track info for signal will be lost !
933   
934    char *all = strstr(opt,"All");
935    char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
936
937    Int_t nmodules;
938    InitModules(size,nmodules); 
939    FillModules(evNumber,bgrev,nmodules,option,filename);
940
941    TBranch *branch;
942    AliITSsimulation* sim;
943
944    TObjArray *branches=gAlice->TreeD()->GetListOfBranches();
945    AliITSgeom *geom = GetITSgeom();
946
947    Int_t id,module;
948    for (id=0;id<fgkNTYPES;id++) {
949         if (!all && !det[id]) continue;
950         branch = (TBranch*)branches->UncheckedAt(id);
951         AliITSDetType *iDetType=DetType(id); 
952         sim = (AliITSsimulation*)iDetType->GetSimulationModel();
953         if (!sim) {
954            Error("HitsToDigits","The simulation class was not instantiated!");
955            exit(1);
956            // or SetDefaultSimulation();
957         }
958         Int_t first = geom->GetStartDet(id);
959         Int_t last = geom->GetLastDet(id);
960         printf("det type %d first, last %d %d \n",id,first,last);
961         for(module=first;module<=last;module++) {
962             AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
963             sim->DigitiseModule(mod,module,evNumber);
964             // fills all branches - wasted disk space
965             gAlice->TreeD()->Fill(); 
966             ResetDigits();
967             // try and fill only the branch 
968             //branch->Fill();
969             //ResetDigits(id);
970         } // loop over modules
971    } // loop over detector types
972
973    ClearModules();
974
975    Int_t nentries=(Int_t)gAlice->TreeD()->GetEntries();
976    printf("nentries in TreeD %d\n",nentries);
977
978    char hname[30];
979    sprintf(hname,"TreeD%d",evNumber);
980    gAlice->TreeD()->Write(hname);
981    // reset tree
982    gAlice->TreeD()->Reset();
983
984 }
985
986
987 //____________________________________________________________________________
988 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt)
989 {
990   // cluster finding and reconstruction of space points
991   
992    char *all = strstr(opt,"All");
993    char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
994
995    static Bool_t first=kTRUE;
996    if (first) {
997        MakeTreeC("C");
998        first=kFALSE;
999    }
1000    TTree *iTC=TreeC();
1001  
1002    TBranch *branch;
1003    AliITSClusterFinder* rec;
1004
1005    TObjArray *branches=gAlice->TreeR()->GetListOfBranches();
1006    AliITSgeom *geom = GetITSgeom();
1007
1008    Int_t id,module;
1009    for (id=0;id<fgkNTYPES;id++) {
1010         if (!all && !det[id]) continue;
1011         branch = (TBranch*)branches->UncheckedAt(id);
1012         AliITSDetType *iDetType=DetType(id); 
1013         rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1014         if (!rec) {
1015            Error("DigitsToRecPoints","The cluster finder class was not instantiated!");
1016            exit(1);
1017            // or SetDefaultClusterFinders();
1018         }
1019         TClonesArray *itsDigits  = this->DigitsAddress(id);
1020
1021         Int_t first = geom->GetStartDet(id);
1022         Int_t last = geom->GetLastDet(id);
1023         for(module=first;module<=last;module++) {
1024               this->ResetDigits();
1025               if (all) gAlice->TreeD()->GetEvent(lastentry+module);
1026               else gAlice->TreeD()->GetEvent(lastentry+(module-first));
1027               Int_t ndigits = itsDigits->GetEntriesFast();
1028               if (ndigits) rec->FindRawClusters();
1029               gAlice->TreeR()->Fill(); 
1030               ResetRecPoints();
1031               iTC->Fill();
1032               ResetClusters();
1033               // try and fill only the branch 
1034               //branch->Fill();
1035               //ResetRecPoints(id);
1036         } // loop over modules
1037    } // loop over detector types
1038
1039
1040    Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1041    Int_t ncentries=(Int_t)iTC->GetEntries();
1042    printf(" nentries ncentries %d %d\n", nentries, ncentries);
1043
1044    char hname[30];
1045    sprintf(hname,"TreeR%d",evNumber);
1046    gAlice->TreeR()->Write(hname);
1047    // reset tree
1048    gAlice->TreeR()->Reset();
1049
1050    sprintf(hname,"TreeC%d",evNumber);
1051    iTC->Write(hname);
1052    iTC->Reset();
1053 }
1054
1055
1056 //____________________________________________________________________________
1057 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
1058 Option_t *option,Option_t *opt,Text_t *filename)
1059 {
1060     // keep galice.root for signal and name differently the file for 
1061     // background when add! otherwise the track info for signal will be lost !
1062   
1063    char *all = strstr(opt,"All");
1064    char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1065
1066    Int_t nmodules;
1067    InitModules(size,nmodules);
1068    FillModules(evNumber,bgrev,nmodules,option,filename);
1069
1070
1071    AliITSsimulation* sim;
1072    AliITSgeom *geom = GetITSgeom();
1073
1074    TRandom *random=new TRandom[9];
1075    random[0].SetSeed(111);
1076    random[1].SetSeed(222);
1077    random[2].SetSeed(333);              
1078    random[3].SetSeed(444);
1079    random[4].SetSeed(555);
1080    random[5].SetSeed(666);              
1081    random[6].SetSeed(777);
1082    random[7].SetSeed(888);
1083    random[8].SetSeed(999);              
1084
1085
1086    Int_t id,module;
1087    for (id=0;id<fgkNTYPES;id++) {
1088         if (!all && !det[id]) continue;
1089         AliITSDetType *iDetType=DetType(id); 
1090         sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1091         if (!sim) {
1092            Error("HitsToFastPoints","The simulation class was not instantiated!");
1093            exit(1);
1094            // or SetDefaultSimulation();
1095         }
1096         Int_t first = geom->GetStartDet(id);
1097         Int_t last = geom->GetLastDet(id);
1098         for(module=first;module<=last;module++) {
1099             AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
1100             sim->CreateFastRecPoints(mod,module,random);
1101             gAlice->TreeR()->Fill(); 
1102             ResetRecPoints();
1103         } // loop over modules
1104    } // loop over detector types
1105
1106
1107    ClearModules();
1108
1109    //Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1110
1111    char hname[30];
1112    sprintf(hname,"TreeR%d",evNumber);
1113    gAlice->TreeR()->Write(hname);
1114    // reset tree
1115    gAlice->TreeR()->Reset();
1116
1117    delete [] random;
1118
1119 }
1120
1121 //____________________________________________________________________________
1122 void AliITS::Streamer(TBuffer &R__b){
1123    // Stream an object of class AliITS.
1124
1125    Int_t i;
1126
1127    if (R__b.IsReading()) {
1128       Version_t R__v = R__b.ReadVersion();
1129       if (R__v) {
1130           AliDetector::Streamer(R__b);
1131           R__b >> fIdN;
1132           R__b.ReadArray(fIdSens); 
1133           for(i=0;i<fIdN;i++) fIdName[i].Streamer(R__b);
1134           R__b >> fITSgeom;
1135           R__b >> fITSmodules;
1136           R__b >> fEuclidOut;
1137           R__b >> fMajorVersion;
1138           R__b >> fMinorVersion;
1139           R__b >> fDetTypes;
1140           R__b >> fDtype;
1141           delete []fNdtype; 
1142           fNdtype = new Int_t[fgkNTYPES];
1143           R__b.ReadFastArray(fNdtype,fgkNTYPES);
1144           R__b >> fCtype;
1145           delete []fNctype; 
1146           fNctype = new Int_t[fgkNTYPES];
1147           R__b.ReadFastArray(fNctype,fgkNTYPES);
1148           R__b >> fRecPoints;
1149           R__b >> fNRecPoints;
1150           R__b >> fTreeC;
1151       } // end if R__v
1152    } else { // writing
1153       R__b.WriteVersion(AliITS::IsA());
1154       AliDetector::Streamer(R__b);
1155       R__b << fIdN;
1156       R__b.WriteArray(fIdSens,fIdN); 
1157       for(i=0;i<fIdN;i++) fIdName[i].Streamer(R__b);
1158       R__b << fITSgeom;
1159       R__b << fITSmodules;
1160       R__b << fEuclidOut;
1161       R__b << fMajorVersion;
1162       R__b << fMinorVersion;
1163       R__b << fDetTypes;
1164       R__b << fDtype;
1165       R__b.WriteFastArray(fNdtype,fgkNTYPES);
1166       R__b << fCtype;
1167       R__b.WriteFastArray(fNctype,fgkNTYPES);
1168       R__b << fRecPoints;
1169       R__b << fNRecPoints;
1170       R__b << fTreeC;
1171    } // end if
1172
1173 }
1174
1175