]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITS.cxx
Reducing dependencies (F.Carminati)
[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 /* $Id$ */
17
18
19 ///////////////////////////////////////////////////////////////////////////////
20 //
21 //      An overview of the basic philosophy of the ITS code development
22 // and analysis is show in the figure below.
23 //Begin_Html
24 /*
25 <img src="picts/ITS/ITS_Analysis_schema.gif">
26 </pre>
27 <br clear=left>
28 <font size=+2 color=red>
29 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
30 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
31 </font>
32 <pre>
33 */
34 //End_Html
35 //
36 //  AliITS. Inner Traking System base class.
37 //  This class contains the base procedures for the Inner Tracking System
38 //
39 //Begin_Html
40 /*
41 <img src="picts/ITS/AliITS_Class_Diagram.gif">
42 </pre>
43 <br clear=left>
44 <font size=+2 color=red>
45 <p>This show the class diagram of the different elements that are part of
46 the AliITS class.
47 </font>
48 <pre>
49 */
50 //End_Html
51 //
52 // Version: 0
53 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
54 //
55 // Version: 1
56 // Modified and documented by Bjorn S. Nilsen
57 // July 11 1999
58 //
59 // Version: 2
60 // Modified and documented by A. Bologna
61 // October 18 1999
62 //
63 // AliITS is the general base class for the ITS. Also see AliDetector for
64 // futher information.
65 //
66 ///////////////////////////////////////////////////////////////////////////////
67
68 #include <Riostream.h>
69 #include <stdlib.h>
70
71 #include <TBranch.h>
72 #include <TClonesArray.h>
73 #include <TFile.h>
74 #include <TMath.h>
75 #include <TObjectTable.h>
76 #include <TROOT.h>
77 #include <TRandom.h>
78 #include <TString.h>
79 #include <TTree.h>
80 #include <TVector.h>
81 #include <TVirtualMC.h>
82
83 #include "AliConfig.h"
84 #include "AliHeader.h"
85 #include "AliITS.h"
86 #include "AliITSClusterFinderSDD.h"
87 #include "AliITSClusterFinderSPD.h"
88 #include "AliITSClusterFinderSSD.h"
89 #include "AliITSDetType.h"
90 #include "AliITSLoader.h"
91 #include "AliITSRawCluster.h"
92 #include "AliITSRecPoint.h"
93 #include "AliITSdigit.h"
94 #include "AliITSgeom.h"
95 #include "AliITShit.h"
96 #include "AliITSmodule.h"
97 #include "AliITSpList.h"
98 #include "AliITSresponseSDD.h"
99 #include "AliITSresponseSPD.h"
100 #include "AliITSresponseSSD.h"
101 #include "AliITSsegmentationSDD.h"
102 #include "AliITSsegmentationSPD.h"
103 #include "AliITSsegmentationSSD.h"
104 #include "AliITSsimulationSDD.h"
105 #include "AliITSsimulationSPD.h"
106 #include "AliITSsimulationSSD.h"
107 #include "AliMC.h"
108
109 ClassImp(AliITS)
110
111 //______________________________________________________________________
112 AliITS::AliITS() : AliDetector() {
113     // Default initializer for ITS
114     //      The default constructor of the AliITS class. In addition to
115     // creating the AliITS class it zeros the variables fIshunt (a member
116     // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
117     // fITSpoints, fIdSens, and fIdName. The AliDetector default constructor
118     // is also called.
119     // Inputs:
120     //      none.
121     // Outputs:
122     //      none.
123     // Return:
124     //      Blank ITS class.
125
126     fIshunt     = 0;   // not zeroed in AliDetector.
127
128     // AliITS variables.
129     fEuclidOut  = 0;
130     fITSgeom    = 0;
131     fITSmodules = 0;
132     fOpt        = "All";
133 //    SetDetectors(); // default to fOpt="All". This variable not written out.
134
135     fIdN        = 0;
136     fIdName     = 0;
137     fIdSens     = 0;
138
139     fNDetTypes  = kNTYPES;
140     fDetTypes   = 0;
141
142     fSDigits    = 0;
143     fNSDigits   = 0;
144
145     fNdtype     = 0;
146     fDtype      = 0;
147
148     fCtype      = 0;
149     fNctype     = 0;
150
151     fRecPoints  = 0;
152     fNRecPoints = 0;
153
154     SetMarkerColor(kRed);
155 }
156 //______________________________________________________________________
157 AliITS::AliITS(const char *name, const char *title):AliDetector(name,title){
158     //     The standard Constructor for the ITS class. In addition to 
159     // creating the AliITS class, it allocates memory for the TClonesArrays 
160     // fHits, fSDigits, fDigits, fITSpoints, and the TObjArray of fCtype 
161     // (clusters). It also zeros the variables
162     // fIshunt (a member of AliDetector class), fEuclidOut, and fIdN, and zeros
163     // the pointers fIdSens and fIdName. To help in displaying hits via the
164     // ROOT macro display.C AliITS also sets the marker color to red. The
165     // variables passes with this constructor, const char *name and *title,
166     // are used by the constructor of AliDetector class. See AliDetector
167     // class for a description of these parameters and its constructor
168     // functions.
169     // Inputs:
170     //      const char *name      Detector name. Should always be "ITS"
171     //      const char *title     Detector title.
172     // Outputs:
173     //      none.
174     // Return:
175     //      An ITS class.
176
177     fIshunt     = 0;  // not zeroed in AliDetector
178     fHits       = new TClonesArray("AliITShit", 1560);//not done in AliDetector
179     gAlice->GetMCApp()->AddHitList(fHits);  // Not done in AliDetector.
180
181     fEuclidOut  = 0;
182     fITSgeom    = 0;
183     fITSmodules = 0;
184     SetDetectors(); // default to fOpt="All". This variable not written out.
185
186     fIdN        = 0;
187     fIdName     = 0;
188     fIdSens     = 0;
189
190     fNDetTypes  = kNTYPES;
191     fDetTypes   = new TObjArray(fNDetTypes);
192
193     fSDigits    = new TClonesArray("AliITSpListItem",1000);
194     fNSDigits   = 0;
195
196     fNdtype     = new Int_t[fNDetTypes];
197     fDtype      = new TObjArray(fNDetTypes);
198
199     fCtype      = new TObjArray(fNDetTypes);
200     fNctype     = new Int_t[fNDetTypes];
201
202     fRecPoints  = new TClonesArray("AliITSRecPoint",1000);
203     fNRecPoints = 0;
204
205     Int_t i;
206     for(i=0;i<fNDetTypes;i++) {
207         fDetTypes->AddAt(new AliITSDetType(),i); 
208         fNdtype[i] = 0;
209         fNctype[i] = 0;
210     } // end for i
211
212     SetMarkerColor(kRed);
213
214 }
215 //______________________________________________________________________
216 AliITS::~AliITS(){
217     // Default destructor for ITS.
218     //     The default destructor of the AliITS class. In addition to deleting
219     // the AliITS class it deletes the memory pointed to by the fHits, fDigits,
220     // fSDigits, fCtype, fITSmodules, fITSgeom, fRecPoints, fIdSens, fIdName, 
221     // fITSpoints, fDetType and it's contents.
222     // Inputs:
223     //      none.
224     // Outputs:
225     //      none.
226     // Return:
227     //      none.
228
229     if (fHits) {
230       fHits->Delete();
231       delete fHits;
232       fHits=0;
233     }
234     if (fSDigits) {
235       fSDigits->Delete();
236       delete fSDigits;
237       fSDigits=0;
238     }
239     if (fDigits) {
240       fDigits->Delete();
241       delete fDigits;
242       fDigits=0;
243     }
244     if (fRecPoints) {
245       fRecPoints->Delete();
246       delete fRecPoints;
247       fRecPoints=0;
248     }
249     delete[] fIdName;  // Array of TStrings
250     delete[] fIdSens;
251     if(fITSmodules) {
252         this->ClearModules();
253         delete fITSmodules;
254         fITSmodules = 0;
255     }// end if fITSmodules!=0
256
257     if(fDtype) {
258         fDtype->Delete();
259         delete fDtype;
260     } // end if fDtype
261     delete [] fNdtype;
262     if (fCtype) {
263         fCtype->Delete();
264         delete fCtype;
265         fCtype = 0;
266     } // end if fCtype
267     delete [] fNctype;
268
269     if (fDetTypes) {
270         fDetTypes->Delete();
271         delete fDetTypes;
272         fDetTypes = 0;
273     } // end if fDetTypes
274
275
276     if (fITSgeom) delete fITSgeom;
277 }
278 //______________________________________________________________________
279 AliITS::AliITS(const AliITS &source) : AliDetector(source){
280     // Copy constructor. This is a function which is not allowed to be
281     // done to the ITS. It exits with an error.
282     // Inputs:
283     //      AliITS &source  An AliITS class.
284     // Outputs:
285     //      none.
286     // Return:
287     //      none.
288
289     if(this==&source) return;
290     Error("Copy constructor",
291           "You are not allowed to make a copy of the AliITS");
292     exit(1);
293 }
294 //______________________________________________________________________
295 AliITS& AliITS::operator=(AliITS &source){
296     // Assignment operator. This is a function which is not allowed to be
297     // done to the ITS. It exits with an error.
298     // Inputs:
299     //      AliITS &source  An AliITS class.
300     // Outputs:
301     //      none.
302     // Return:
303     //      none.
304
305     if(this==&source) return *this;
306     Error("operator=","You are not allowed to make a copy of the AliITS");
307     exit(1);
308     return *this; //fake return
309 }
310 //______________________________________________________________________
311 Int_t AliITS::DistancetoPrimitive(Int_t,Int_t) const{
312     // Distance from mouse to ITS on the screen. Dummy routine
313     //     A dummy routine used by the ROOT macro display.C to allow for the
314     // use of the mouse (pointing device) in the macro. In general this should
315     // never be called. If it is it returns the number 9999 for any value of
316     // x and y.
317     // Inputs:
318     //      Int_t     Dummy screen coordinate.
319     //      Int_t     Dummy screen coordinate.
320     // Outputs:
321     //      none.
322     // Return:
323     //      Int_t     Dummy = 9999 distance to ITS.
324
325     return 9999;
326 }
327 //______________________________________________________________________
328 void AliITS::Init(){
329     // Initializer ITS after it has been built
330     //     This routine initializes the AliITS class. It is intended to be
331     // called from the Init function in AliITSv?. Besides displaying a banner
332     // indicating that it has been called it initializes the array fIdSens
333     // and sets the default segmentation, response, digit and raw cluster
334     // classes therefore it should be called after a call to CreateGeometry.
335     // Inputs:
336     //      none.
337     // Outputs:
338     //      none.
339     // Return:
340     //      none.
341     Int_t i;
342
343     SetDefaults();
344     // Array of TStrings
345     for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
346 }
347 //______________________________________________________________________
348 void AliITS::SetDefaults(){
349     // sets the default segmentation, response, digit and raw cluster classes.
350     // Inputs:
351     //      none.
352     // Outputs:
353     //      none.
354     // Return:
355     //      none.
356
357     if(fDebug) printf("%s: SetDefaults\n",ClassName());
358
359     AliITSDetType *iDetType;
360
361     //SPD
362     iDetType=DetType(0); 
363     if (!iDetType->GetSegmentationModel()) {
364         AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(fITSgeom);
365         SetSegmentationModel(0,seg0); 
366     } // end if
367     if (!iDetType->GetResponseModel()) {
368         SetResponseModel(0,new AliITSresponseSPD()); 
369     } // end if
370     // set digit and raw cluster classes to be used
371
372     const char *kData0=(iDetType->GetResponseModel())->DataType();
373     if (strstr(kData0,"real")) {
374         iDetType->ClassNames("AliITSdigit","AliITSRawClusterSPD");
375     } else iDetType->ClassNames("AliITSdigitSPD","AliITSRawClusterSPD");
376
377     // SDD
378     iDetType=DetType(1); 
379     if (!iDetType->GetResponseModel()) {
380         SetResponseModel(1,new AliITSresponseSDD("simulated")); 
381     } // end if
382     AliITSresponse *resp1=iDetType->GetResponseModel();
383     if (!iDetType->GetSegmentationModel()) {
384         AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(fITSgeom,resp1);
385         SetSegmentationModel(1,seg1); 
386     } // end if
387     const char *kData1=(iDetType->GetResponseModel())->DataType();
388     const char *kopt=iDetType->GetResponseModel()->ZeroSuppOption();
389     if((!strstr(kopt,"2D"))&&(!strstr(kopt,"1D")) || strstr(kData1,"real") ){
390         iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
391     } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
392
393     // SSD
394     iDetType=DetType(2); 
395     if (!iDetType->GetSegmentationModel()) {
396         AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(fITSgeom);
397         SetSegmentationModel(2,seg2); 
398     } // end if
399     if (!iDetType->GetResponseModel()) {
400         SetResponseModel(2,new AliITSresponseSSD("simulated"));
401     } // end if
402     const char *kData2=(iDetType->GetResponseModel())->DataType();
403     if (strstr(kData2,"real")) {
404         iDetType->ClassNames("AliITSdigit","AliITSRawClusterSSD");
405     } else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
406
407     if (kNTYPES>3) {
408         Warning("SetDefaults",
409                 "Only the three basic detector types are initialized!");
410     }  // end if
411 }
412 //______________________________________________________________________
413 void AliITS::SetDefaultSimulation(){
414     // sets the default simulation.
415     // Inputs:
416     //      none.
417     // Outputs:
418     //      none.
419     // Return:
420     //      none.
421
422     AliITSDetType *iDetType;
423     AliITSsimulation *sim;
424     iDetType=DetType(0);
425     sim = iDetType->GetSimulationModel();
426     if (!sim) {
427         AliITSsegmentation *seg0=
428             (AliITSsegmentation*)iDetType->GetSegmentationModel();
429         AliITSresponse *res0 = (AliITSresponse*)iDetType->GetResponseModel();
430         AliITSsimulationSPD *sim0=new AliITSsimulationSPD(seg0,res0);
431         SetSimulationModel(0,sim0);
432     }else{ // simulation exists, make sure it is set up properly.
433         ((AliITSsimulationSPD*)sim)->Init(
434             (AliITSsegmentationSPD*) iDetType->GetSegmentationModel(),
435             (AliITSresponseSPD*) iDetType->GetResponseModel());
436 //        if(sim->GetResponseModel()==0) sim->SetResponseModel(
437 //            (AliITSresponse*)iDetType->GetResponseModel());
438 //        if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
439 //            (AliITSsegmentation*)iDetType->GetSegmentationModel());
440     } // end if
441     iDetType=DetType(1);
442     sim = iDetType->GetSimulationModel();
443     if (!sim) {
444         AliITSsegmentation *seg1=
445             (AliITSsegmentation*)iDetType->GetSegmentationModel();
446         AliITSresponse *res1 = (AliITSresponse*)iDetType->GetResponseModel();
447         AliITSsimulationSDD *sim1=new AliITSsimulationSDD(seg1,res1);
448         SetSimulationModel(1,sim1);
449     }else{ // simulation exists, make sure it is set up properly.
450         ((AliITSsimulationSDD*)sim)->Init(
451             (AliITSsegmentationSDD*) iDetType->GetSegmentationModel(),
452             (AliITSresponseSDD*) iDetType->GetResponseModel());
453 //        if(sim->GetResponseModel()==0) sim->SetResponseModel(
454 //            (AliITSresponse*)iDetType->GetResponseModel());
455 //        if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
456 //            (AliITSsegmentation*)iDetType->GetSegmentationModel());
457     } //end if
458     iDetType=DetType(2);
459     sim = iDetType->GetSimulationModel();
460     if (!sim) {
461         AliITSsegmentation *seg2=
462             (AliITSsegmentation*)iDetType->GetSegmentationModel();
463         AliITSresponse *res2 = (AliITSresponse*)iDetType->GetResponseModel();
464         AliITSsimulationSSD *sim2=new AliITSsimulationSSD(seg2,res2);
465         SetSimulationModel(2,sim2);
466     }else{ // simulation exists, make sure it is set up properly.
467         ((AliITSsimulationSSD*)sim)->Init(
468             (AliITSsegmentationSSD*) iDetType->GetSegmentationModel(),
469             (AliITSresponseSSD*) iDetType->GetResponseModel());
470 //        if(sim->GetResponseModel()==0) sim->SetResponseModel(
471 //            (AliITSresponse*)iDetType->GetResponseModel());
472 //        if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
473 //            (AliITSsegmentation*)iDetType->GetSegmentationModel());
474     } // end if
475 }
476 //______________________________________________________________________
477 void AliITS::SetDefaultClusterFinders(){
478     // Sets the default cluster finders. Used in finding RecPoints.
479     // Inputs:
480     //      none.
481     // Outputs:
482     //      none.
483     // Return:
484     //      none.
485
486     MakeTreeC();
487     AliITSDetType *iDetType;
488
489     // SPD
490     iDetType=DetType(0);
491     if (!iDetType->GetReconstructionModel()) {
492         AliITSsegmentation *seg0 =
493             (AliITSsegmentation*)iDetType->GetSegmentationModel();
494         TClonesArray *dig0=DigitsAddress(0);
495         TClonesArray *recp0=ClustersAddress(0);
496         AliITSClusterFinderSPD *rec0 = new AliITSClusterFinderSPD(seg0,dig0,
497                                                                   recp0);
498         SetReconstructionModel(0,rec0);
499     } // end if
500
501     // SDD
502     iDetType=DetType(1);
503     if (!iDetType->GetReconstructionModel()) {
504         AliITSsegmentation *seg1 =
505             (AliITSsegmentation*)iDetType->GetSegmentationModel();
506         AliITSresponse *res1 = (AliITSresponse*)iDetType->GetResponseModel();
507         TClonesArray *dig1=DigitsAddress(1);
508         TClonesArray *recp1=ClustersAddress(1);
509         AliITSClusterFinderSDD *rec1 =
510             new AliITSClusterFinderSDD(seg1,res1,dig1,recp1);
511       SetReconstructionModel(1,rec1);
512     } // end if
513
514     // SSD
515     iDetType=DetType(2);
516     if (!iDetType->GetReconstructionModel()) {
517         AliITSsegmentation *seg2=
518             (AliITSsegmentation*)iDetType->GetSegmentationModel();
519         TClonesArray *dig2=DigitsAddress(2);
520         AliITSClusterFinderSSD *rec2= new AliITSClusterFinderSSD(seg2,dig2);
521         SetReconstructionModel(2,rec2);
522     } // end if
523 }
524 //______________________________________________________________________
525 void AliITS::MakeBranch(Option_t* option){
526     // Creates Tree branches for the ITS.
527     // Inputs:
528     //      Option_t *option    String of Tree types S,D, and/or R.
529     //      const char *file    String of the file name where these branches
530     //                          are to be stored. If blank then these branches
531     //                          are written to the same tree as the Hits were
532     //                          read from.
533     // Outputs:
534     //      none.
535     // Return:
536     //      none.
537     Bool_t cH = (strstr(option,"H")!=0);
538     Bool_t cS = (strstr(option,"S")!=0);
539     Bool_t cD = (strstr(option,"D")!=0);
540     Bool_t cR = (strstr(option,"R")!=0);
541     Bool_t cRF = (strstr(option,"RF")!=0);
542     
543     if(cRF)cR = kFALSE;
544     if(cH && (fHits == 0x0)) fHits  = new TClonesArray("AliITShit", 1560);
545     
546     AliDetector::MakeBranch(option);
547
548     if(cS) MakeBranchS(0);
549     if(cD) MakeBranchD(0);
550     if(cR) MakeBranchR(0);
551     if(cRF) MakeBranchRF(0);
552 }
553 //______________________________________________________________________
554 void AliITS::SetTreeAddress(){
555     // Set branch address for the Trees.
556     // Inputs:
557     //      none.
558     // Outputs:
559     //      none.
560     // Return:
561     //      none.
562     TTree *treeS = fLoader->TreeS();
563     TTree *treeD = fLoader->TreeD();
564     TTree *treeR = fLoader->TreeR();
565     if (fLoader->TreeH() && (fHits == 0x0)) fHits = new TClonesArray("AliITShit", 1560);
566       
567     AliDetector::SetTreeAddress();
568
569     SetTreeAddressS(treeS);
570     SetTreeAddressD(treeD);
571     SetTreeAddressR(treeR);
572 }
573 //______________________________________________________________________
574 AliITSDetType* AliITS::DetType(Int_t id){
575     // Return pointer to id detector type.
576     // Inputs:
577     //      Int_t id   detector id number.
578     // Outputs:
579     //      none.
580     // Return:
581     //      returned, a pointer to a AliITSDetType.
582
583     return ((AliITSDetType*) fDetTypes->At(id));
584 }
585 //______________________________________________________________________
586 void AliITS::SetResponseModel(Int_t id, AliITSresponse *response){
587     // Set the response model for the id detector type.
588     // Inputs:
589     //      Int_t id        detector id number.
590     //      AliITSresponse* a pointer containing an instance of AliITSresponse
591     //                      to be stored/owned b y AliITSDetType.
592     // Outputs:
593     //      none.
594     // Return:
595     //      none.
596
597     ((AliITSDetType*) fDetTypes->At(id))->ResponseModel(response);
598 }
599 //______________________________________________________________________
600 void AliITS::SetSegmentationModel(Int_t id, AliITSsegmentation *seg){
601     // Set the segmentation model for the id detector type.
602     // Inputs:
603     //      Int_t id            detector id number.
604     //      AliITSsegmentation* a pointer containing an instance of 
605     //                          AliITSsegmentation to be stored/owned b y 
606     //                          AliITSDetType.
607     // Outputs:
608     //      none.
609     // Return:
610     //      none.
611
612     ((AliITSDetType*) fDetTypes->At(id))->SegmentationModel(seg);
613 }
614 //______________________________________________________________________
615 void AliITS::SetSimulationModel(Int_t id, AliITSsimulation *sim){
616     // Set the simulation model for the id detector type.
617     // Inputs:
618     //      Int_t id        detector id number.
619     //      AliITSresponse* a pointer containing an instance of AliITSresponse
620     //                      to be stored/owned b y AliITSDetType.
621     // Outputs:
622     //      none.
623     // Return:
624     //      none.
625
626    ((AliITSDetType*) fDetTypes->At(id))->SimulationModel(sim);
627
628 }
629 //______________________________________________________________________
630 void AliITS::SetReconstructionModel(Int_t id, AliITSClusterFinder *reconst){
631     // Set the cluster finder model for the id detector type.
632     // Inputs:
633     //      Int_t id             detector id number.
634     //      AliITSClusterFinder* a pointer containing an instance of 
635     //                           AliITSClusterFinder to be stored/owned b y 
636     //                           AliITSDetType.
637     // Outputs:
638     //      none.
639     // Return:
640     //      none.
641
642     ((AliITSDetType*) fDetTypes->At(id))->ReconstructionModel(reconst);
643 }
644 //______________________________________________________________________
645 void AliITS::SetClasses(Int_t id, const char *digit, const char *cluster){
646     // Set the digit and cluster classes name to be used for the id detector
647     // type.
648     // Inputs:
649     //      Int_t id            detector id number.
650     //      const char *digit   Digit class name for detector id.
651     //      const char *cluster Cluster class name for detector id.
652     // Outputs:
653     //      none.
654     // Return:
655     //      none.
656
657     ((AliITSDetType*) fDetTypes->At(id))->ClassNames(digit,cluster);
658 }
659 //______________________________________________________________________
660 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
661     // Add an ITS hit
662     //     The function to add information to the AliITShit class. See the
663     // AliITShit class for a full description. This function allocates the
664     // necessary new space for the hit information and passes the variable
665     // track, and the pointers *vol and *hits to the AliITShit constructor
666     // function.
667     // Inputs:
668     //      Int_t   track   Track number which produced this hit.
669     //      Int_t   *vol    Array of Integer Hit information. See AliITShit.h
670     //      Float_t *hits   Array of Floating Hit information.  see AliITShit.h
671     // Outputs:
672     //      none.
673     // Return:
674     //      none.
675
676     TClonesArray &lhits = *fHits;
677     new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
678 }
679 //______________________________________________________________________
680 void AliITS::InitModules(Int_t size,Int_t &nmodules){
681     // Initialize the modules array.
682     // Inputs:
683     //      Int_t size  Size of array of the number of modules to be
684     //                  created. If size <=0 then the number of modules
685     //                  is gotten from AliITSgeom class kept in fITSgeom.
686     // Outputs:
687     //      Int_t &nmodules The number of modules existing.
688     // Return:
689     //      none.
690
691     if(fITSmodules){ 
692         fITSmodules->Delete();
693         delete fITSmodules;
694     } // end fir fITSmoudles
695
696     Int_t nl,indexMAX,index;
697
698     if(size<=0){ // default to using data stored in AliITSgeom
699         if(fITSgeom==0) {
700             Error("InitModules","fITSgeom not defined");
701             return;
702         } // end if fITSgeom==0
703         nl = fITSgeom->GetNlayers();
704         indexMAX = fITSgeom->GetModuleIndex(nl,fITSgeom->GetNladders(nl),
705                                             fITSgeom->GetNdetectors(nl))+1;
706         nmodules = indexMAX;
707         fITSmodules = new TObjArray(indexMAX);
708         for(index=0;index<indexMAX;index++){
709             fITSmodules->AddAt( new AliITSmodule(index),index);
710         } // end for index
711     }else{
712         fITSmodules = new TObjArray(size);
713         for(index=0;index<size;index++) {
714             fITSmodules->AddAt( new AliITSmodule(index),index);
715         } // end for index
716
717         nmodules = size;
718     } // end i size<=0
719 }
720 //______________________________________________________________________
721 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t nmodules,
722                          Option_t *option,Text_t *filename){
723     // fill the modules with the sorted by module hits; add hits from
724     // background if option=Add.
725     // Inputs:
726     //      Int_t evnt       Event to be processed.
727     //      Int_t bgrev      Background Hit tree number.
728     //      Int_t nmodules   Not used.
729     //      Option_t *option String indicating if merging hits or not. To
730     //                       merge hits set equal to "Add". Otherwise no
731     //                       background hits are considered.
732     //      Test_t *filename File name containing the background hits..
733     // Outputs:
734     //      none.
735     // Return:
736     //      none.
737     static TTree *trH1;                 //Tree with background hits
738     static Bool_t first=kTRUE;
739     static TFile *file;
740     const char *addBgr = strstr(option,"Add");
741
742     evnt = nmodules; // Dummy use of variables to remove warnings
743     if (addBgr ) {
744         if(first) {
745             file=new TFile(filename);
746         } // end if first
747         first=kFALSE;
748         file->cd();
749         file->ls();
750         // Get Hits Tree header from file
751         if(trH1) delete trH1;
752         trH1=0;
753
754         char treeName[20];
755         sprintf(treeName,"TreeH%d",bgrev);
756         trH1 = (TTree*)gDirectory->Get(treeName);
757         if (!trH1) {
758             Error("FillModules","cannot find Hits Tree for event:%d",bgrev);
759         } // end if !trH1
760         // Set branch addresses
761     } // end if addBgr
762
763     FillModules(fLoader->TreeH(),0); // fill from this file's tree.
764     
765     if (addBgr ) {
766         FillModules(trH1,10000000); // Default mask 10M.
767         TTree *fAli=fLoader->GetRunLoader()->TreeK();
768         TFile *fileAli=0;
769         if (fAli) fileAli =fAli->GetCurrentFile();
770         fileAli->cd();
771     } // end if add
772 }
773 //______________________________________________________________________
774 void AliITS::FillModules(TTree *treeH, Int_t mask) {
775     // fill the modules with the sorted by module hits; 
776     // can be called many times to do a merging
777     // Inputs:
778     //      TTree *treeH  The tree containing the hits to be copied into
779     //                    the modules.
780     //      Int_t mask    The track number mask to indecate which file
781     //                    this hits came from.
782     // Outputs:
783     //      none.
784     // Return:
785     //      none.
786
787     if (treeH == 0x0)
788      {
789        Error("FillModules","Tree is NULL");
790      }
791     Int_t lay,lad,det,index;
792     AliITShit *itsHit=0;
793     AliITSmodule *mod=0;
794     char branchname[20];
795     sprintf(branchname,"%s",GetName());
796     TBranch *branch = treeH->GetBranch(branchname);
797     if (!branch) {
798         Error("FillModules","%s branch in TreeH not found",branchname);
799         return;
800     } // end if !branch
801     branch->SetAddress(&fHits);
802     Int_t nTracks =(Int_t) treeH->GetEntries();
803     Int_t iPrimTrack,h;
804     for(iPrimTrack=0; iPrimTrack<nTracks; iPrimTrack++){
805         ResetHits();
806         Int_t nBytes = treeH->GetEvent(iPrimTrack);
807         if (nBytes <= 0) continue;
808         Int_t nHits = fHits->GetEntriesFast();
809         for(h=0; h<nHits; h++){
810             itsHit = (AliITShit *)fHits->UncheckedAt(h);
811             itsHit->GetDetectorID(lay,lad,det);
812             if (fITSgeom) {
813                 index = fITSgeom->GetModuleIndex(lay,lad,det);
814             } else {
815                 index=det-1; // This should not be used.
816             } // end if [You must have fITSgeom for this to work!]
817             mod = GetModule(index);
818             itsHit->SetTrack(itsHit->GetTrack()+mask); // Set track mask.
819             mod->AddHit(itsHit,iPrimTrack,h);
820         } // end loop over hits 
821     } // end loop over tracks
822 }
823 //______________________________________________________________________
824 void AliITS::ClearModules(){
825     // Clear the modules TObjArray.
826     // Inputs:
827     //      none.
828     // Outputs:
829     //      none.
830
831     if(fITSmodules) fITSmodules->Delete();
832 }
833 //______________________________________________________________________
834 void AliITS::MakeBranchS(const char *fl){
835     // Creates Tree Branch for the ITS summable digits.
836     // Inputs:
837     //      cont char *fl  File name where SDigits branch is to be written
838     //                     to. If blank it write the SDigits to the same
839     //                     file in which the Hits were found.
840     // Outputs:
841     //      none.
842     // Return:
843     //      none.
844     Int_t buffersize = 4000;
845     char branchname[30];
846
847     // only one branch for SDigits.
848     sprintf(branchname,"%s",GetName());
849     
850
851     if(fLoader->TreeS()){
852         if (fSDigits == 0x0)  fSDigits  = new TClonesArray("AliITSpListItem",1000);
853         MakeBranchInTree(fLoader->TreeS(),branchname,&fSDigits,buffersize,fl);
854     } // end if
855 }
856 //______________________________________________________________________
857 void AliITS::SetTreeAddressS(TTree *treeS){
858     // Set branch address for the ITS summable digits Trees.
859     // Inputs:
860     //      TTree *treeS   Tree containing the SDigits.
861     // Outputs:
862     //      none.
863     // Return:
864     //      none.
865     char branchname[30];
866
867     if(!treeS) return;
868     if (fSDigits == 0x0)  fSDigits = new TClonesArray("AliITSpListItem",1000);
869     TBranch *branch;
870     sprintf(branchname,"%s",GetName());
871     branch = treeS->GetBranch(branchname);
872     if (branch) branch->SetAddress(&fSDigits);
873 }
874 //______________________________________________________________________
875 void AliITS::MakeBranchInTreeD(TTree *treeD,const char *file){
876     // Creates Tree branches for the ITS.
877     // Inputs:
878     //      TTree     *treeD Pointer to the Digits Tree.
879     //      cont char *file  File name where Digits branch is to be written
880     //                       to. If blank it write the SDigits to the same
881     //                       file in which the Hits were found.
882     // Outputs:
883     //      none.
884     // Return:
885     //      none.
886     Int_t buffersize = 4000;
887     char branchname[30];
888
889     sprintf(branchname,"%s",GetName());
890     // one branch for digits per type of detector
891     const char *det[3] = {"SPD","SDD","SSD"};
892     char digclass[40];
893     char clclass[40];
894     Int_t i;
895     for (i=0; i<kNTYPES ;i++) {
896         DetType(i)->GetClassNames(digclass,clclass);
897         // digits
898         if (fDtype == 0x0) fDtype = new TObjArray(fNDetTypes);
899         if(!(fDtype->At(i))) fDtype->AddAt(new TClonesArray(digclass,1000),i);
900         else ResetDigits(i);
901     } // end for i
902     for (i=0; i<kNTYPES ;i++) {
903         if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
904         else  sprintf(branchname,"%sDigits%d",GetName(),i+1);      
905         if (fDtype && treeD) {
906             MakeBranchInTree(treeD, branchname, &((*fDtype)[i]),buffersize,file);
907         } // end if
908     } // end for i
909 }
910 //______________________________________________________________________
911 void AliITS::SetTreeAddressD(TTree *treeD){
912     // Set branch address for the Trees.
913     // Inputs:
914     //      TTree *treeD   Tree containing the Digits.
915     // Outputs:
916     //      none.
917     // Return:
918     //      none.
919     char branchname[30];
920     const char *det[3] = {"SPD","SDD","SSD"};
921     TBranch *branch;
922     char digclass[40];
923     char clclass[40];
924     Int_t i;
925
926     if(!treeD) return;
927     for (i=0; i<kNTYPES; i++) {
928         DetType(i)->GetClassNames(digclass,clclass);
929         // digits
930         if (fDtype == 0x0) fDtype = new TObjArray(fNDetTypes);
931         if(!(fDtype->At(i))) fDtype->AddAt(new TClonesArray(digclass,1000),i);
932         else ResetDigits(i);
933         if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
934         else  sprintf(branchname,"%sDigits%d",GetName(),i+1);
935         if (fDtype) {
936             branch = treeD->GetBranch(branchname);
937             if (branch) branch->SetAddress(&((*fDtype)[i]));
938         } // end if fDtype
939     } // end for i
940 }
941 //______________________________________________________________________
942 void AliITS::Hits2SDigits(){
943     // Standard Hits to summable Digits function.
944     // Inputs:
945     //      none.
946     // Outputs:
947     //      none.
948
949 //    return; // Using Hits in place of the larger sDigits.
950     AliRunLoader* rl = fLoader->GetRunLoader(); 
951     AliHeader *header=rl->GetHeader(); // Get event number from this file.
952     if (header == 0x0)
953      {
954        rl->LoadHeader();
955        header=rl->GetHeader();
956        if (header == 0x0) return;
957      }
958     // Do the Hits to Digits operation. Use Standard input values.
959     // Event number from file, no background hit merging , use size from
960     // AliITSgeom class, option="All", input from this file only.
961     HitsToSDigits(header->GetEvent(),0,-1," ",fOpt," ");
962     
963 }
964 //______________________________________________________________________
965 void AliITS::Hits2PreDigits(){
966     // Standard Hits to summable Digits function.
967     // Inputs:
968     //      none.
969     // Outputs:
970     //      none.
971
972     AliHeader *header=fLoader->GetRunLoader()->GetHeader(); // Get event number from this file.
973     // Do the Hits to Digits operation. Use Standard input values.
974     // Event number from file, no background hit merging , use size from
975     // AliITSgeom class, option="All", input from this file only.
976     HitsToPreDigits(header->GetEvent(),0,-1," ",fOpt," ");
977 }
978 //______________________________________________________________________
979 void AliITS::SDigitsToDigits(Option_t *opt){
980     // Standard Summable digits to Digits function.
981     // Inputs:
982     //      none.
983     // Outputs:
984     //      none.
985     char name[20] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
986
987     if(!GetITSgeom()) return; // need transformations to do digitization.
988     AliITSgeom *geom = GetITSgeom();
989
990     const char *all = strstr(opt,"All");
991     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
992                           strstr(opt,"SSD")};
993     if( !det[0] && !det[1] && !det[2] ) all = "All";
994     else all = 0;
995     static Bool_t setDef=kTRUE;
996     if (setDef) SetDefaultSimulation();
997     setDef=kFALSE;
998
999     AliITSsimulation *sim      = 0;
1000     AliITSDetType    *iDetType = 0;
1001     TTree            *trees    = fLoader->TreeS();
1002     if( !(trees && this->GetSDigits()) ){
1003         Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
1004         return;
1005     } // end if
1006     sprintf( name, "%s", this->GetName() );
1007     TBranch *brchSDigits = trees->GetBranch( name );
1008     
1009     Int_t id,module;
1010     for(module=0;module<geom->GetIndexMax();module++){
1011         id       = geom->GetModuleType(module);
1012         if (!all && !det[id]) continue;
1013         iDetType = DetType(id);
1014         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1015         if (!sim) {
1016             Error("SDigit2Digits",
1017                   "The simulation class was not instanciated!");
1018             exit(1);
1019         } // end if !sim
1020         sim->InitSimulationModule(module,gAlice->GetEvNumber());
1021 //
1022         // add summable digits to module
1023         this->GetSDigits()->Clear();
1024         brchSDigits->GetEvent(module);
1025         sim->AddSDigitsToModule(GetSDigits(),0);
1026 //
1027         // Digitise current module sum(SDigits)->Digits
1028         sim->FinishSDigitiseModule();
1029
1030         // fills all branches - wasted disk space
1031         fLoader->TreeD()->Fill();
1032         this->ResetDigits();
1033     } // end for module
1034
1035     fLoader->TreeD()->GetEntries();
1036
1037     fLoader->TreeD()->AutoSave();
1038     // reset tree
1039     fLoader->TreeD()->Reset();
1040     
1041 }
1042 //______________________________________________________________________
1043 void AliITS::Hits2Digits(){
1044     // Standard Hits to Digits function.
1045     // Inputs:
1046     //      none.
1047     // Outputs:
1048     //      none.
1049
1050     AliHeader *header=fLoader->GetRunLoader()->GetHeader(); // Get event number from this file.
1051     // Do the Hits to Digits operation. Use Standard input values.
1052     // Event number from file, no background hit merging , use size from
1053     // AliITSgeom class, option="All", input from this file only.
1054     HitsToDigits(header->GetEvent(),0,-1," ",fOpt," ");
1055 }
1056 //______________________________________________________________________
1057 void AliITS::HitsToSDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1058                           Option_t *option, Option_t *opt,Text_t *filename){
1059     // keep galice.root for signal and name differently the file for 
1060     // background when add! otherwise the track info for signal will be lost !
1061     // the condition below will disappear when the geom class will be
1062     // initialized for all versions - for the moment it is only for v5 !
1063     // 7 is the SDD beam test version. Dummy routine. Hits are ITS's Summable
1064     // Digits.
1065     // Inputs:
1066     //      Int_t evnt       Event to be processed.
1067     //      Int_t bgrev      Background Hit tree number.
1068     //      Int_t nmodules   Not used.
1069     //      Option_t *option String indicating if merging hits or not. To
1070     //                       merge hits set equal to "Add". Otherwise no
1071     //                       background hits are considered.
1072     //      Test_t *filename File name containing the background hits..
1073     // Outputs:
1074     //      none.
1075     // Return:
1076     //      none.
1077 //    return; // using Hits instead of the larger sdigits.
1078
1079     HitsToPreDigits(evNumber,bgrev,size,option,opt,filename);
1080 }
1081 //______________________________________________________________________
1082 void AliITS::HitsToPreDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1083                           Option_t *option, Option_t *opt,Text_t *filename){
1084     //   Keep galice.root for signal and name differently the file for 
1085     // background when add! otherwise the track info for signal will be lost !
1086     // the condition below will disappear when the geom class will be
1087     // initialized for all versions - for the moment it is only for v5 !
1088     // 7 is the SDD beam test version.
1089     // Inputs:
1090     //      Int_t evnt       Event to be processed.
1091     //      Int_t bgrev      Background Hit tree number.
1092     //      Int_t nmodules   Not used.
1093     //      Option_t *option String indicating if merging hits or not. To
1094     //                       merge hits set equal to "Add". Otherwise no
1095     //                       background hits are considered.
1096     //      Test_t *filename File name containing the background hits..
1097     // Outputs:
1098     //      none.
1099     // Return:
1100     //      none.
1101
1102     if(!GetITSgeom()) return; // need transformations to do digitization.
1103     AliITSgeom *geom = GetITSgeom();
1104
1105     const char *all = strstr(opt,"All");
1106     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1107                           strstr(opt,"SSD")};
1108     static Bool_t setDef=kTRUE;
1109     if (setDef) SetDefaultSimulation();
1110     setDef=kFALSE;
1111
1112     Int_t nmodules;
1113     InitModules(size,nmodules);
1114     FillModules(evNumber,bgrev,nmodules,option,filename);
1115
1116     AliITSsimulation *sim      = 0;
1117     AliITSDetType    *iDetType = 0;
1118     AliITSmodule     *mod      = 0;
1119     Int_t id,module;
1120     for(module=0;module<geom->GetIndexMax();module++){
1121         id       = geom->GetModuleType(module);
1122         if (!all && !det[id]) continue;
1123         iDetType = DetType(id);
1124         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1125         if (!sim) {
1126             Error("HitsToSDigits",
1127                   "The simulation class was not instanciated!");
1128             exit(1);
1129         } // end if !sim
1130         mod      = (AliITSmodule *)fITSmodules->At(module);
1131         sim->SDigitiseModule(mod,module,evNumber);
1132         // fills all branches - wasted disk space
1133         fLoader->TreeS()->Fill(); 
1134         ResetSDigits();
1135     } // end for module
1136
1137     ClearModules();
1138
1139     fLoader->TreeS()->GetEntries();
1140     fLoader->TreeS()->AutoSave();
1141     fLoader->WriteSDigits("OVERWRITE");
1142     // reset tree
1143     fLoader->TreeS()->Reset();
1144 }
1145 //______________________________________________________________________
1146 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1147                           Option_t *option, Option_t *opt,Text_t *filename){
1148     //   Keep galice.root for signal and name differently the file for 
1149     // background when add! otherwise the track info for signal will be lost !
1150     // the condition below will disappear when the geom class will be
1151     // initialized for all versions - for the moment it is only for v5 !
1152     // 7 is the SDD beam test version.
1153     // Inputs:
1154     //      Int_t evnt       Event to be processed.
1155     //      Int_t bgrev      Background Hit tree number.
1156     //      Int_t nmodules   Not used.
1157     //      Option_t *option String indicating if merging hits or not. To
1158     //                       merge hits set equal to "Add". Otherwise no
1159     //                       background hits are considered.
1160     //      Test_t *filename File name containing the background hits..
1161     // Outputs:
1162     //      none.
1163     // Return:
1164     //      none.
1165
1166     if(!GetITSgeom()) return; // need transformations to do digitization.
1167     AliITSgeom *geom = GetITSgeom();
1168
1169     const char *all = strstr(opt,"All");
1170     const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1171                           strstr(opt,"SSD")};
1172     static Bool_t setDef=kTRUE;
1173     if (setDef) SetDefaultSimulation();
1174     setDef=kFALSE;
1175
1176     Int_t nmodules;
1177     InitModules(size,nmodules);
1178     FillModules(evNumber,bgrev,nmodules,option,filename);
1179
1180     AliITSsimulation *sim      = 0;
1181     AliITSDetType    *iDetType = 0;
1182     AliITSmodule     *mod      = 0;
1183     Int_t id,module;
1184     for(module=0;module<geom->GetIndexMax();module++){
1185         id       = geom->GetModuleType(module);
1186         if (!all && !det[id]) continue;
1187         iDetType = DetType(id);
1188         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1189         if (!sim) {
1190             Error("HitsToDigits",
1191                   "The simulation class was not instanciated!");
1192             exit(1);
1193         } // end if !sim
1194         mod      = (AliITSmodule *)fITSmodules->At(module);
1195         sim->DigitiseModule(mod,module,evNumber);
1196         // fills all branches - wasted disk space
1197         fLoader->TreeD()->Fill(); 
1198         ResetDigits();
1199     } // end for module
1200
1201     ClearModules();
1202
1203     fLoader->TreeD()->GetEntries();
1204     fLoader->TreeD()->AutoSave();
1205     // reset tree
1206     fLoader->TreeD()->Reset();
1207 }
1208 //______________________________________________________________________
1209 void AliITS::ResetSDigits(){
1210     // Reset the Summable Digits array.
1211     // Inputs:
1212     //      none.
1213     // Outputs:
1214     //      none.
1215
1216     if (fSDigits) fSDigits->Clear();
1217     fNSDigits = 0;
1218 }
1219 //______________________________________________________________________
1220 void AliITS::ResetDigits(){
1221     // Reset number of digits and the digits array for the ITS detector.
1222     // Inputs:
1223     //      none.
1224     // Outputs:
1225     //      none.
1226
1227     if (!fDtype) return;
1228
1229     Int_t i;
1230     for (i=0;i<kNTYPES;i++ ) {
1231         if (fDtype->At(i))    ((TClonesArray*)fDtype->At(i))->Clear();
1232         if (fNdtype)  fNdtype[i]=0;
1233     } // end for i
1234 }
1235 //______________________________________________________________________
1236 void AliITS::ResetDigits(Int_t i){
1237     // Reset number of digits and the digits array for this branch.
1238     // Inputs:
1239     //      none.
1240     // Outputs:
1241     //      none.
1242
1243     if (fDtype->At(i))    ((TClonesArray*)fDtype->At(i))->Clear();
1244     if (fNdtype)  fNdtype[i]=0;
1245 }
1246 //______________________________________________________________________
1247 void AliITS::AddSumDigit(AliITSpListItem &sdig){
1248     // Adds the a module full of summable digits to the summable digits tree.
1249     // Inputs:
1250     //      AliITSpListItem &sdig   SDigit to be added to SDigits tree.
1251     // Outputs:
1252     //      none.
1253     // Return:
1254     //      none.
1255
1256     TClonesArray &lsdig = *fSDigits;
1257     new(lsdig[fNSDigits++]) AliITSpListItem(sdig);
1258 }
1259 //______________________________________________________________________
1260 void AliITS::AddRealDigit(Int_t id, Int_t *digits){
1261     //   Add a real digit - as coming from data.
1262     // Inputs:
1263     //      Int_t id        Detector type number.
1264     //      Int_t *digits   Integer array containing the digits info. See 
1265     //                      AliITSdigit.h
1266     // Outputs:
1267     //      none.
1268     // Return:
1269     //      none.
1270
1271     TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1272     new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
1273 }
1274 //______________________________________________________________________
1275 void AliITS::AddSimDigit(Int_t id, AliITSdigit *d){
1276     //    Add a simulated digit.
1277     // Inputs:
1278     //      Int_t id        Detector type number.
1279     //      AliITSdigit *d  Digit to be added to the Digits Tree. See 
1280     //                      AliITSdigit.h
1281     // Outputs:
1282     //      none.
1283     // Return:
1284     //      none.
1285
1286     TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1287
1288     switch(id){
1289     case 0:
1290         new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
1291         break;
1292     case 1:
1293         new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
1294         break;
1295     case 2:
1296         new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
1297         break;
1298     } // end switch id
1299 }
1300 //______________________________________________________________________
1301 void AliITS::AddSimDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,
1302                          Int_t *hits,Float_t *charges){
1303     //   Add a simulated digit to the list.
1304     // Inputs:
1305     //      Int_t id        Detector type number.
1306     //      Float_t phys    Physics indicator. See AliITSdigits.h
1307     //      Int_t *digits   Integer array containing the digits info. See 
1308     //                      AliITSdigit.h
1309     //      Int_t *tracks   Integer array [AliITSdigitS?D::GetNTracks()] 
1310     //                      containing the track numbers that contributed to
1311     //                      this digit.
1312     //      Int_t *hits     Integer array [AliITSdigitS?D::GetNTracks()]
1313     //                      containing the hit numbers, from AliITSmodule, that
1314     //                      contributed to this digit.
1315     //      Float_t *charge Floating point array of the signals contributed
1316     //                      to this digit by each track.
1317     // Outputs:
1318     //      none.
1319     // Return:
1320     //      none.
1321
1322     TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1323     switch(id){
1324     case 0:
1325         new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks,hits);
1326         break;
1327     case 1:
1328         new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,
1329                                                    hits,charges);
1330         break;
1331     case 2:
1332         new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks,hits);
1333         break;
1334     } // end switch id
1335 }
1336 //______________________________________________________________________
1337 void AliITS::MakeTreeC(Option_t *option){
1338   //   Create a separate tree to store the clusters.
1339   // Inputs:
1340   //      Option_t *option  string which must contain "C" otherwise
1341   //                        no Cluster Tree is created.
1342   // Outputs:
1343   //      none.
1344   // Return:
1345   //      none.
1346
1347   AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;    
1348     
1349   if (pITSLoader == 0x0) {
1350     Error("MakeTreeC","fLoader == 0x0 option=%s",option);
1351     return;
1352   }
1353   if (pITSLoader->TreeC() == 0x0) pITSLoader->MakeTree("C");
1354   MakeBranchC();
1355 }
1356
1357 void AliITS::MakeBranchC()
1358 {
1359 //Makes barnches in treeC
1360   AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;    
1361   if (pITSLoader == 0x0) 
1362    {
1363     Error("MakeTreeC","fLoader == 0x0");
1364     return;
1365    }
1366   TTree * lTC = pITSLoader->TreeC();
1367   if (lTC == 0x0)
1368    {
1369      Error("MakeTreeC","Can not get TreeC from Loader");
1370      return;
1371    }
1372
1373   Int_t buffersize = 4000;
1374   char branchname[30];
1375   const char *det[3] = {"SPD","SDD","SSD"};
1376   char digclass[40];
1377   char clclass[40];
1378
1379     // one branch for Clusters per type of detector
1380   Int_t i;   
1381   for (i=0; i<kNTYPES ;i++) 
1382     {
1383         AliITSDetType *iDetType=DetType(i); 
1384         iDetType->GetClassNames(digclass,clclass);
1385         // clusters
1386         if (fCtype == 0x0) fCtype  = new TObjArray(fNDetTypes);
1387         if(!ClustersAddress(i))
1388          {
1389           fCtype->AddAt(new TClonesArray(clclass,1000),i);
1390          }
1391         if (kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
1392         else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
1393         if (fCtype  && lTC) 
1394          {
1395            if (lTC->GetBranch(branchname))
1396             {
1397               Warning("MakeBranchC","Branch %s alread exists in TreeC",branchname);
1398             }
1399            else
1400             {
1401               Info("MakeBranchC","Creating branch %s in TreeC",branchname);
1402               lTC->Branch(branchname,&((*fCtype)[i]), buffersize);
1403             }
1404          } // end if fCtype && lTC
1405   } // end for i
1406 }
1407
1408 //______________________________________________________________________
1409 void AliITS::GetTreeC(Int_t event){
1410   //    Get the clusters tree for this event and set the branch address.
1411   // Inputs:
1412   //      Int_t event    Event number for the cluster tree.
1413   // Outputs:
1414   //      none.
1415   // Return:
1416   //      none.
1417   char branchname[30];
1418   const char *det[3] = {"SPD","SDD","SSD"};
1419
1420   AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;
1421   TTree * lTC = pITSLoader->TreeC();
1422
1423   ResetClusters();
1424   if (lTC) {
1425     pITSLoader->CleanRawClusters();
1426   } // end if TreeC()
1427
1428
1429   TBranch *branch;
1430
1431   if (lTC) {
1432     Int_t i;
1433         char digclass[40];
1434         char clclass[40];
1435         for (i=0; i<kNTYPES; i++) {
1436       AliITSDetType *iDetType=DetType(i); 
1437       iDetType->GetClassNames(digclass,clclass);
1438       // clusters
1439       if (fCtype == 0x0) fCtype  = new TObjArray(fNDetTypes);
1440       if(!fCtype->At(i)) fCtype->AddAt(new TClonesArray(clclass,1000),i);
1441       if(kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
1442       else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
1443       if (fCtype) {
1444                 branch = lTC->GetBranch(branchname);
1445         if (branch) branch->SetAddress(&((*fCtype)[i]));
1446       } // end if fCtype
1447         } // end for i
1448   } else {
1449         Error("GetTreeC","cannot find Clusters Tree for event:%d",event);
1450   } // end if lTC
1451 }
1452 //______________________________________________________________________
1453 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c){
1454     //   Add a cluster to the list.
1455     // Inputs:
1456     //      Int_t id             Detector type number.
1457     //      AliITSRawCluster *c  Cluster class to be added to the tree of
1458     //                           clusters.
1459     // Outputs:
1460     //      none.
1461     // Return:
1462     //      none.
1463
1464     TClonesArray &lc = *((TClonesArray*)fCtype->At(id));
1465
1466     switch(id){
1467     case 0:
1468         new(lc[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
1469         break;
1470     case 1:
1471         new(lc[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
1472         break;
1473     case 2:
1474         new(lc[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
1475         break;
1476     } // end switch id
1477 }
1478 //______________________________________________________________________
1479 void AliITS::ResetClusters(){
1480     // Reset number of clusters and the clusters array for ITS.
1481     // Inputs:
1482     //      none.
1483     // Outputs:
1484     //      none.
1485
1486     Int_t i;
1487     for (i=0;i<kNTYPES;i++ ) ResetClusters(i);
1488 }
1489 //______________________________________________________________________
1490 void AliITS::ResetClusters(Int_t i){
1491     //    Reset number of clusters and the clusters array for this branch.
1492     // Inputs:
1493     //      Int_t i        Detector type number.
1494     // Outputs:
1495     //      none.
1496     // Return:
1497     //      none.
1498
1499     if (fCtype->At(i))    ((TClonesArray*)fCtype->At(i))->Clear();
1500     if (fNctype)  fNctype[i]=0;
1501 }
1502 //______________________________________________________________________
1503 void AliITS::MakeBranchR(const char *file, Option_t *opt){
1504     // Creates Tree branches for the ITS Reconstructed points.
1505     // Inputs:
1506     //      cont char *file  File name where RecPoints branch is to be written
1507     //                       to. If blank it write the SDigits to the same
1508     //                       file in which the Hits were found.
1509     // Outputs:
1510     //      none.
1511     // Return:
1512     //      none.
1513     Int_t buffsz = 4000;
1514     char branchname[30];
1515
1516     // only one branch for rec points for all detector types
1517     Bool_t oFast= (strstr(opt,"Fast")!=0);
1518     if(oFast){
1519       sprintf(branchname,"%sRecPointsF",GetName());
1520     } else {
1521       sprintf(branchname,"%sRecPoints",GetName());
1522     }
1523     
1524     
1525     if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1526     if (fLoader->TreeR()) {
1527         if (fRecPoints == 0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1528         MakeBranchInTree(fLoader->TreeR(),branchname,&fRecPoints,buffsz,file);
1529     } // end if
1530 }
1531 //______________________________________________________________________
1532 void AliITS::SetTreeAddressR(TTree *treeR){
1533     // Set branch address for the Reconstructed points Trees.
1534     // Inputs:
1535     //      TTree *treeR   Tree containing the RecPoints.
1536     // Outputs:
1537     //      none.
1538     // Return:
1539     //      none.
1540     char branchname[30];
1541
1542     if(!treeR) return;
1543     if (fRecPoints == 0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1544     TBranch *branch;
1545     sprintf(branchname,"%sRecPoints",GetName());
1546     branch = treeR->GetBranch(branchname);
1547     if (branch) {
1548       branch->SetAddress(&fRecPoints);
1549     }
1550     else {
1551       sprintf(branchname,"%sRecPointsF",GetName());
1552       branch = treeR->GetBranch(branchname);
1553       if (branch) {
1554         branch->SetAddress(&fRecPoints);
1555       }
1556     }
1557 }
1558 //______________________________________________________________________
1559 void AliITS::AddRecPoint(const AliITSRecPoint &r){
1560     // Add a reconstructed space point to the list
1561     // Inputs:
1562     //      const AliITSRecPoint &r RecPoint class to be added to the tree
1563     //                              of reconstructed points TreeR.
1564     // Outputs:
1565     //      none.
1566     // Return:
1567     //      none.
1568
1569     TClonesArray &lrecp = *fRecPoints;
1570     new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
1571 }
1572 //______________________________________________________________________
1573 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
1574                                   Option_t *opt0,Option_t *opt1,Text_t *flnm){
1575     // keep galice.root for signal and name differently the file for 
1576     // background when add! otherwise the track info for signal will be lost !
1577     // the condition below will disappear when the geom class will be
1578     // initialized for all versions - for the moment it is only for v5 !
1579     // Inputs:
1580     //      Int_t evnt       Event to be processed.
1581     //      Int_t bgrev      Background Hit tree number.
1582     //      Int_t size       Size used by InitModules. See InitModules.
1583     //      Option_t *opt0   Option passed to FillModules. See FillModules.
1584     //      Option_t *opt1   String indicating if merging hits or not. To
1585     //                       merge hits set equal to "Add". Otherwise no
1586     //                       background hits are considered.
1587     //      Test_t *flnm     File name containing the background hits..
1588     // Outputs:
1589     //      none.
1590     // Return:
1591     //      none.
1592
1593     if(!GetITSgeom()) return;
1594     AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
1595     AliITSgeom *geom = GetITSgeom();
1596
1597     const char *all = strstr(opt1,"All");
1598     const char *det[3] ={strstr(opt1,"SPD"),strstr(opt1,"SDD"),
1599                          strstr(opt1,"SSD")};
1600     Int_t nmodules;
1601     InitModules(size,nmodules);
1602     FillModules(evNumber,bgrev,nmodules,opt0,flnm);
1603
1604     AliITSsimulation *sim      = 0;
1605     AliITSDetType    *iDetType = 0;
1606     AliITSmodule     *mod      = 0;
1607     Int_t id,module;
1608
1609     //m.b. : this change is nothing but a nice way to make sure
1610     //the CPU goes up !
1611     
1612     cout<<"HitsToFastRecPoints: N mod = "<<geom->GetIndexMax()<<endl;
1613     
1614     for(module=0;module<geom->GetIndexMax();module++)
1615      {
1616         id       = geom->GetModuleType(module);
1617         if (!all && !det[id]) continue;
1618         iDetType = DetType(id);
1619         sim      = (AliITSsimulation*)iDetType->GetSimulationModel();
1620         if (!sim) 
1621          {
1622            Error("HitsToFastPoints","The simulation class was not instanciated!");
1623            exit(1);
1624          } // end if !sim
1625         mod      = (AliITSmodule *)fITSmodules->At(module);
1626         sim->CreateFastRecPoints(mod,module,gRandom);
1627         cout<<module<<"\r";fflush(0);
1628         //gAlice->TreeR()->Fill();
1629         TTree *lTR = pITSloader->TreeR();
1630         TBranch *br=lTR->GetBranch("ITSRecPointsF");
1631         br->Fill();
1632         ResetRecPoints();
1633     } // end for module
1634
1635     ClearModules();
1636     
1637     fLoader->WriteRecPoints("OVERWRITE");
1638 }
1639 //______________________________________________________________________
1640 void AliITS::Digits2Reco(){
1641     // Find clusters and reconstruct space points.
1642     // Inputs:
1643     //      none.
1644     // Outputs:
1645     //      none.
1646
1647     AliHeader *header=fLoader->GetRunLoader()->GetHeader();
1648     // to Digits to RecPoints for event in file, all digits in file, and
1649     // all ITS detectors.
1650     DigitsToRecPoints(header->GetEvent(),0,fOpt);
1651 }
1652 //______________________________________________________________________
1653 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt){
1654   // cluster finding and reconstruction of space points
1655   // the condition below will disappear when the geom class will be
1656   // initialized for all versions - for the moment it is only for v5 !
1657   // 7 is the SDD beam test version
1658   // Inputs:
1659   //      Int_t evNumber   Event number to be processed.
1660   //      Int_t lastentry  Offset for module when not all of the modules
1661   //                       are processed.
1662   //      Option_t *opt    String indicating which ITS sub-detectors should
1663   //                       be processed. If ="All" then all of the ITS
1664   //                       sub detectors are processed.
1665   // Outputs:
1666   //      none.
1667   // Return:
1668   //      none.
1669
1670   if(!GetITSgeom()) return;
1671   AliITSgeom *geom = GetITSgeom();
1672     
1673   const char *all = strstr(opt,"All");
1674   const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1675                         strstr(opt,"SSD")};
1676   static Bool_t setRec=kTRUE;
1677   if (setRec) SetDefaultClusterFinders();
1678   setRec=kFALSE;
1679
1680   AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
1681   TTree *treeC=pITSloader->TreeC();
1682   AliITSClusterFinder *rec     = 0;
1683   AliITSDetType      *iDetType = 0;
1684   Int_t id,module,first=0;
1685   for(module=0;module<geom->GetIndexMax();module++){
1686       id       = geom->GetModuleType(module);
1687       if (!all && !det[id]) continue;
1688       if(det[id]) first = geom->GetStartDet(id);
1689       iDetType = DetType(id);
1690       rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1691       TClonesArray *itsDigits  = this->DigitsAddress(id);
1692       if (!rec) {
1693           Error("DigitsToRecPoints",
1694                 "The reconstruction class was not instanciated! event=%d",
1695                 evNumber);
1696           exit(1);
1697       } // end if !rec
1698       this->ResetDigits();
1699       TTree *lTD = pITSloader->TreeD();
1700       if (all) {
1701           lTD->GetEvent(lastentry+module);
1702       }else {
1703           lTD->GetEvent(lastentry+(module-first));
1704       }
1705       Int_t ndigits = itsDigits->GetEntriesFast();
1706       if (ndigits) rec->FindRawClusters(module);
1707       pITSloader->TreeR()->Fill(); 
1708       ResetRecPoints();
1709       treeC->Fill();
1710       ResetClusters();
1711   } // end for module
1712
1713
1714   pITSloader->WriteRecPoints("OVERWRITE");
1715   pITSloader->WriteRawClusters("OVERWRITE");
1716 }
1717 //______________________________________________________________________
1718 void AliITS::ResetRecPoints(){
1719     // Reset number of rec points and the rec points array.
1720     // Inputs:
1721     //      none.
1722     // Outputs:
1723     //      none.
1724
1725     if (fRecPoints) fRecPoints->Clear();
1726     fNRecPoints = 0;
1727 }
1728 //______________________________________________________________________
1729 AliLoader* AliITS::MakeLoader(const char* topfoldername)
1730
1731   //builds ITSgetter (AliLoader type)
1732   //if detector wants to use castomized getter, it must overload this method
1733
1734   Info("MakeLoader","Creating AliITSLoader. Top folder is %s.",topfoldername);
1735   fLoader = new AliITSLoader(GetName(),topfoldername);
1736   return fLoader;
1737 }
1738
1739