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