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