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