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