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