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