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