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