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