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