]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITS.cxx
stdlib.h included (HP,Sun)
[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.39  2001/03/12 17:45:32  hristov
19 Changes needed on Sun with CC 5.0
20
21 Revision 1.38  2001/03/07 14:04:51  barbera
22 Some vector dimensions increased to cope with full events
23
24 Revision 1.37  2001/03/07 12:36:35  barbera
25 A change added in the tracking part to manage delta rays
26
27 Revision 1.36  2001/03/02 19:44:11  barbera
28  modified to taking into account new version tracking v1
29
30 Revision 1.35  2001/02/28 18:16:46  mariana
31 Make the code compatible with the new AliRun
32
33 Revision 1.34  2001/02/11 15:51:39  mariana
34 Set protection in MakeBranch
35
36 Revision 1.33  2001/02/10 22:26:39  mariana
37 Move the initialization of the containers for raw clusters in MakeTreeC()
38
39 Revision 1.32  2001/02/08 23:55:31  nilsen
40 Removed fMajor/MinorVersion variables in favor of variables in derived classes.
41 Set arrays char *det[3] = {"SPD","SDD","SSD"} as const.
42
43 Revision 1.31  2001/02/02 23:57:28  nilsen
44 Added include file that are no londer included in AliITSgeom.h
45
46 Revision 1.30  2001/01/30 09:23:13  hristov
47 Streamers removed (R.Brun)
48
49 Revision 1.29  2001/01/26 20:01:09  hristov
50 Major upgrade of AliRoot code
51
52 Revision 1.28  2000/12/18 14:02:00  barbera
53 new version of the ITS tracking to take into account the new TPC track parametrization
54
55 Revision 1.27  2000/12/08 13:49:27  barbera
56 Hidden declaration in a for loop removed to be compliant with HP-UX compiler
57
58 Revision 1.26  2000/11/27 13:12:13  barbera
59 New version containing the files for tracking
60
61 Revision 1.25  2000/11/12 22:38:05  barbera
62 Added header file for the SPD Bari model
63
64 Revision 1.24  2000/10/09 22:18:12  barbera
65 Bug fixes from MAriana to le AliITStest.C run correctly
66
67 Revision 1.23  2000/10/05 20:47:42  nilsen
68 fixed dependencies of include files. Tryed but failed to get a root automaticly
69 generates streamer function to work. Modified SetDefaults.
70
71 Revision 1.9.2.15  2000/10/04 16:56:40  nilsen
72 Needed to include stdlib.h
73
74 =======
75 Revision 1.22  2000/10/04 19:45:52  barbera
76 Corrected by F. Carminati for v3.04
77
78 Revision 1.21  2000/10/02 21:28:08  fca
79 Removal of useless dependecies via forward declarations
80
81 Revision 1.20  2000/10/02 16:31:39  barbera
82 General code clean-up
83
84 Revision 1.9.2.14  2000/10/02 15:43:51  barbera
85 General code clean-up (e.g., printf -> cout)
86
87 Revision 1.19  2000/09/22 12:13:25  nilsen
88 Patches and updates for fixes to this and other routines.
89
90 Revision 1.18  2000/07/12 05:32:20  fca
91 Correcting several syntax problem with static members
92
93 Revision 1.17  2000/07/10 16:07:18  fca
94 Release version of ITS code
95
96 Revision 1.9.2.3  2000/02/02 13:42:09  barbera
97 fixed AliITS.cxx for new AliRun structure. Added ITS hits list to list of hits which will have their track numbers updated
98
99 Revision 1.9.2.2  2000/01/23 03:03:13  nilsen
100 //fixed FillModule. Removed fi(fabs(xl)<dx....
101
102 Revision 1.9.2.1  2000/01/12 19:03:32  nilsen
103 This is the version of the files after the merging done in December 1999.
104 See the ReadMe110100.txt file for details
105
106 Revision 1.9  1999/11/14 14:33:25  fca
107 Correct problems with distructors and pointers, thanks to I.Hrivnacova
108
109 Revision 1.8  1999/09/29 09:24:19  fca
110 Introduction of the Copyright and cvs Log
111
112 */
113
114 ///////////////////////////////////////////////////////////////////////////////
115 //
116 //      An overview of the basic philosophy of the ITS code development
117 // and analysis is show in the figure below.
118 //Begin_Html
119 /*
120 <img src="picts/ITS/ITS_Analysis_schema.gif">
121 </pre>
122 <br clear=left>
123 <font size=+2 color=red>
124 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
125 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
126 </font>
127 <pre>
128 */
129 //End_Html
130 //
131 //  AliITS. Inner Traking System base class.
132 //  This class contains the base procedures for the Inner Tracking System
133 //
134 //Begin_Html
135 /*
136 <img src="picts/ITS/AliITS_Class_Diagram.gif">
137 </pre>
138 <br clear=left>
139 <font size=+2 color=red>
140 <p>This show the class diagram of the different elements that are part of
141 the AliITS class.
142 </font>
143 <pre>
144 */
145 //End_Html
146 //
147 // Version: 0
148 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
149 //
150 // Version: 1
151 // Modified and documented by Bjorn S. Nilsen
152 // July 11 1999
153 //
154 // Version: 2
155 // Modified and documented by A. Bologna
156 // October 18 1999
157 //
158 // AliITS is the general base class for the ITS. Also see AliDetector for
159 // futher information.
160 //
161 ///////////////////////////////////////////////////////////////////////////////
162 #include <iostream.h>
163 #include <iomanip.h>
164 #include <fstream.h>
165 #include <stdlib.h>
166 #include <TMath.h>
167 #include <TRandom.h>
168 #include <TBranch.h>
169 #include <TVector.h>
170 #include <TClonesArray.h>
171 #include <TROOT.h>
172 #include <TObjectTable.h>
173 #include <TFile.h>
174 #include <TTree.h>
175 #include <TString.h>
176 #include <TParticle.h>
177
178
179 #include "AliRun.h"
180 #include "AliITS.h"
181 #include "AliITSMap.h"
182 #include "AliITSDetType.h"
183 #include "AliITSClusterFinder.h"
184 //#include "AliITSsimulation.h"
185 #include "AliITSsimulationSPD.h"
186 #include "AliITSsimulationSDD.h"
187 #include "AliITSsimulationSSD.h"
188 #include "AliITSresponse.h"
189 #include "AliITSsegmentationSPD.h"
190 #include "AliITSresponseSPD.h"
191 #include "AliITSresponseSPDbari.h"
192 #include "AliITSsegmentationSDD.h"
193 #include "AliITSresponseSDD.h"
194 #include "AliITSsegmentationSSD.h"
195 #include "AliITSresponseSSD.h"
196 #include "AliITShit.h"
197 #include "AliITSgeom.h"
198 #include "AliITSdigit.h"
199 #include "AliITSmodule.h"
200 #include "AliITSRecPoint.h"
201 #include "AliITSRawCluster.h"
202 #include "AliMC.h"
203 #include "stdlib.h"
204
205 #include "AliITStrack.h"
206 #include "AliITSiotrack.h"
207 #include "AliITStracking.h"
208 #include "AliITSRad.h"   
209 #include "../TPC/AliTPC.h"
210 #include "../TPC/AliTPCParam.h"
211
212
213 ClassImp(AliITS)
214  
215 //_____________________________________________________________________________
216 AliITS::AliITS() : AliDetector() {
217   //
218   // Default initialiser for ITS
219   //     The default constructor of the AliITS class. In addition to
220   // creating the AliITS class it zeros the variables fIshunt (a member
221   // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
222   // fITSpoints, fIdSens, and fIdName. The AliDetector default constructor
223   // is also called.
224   //
225
226
227   fIshunt     = 0;
228   fEuclidOut  = 0;
229
230   fNDetTypes = kNTYPES;
231   fIdN        = 0;
232   fIdName     = 0;
233   fIdSens     = 0;
234   fITSmodules = 0;
235   //
236   fDetTypes   = 0;
237   //
238   fDtype  = 0;
239   fNdtype = 0;
240   fCtype  = 0;
241   fNctype = 0;
242   fRecPoints = 0;
243   fNRecPoints = 0;
244   fTreeC = 0;
245   //
246   fITSgeom=0;
247 }
248
249 //_____________________________________________________________________________
250 AliITS::AliITS(const char *name, const char *title):AliDetector(name,title){
251   //
252   // Default initialiser for ITS
253   //     The constructor of the AliITS class. In addition to creating the
254   // AliITS class, it allocates memory for the TClonesArrays fHits and
255   // fDigits, and for the TObjArray fITSpoints. It also zeros the variables
256   // fIshunt (a member of AliDetector class), fEuclidOut, and fIdN, and zeros
257   // the pointers fIdSens and fIdName. To help in displaying hits via the ROOT
258   // macro display.C AliITS also sets the marker color to red. The variables
259   // passes with this constructor, const char *name and *title, are used by
260   // the constructor of AliDetector class. See AliDetector class for a
261   // description of these parameters and its constructor functions.
262   //
263
264
265   fHits       = new TClonesArray("AliITShit", 1560);
266   gAlice->AddHitList(fHits);
267
268   fNDetTypes = kNTYPES;
269
270   fNdtype = new Int_t[kNTYPES];
271   fDtype = new TObjArray(kNTYPES);
272
273   fNctype = new Int_t[kNTYPES];
274   fCtype = new TObjArray(kNTYPES);
275
276
277   fRecPoints = 0;
278   fNRecPoints = 0;
279
280   fTreeC = 0;
281
282   fITSmodules = 0; 
283
284   fIshunt     = 0;
285   fEuclidOut  = 0;
286   fIdN        = 0;
287   fIdName     = 0;
288   fIdSens     = 0;
289  
290   fDetTypes = new TObjArray(kNTYPES);  
291
292   Int_t i;
293   for(i=0;i<kNTYPES;i++) {
294     (*fDetTypes)[i]=new AliITSDetType(); 
295     fNdtype[i]=0;
296     fNctype[i]=0;
297    }
298   //
299
300   SetMarkerColor(kRed);
301
302   fITSgeom=0;
303 }
304 //___________________________________________________________________________
305 AliITS::AliITS(AliITS &source){
306   // copy constructor
307   if(this==&source) return;
308   Error("AliITS::Copy constructor",
309         "You are not allowed to make a copy of the AliITS");
310   exit(1);
311 }
312 //____________________________________________________________________________
313 AliITS& AliITS::operator=(AliITS &source){
314   // assignment operator
315   if(this==&source) return *this;
316   Error("AliITS::operator=",
317         "You are not allowed to make a copy of the AliITS");
318   exit(1);
319   return *this; //fake return
320 }
321 //____________________________________________________________________________
322 void AliITS::ClearModules(){
323   //clear the modules TObjArray
324
325   if(fITSmodules) fITSmodules->Delete();
326
327 }
328 //_____________________________________________________________________________
329 AliITS::~AliITS(){
330   //
331   // Default distructor for ITS
332   //     The default destructor of the AliITS class. In addition to deleting
333   // the AliITS class it deletes the memory pointed to by the fHits, fDigits,
334   // fIdSens, fIdName, and fITSpoints.
335   //
336
337
338   delete fHits;
339   delete fDigits;
340   delete fRecPoints;
341 //  delete fIdName;        // TObjArray of TObjStrings
342   if(fIdName!=0) delete[] fIdName;  // Array of TStrings
343   if(fIdSens!=0) delete[] fIdSens;
344   if(fITSmodules!=0) {
345       this->ClearModules();
346       delete fITSmodules;
347   }// end if fITSmodules!=0
348
349   //
350   if(fDtype) {
351     fDtype->Delete();
352     delete fDtype;
353   }
354   delete [] fNdtype;
355   if (fCtype) {
356     fCtype->Delete();
357     delete fCtype;
358   }
359   delete [] fNctype;
360   //
361
362   if (fDetTypes) {
363     fDetTypes->Delete();
364     delete fDetTypes;
365   }
366
367   if (fTreeC) delete fTreeC;
368
369   if (fITSgeom) delete fITSgeom;
370
371 }
372
373 //___________________________________________
374 AliITSDetType* AliITS::DetType(Int_t id)
375 {
376   //return pointer to id detector type
377     return ((AliITSDetType*) (*fDetTypes)[id]);
378
379 }
380 //___________________________________________
381 void AliITS::SetClasses(Int_t id, const char *digit, const char *cluster)
382 {
383   //set the digit and cluster classes to be used for the id detector type
384     ((AliITSDetType*) (*fDetTypes)[id])->ClassNames(digit,cluster);
385
386 }
387 //___________________________________________
388 void AliITS::SetResponseModel(Int_t id, AliITSresponse *response)
389 {
390   //set the response model for the id detector type
391
392     ((AliITSDetType*) (*fDetTypes)[id])->ResponseModel(response);
393
394 }
395
396 //___________________________________________
397 void AliITS::SetSegmentationModel(Int_t id, AliITSsegmentation *seg)
398 {
399   //set the segmentation model for the id detector type
400
401     ((AliITSDetType*) (*fDetTypes)[id])->SegmentationModel(seg);
402
403 }
404
405 //___________________________________________
406 void AliITS::SetSimulationModel(Int_t id, AliITSsimulation *sim)
407 {
408   //set the simulation model for the id detector type
409
410    ((AliITSDetType*) (*fDetTypes)[id])->SimulationModel(sim);
411
412 }
413 //___________________________________________
414 void AliITS::SetReconstructionModel(Int_t id, AliITSClusterFinder *reconst)
415 {
416   //set the cluster finder model for the id detector type
417
418    ((AliITSDetType*) (*fDetTypes)[id])->ReconstructionModel(reconst);
419
420 }
421
422 //_____________________________________________________________________________
423 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
424   //
425   // Add an ITS hit
426   //     The function to add information to the AliITShit class. See the
427   // AliITShit class for a full description. This function allocates the
428   // necessary new space for the hit information and passes the variable
429   // track, and the pointers *vol and *hits to the AliITShit constructor
430   // function.
431   //
432   TClonesArray &lhits = *fHits;
433   new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
434 }
435 //_____________________________________________________________________________
436 void AliITS::AddRealDigit(Int_t id, Int_t *digits) 
437 {
438   // add a real digit - as coming from data
439
440   TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
441   new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
442
443 }
444 //_____________________________________________________________________________
445 void AliITS::AddSimDigit(Int_t id, AliITSdigit *d) 
446 {
447
448   // add a simulated digit
449
450   TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
451
452   switch(id)
453   {
454   case 0:
455      new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
456      break;
457   case 1:
458      new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
459      break;
460   case 2:
461      new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
462      break;
463   }
464
465 }
466
467 //_____________________________________________________________________________
468 void AliITS::AddSimDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,Int_t *hits,Float_t *charges){
469
470   // add a simulated digit to the list
471
472   TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
473   switch(id)
474   {
475   case 0:
476      new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks,hits);
477      break;
478   case 1:
479      new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,hits,charges);
480      break;
481   case 2:
482      new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks,hits);
483      break;
484   }
485  
486 }
487
488 //_____________________________________________________________________________
489 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c) 
490 {
491
492   // add a cluster to the list
493
494   TClonesArray &lcl = *((TClonesArray*)(*fCtype)[id]);
495
496   switch(id)
497   {
498   case 0:
499      new(lcl[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
500      break;
501   case 1:
502      new(lcl[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
503      break;
504   case 2:
505      new(lcl[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
506      break;
507   }
508
509 }
510
511
512 //_____________________________________________________________________________
513 void AliITS::AddRecPoint(const AliITSRecPoint &r)
514 {
515   //
516   // Add a reconstructed space point to the list
517   //
518   TClonesArray &lrecp = *fRecPoints;
519   new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
520 }
521  
522
523 //____________________________________________
524 void AliITS::ResetDigits()
525 {
526     //
527     // Reset number of digits and the digits array for the ITS detector
528     //
529
530     if (!fDtype) return;
531
532     Int_t i;
533     for (i=0;i<kNTYPES;i++ ) {
534         if ((*fDtype)[i])    ((TClonesArray*)(*fDtype)[i])->Clear();
535         if (fNdtype)  fNdtype[i]=0;
536     }
537 }
538
539 //____________________________________________
540 void AliITS::ResetDigits(Int_t i)
541 {
542     //
543     // Reset number of digits and the digits array for this branch
544     //
545   if ((*fDtype)[i])    ((TClonesArray*)(*fDtype)[i])->Clear();
546   if (fNdtype)  fNdtype[i]=0;
547 }
548
549
550 //____________________________________________
551 void AliITS::ResetClusters()
552 {
553     //
554     // Reset number of clusters and the clusters array for ITS
555     //
556
557     Int_t i;
558     for (i=0;i<kNTYPES;i++ ) {
559         if ((*fCtype)[i])    ((TClonesArray*)(*fCtype)[i])->Clear();
560         if (fNctype)  fNctype[i]=0;
561     }
562
563 }
564
565 //____________________________________________
566 void AliITS::ResetClusters(Int_t i)
567 {
568     //
569     // Reset number of clusters and the clusters array for this branch
570     //
571         if ((*fCtype)[i])    ((TClonesArray*)(*fCtype)[i])->Clear();
572         if (fNctype)  fNctype[i]=0;
573
574 }
575
576
577 //____________________________________________
578 void AliITS::ResetRecPoints()
579 {
580     //
581     // Reset number of rec points and the rec points array 
582     //
583     if (fRecPoints) fRecPoints->Clear();
584     fNRecPoints = 0;
585
586 }
587
588 //_____________________________________________________________________________
589 Int_t AliITS::DistancetoPrimitive(Int_t , Int_t ){
590   //
591   // Distance from mouse to ITS on the screen. Dummy routine
592   //     A dummy routine used by the ROOT macro display.C to allow for the
593   // use of the mouse (pointing device) in the macro. In general this should
594   // never be called. If it is it returns the number 9999 for any value of
595   // x and y.
596   //
597   return 9999;
598 }
599
600 //_____________________________________________________________________________
601 void AliITS::Init(){
602   //
603   // Initialise ITS after it has been built
604   //     This routine initializes the AliITS class. It is intended to be called
605   // from the Init function in AliITSv?. Besides displaying a banner
606   // indicating that it has been called it initializes the array fIdSens
607   // and sets the default segmentation, response, digit and raw cluster classes
608   // Therefore it should be called after a call to CreateGeometry.
609   //
610   Int_t i;
611
612 //
613   SetDefaults();
614 // Array of TStrings
615   for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
616 //
617 }
618
619 //_____________________________________________________________________________
620 void AliITS::SetDefaults()
621 {
622   // sets the default segmentation, response, digit and raw cluster classes
623
624   printf("SetDefaults\n");
625
626   AliITSDetType *iDetType;
627
628
629   //SPD 
630
631   iDetType=DetType(0); 
632   if (!iDetType->GetSegmentationModel()) {
633     AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(fITSgeom);
634     SetSegmentationModel(0,seg0); 
635   }
636   if (!iDetType->GetResponseModel()) {
637      SetResponseModel(0,new AliITSresponseSPD()); 
638   }
639   // set digit and raw cluster classes to be used
640   
641   const char *kData0=(iDetType->GetResponseModel())->DataType();
642   if (strstr(kData0,"real")) {
643       iDetType->ClassNames("AliITSdigit","AliITSRawClusterSPD");
644   } else iDetType->ClassNames("AliITSdigitSPD","AliITSRawClusterSPD");
645
646   // SDD                                          //
647   iDetType=DetType(1); 
648   if (!iDetType->GetResponseModel()) {
649     SetResponseModel(1,new AliITSresponseSDD()); 
650   }
651   AliITSresponse *resp1=iDetType->GetResponseModel();
652   if (!iDetType->GetSegmentationModel()) {
653     AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(fITSgeom,resp1);
654     SetSegmentationModel(1,seg1); 
655   }
656   const char *kData1=(iDetType->GetResponseModel())->DataType();
657   const char *kopt=iDetType->GetResponseModel()->ZeroSuppOption();
658   if ((!strstr(kopt,"2D")) && (!strstr(kopt,"1D")) || strstr(kData1,"real") ) {
659       iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
660   } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
661
662   // SSD
663   iDetType=DetType(2); 
664   if (!iDetType->GetSegmentationModel()) {
665     AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(fITSgeom);
666     SetSegmentationModel(2,seg2); 
667   }
668   if (!iDetType->GetResponseModel()) {
669     SetResponseModel(2,new AliITSresponseSSD()); 
670   }
671   const char *kData2=(iDetType->GetResponseModel())->DataType();
672   if (strstr(kData2,"real")) {
673       iDetType->ClassNames("AliITSdigit","AliITSRawClusterSSD");
674   } else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
675
676   if (kNTYPES>3) {
677     Warning("SetDefaults","Only the three basic detector types are initialised!");
678   } 
679
680 }
681
682
683
684 //_____________________________________________________________________________
685 void AliITS::SetDefaultSimulation()
686 {
687   // to be written
688
689 }
690 //_____________________________________________________________________________
691 void AliITS::SetDefaultClusterFinders()
692 {
693   // to be written
694
695 }
696 //_____________________________________________________________________________
697
698 void AliITS::MakeTreeC(Option_t *option)
699 {
700   // create a separate tree to store the clusters
701
702   cout << "AliITS::MakeTreeC" << endl;
703
704      const char *optC = strstr(option,"C");
705      if (optC && !fTreeC) fTreeC = new TTree("TC","Clusters in ITS");
706      else return;
707
708      Int_t buffersize = 4000;
709      char branchname[30];
710
711      const char *det[3] = {"SPD","SDD","SSD"};
712
713      char digclass[40];
714      char clclass[40];
715
716      // one branch for Clusters per type of detector
717      Int_t i;
718      for (i=0; i<kNTYPES ;i++) {
719        AliITSDetType *iDetType=DetType(i); 
720        iDetType->GetClassNames(digclass,clclass);
721        // clusters
722        (*fCtype)[i] = new TClonesArray(clclass,10000); 
723         if (kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
724         else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
725         if (fCtype   && fTreeC) {
726            TreeC()->Branch(branchname,&((*fCtype)[i]), buffersize);
727            cout << "Making Branch " << branchname;
728            cout << " for Clusters of detector type " << i+1 << endl;
729         }       
730      }
731
732 }
733
734 //_____________________________________________________________________________
735 void AliITS::GetTreeC(Int_t event)
736 {
737
738   cout << "AliITS::GetTreeC" << endl;
739
740   // get the clusters tree for this event and set the branch address
741     char treeName[20];
742     char branchname[30];
743
744     const char *det[3] = {"SPD","SDD","SSD"};
745
746     ResetClusters();
747     if (fTreeC) {
748           delete fTreeC;
749     }
750
751     sprintf(treeName,"TreeC%d",event);
752     fTreeC = (TTree*)gDirectory->Get(treeName);
753
754     TBranch *branch;
755     if (fTreeC) {
756         Int_t i;
757         for (i=0; i<kNTYPES; i++) {
758            if (kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
759            else  sprintf(branchname,"%sClusters%d",GetName(),i+1);
760            if (fCtype) {
761                 branch = fTreeC->GetBranch(branchname);
762                 if (branch) branch->SetAddress(&((*fCtype)[i]));
763             }
764         }
765     } else {
766         Error("AliITS::GetTreeC",
767                 "cannot find Clusters Tree for event:%d\n",event);
768     }
769
770 }
771 //_____________________________________________________________________________
772 void AliITS::MakeBranch(Option_t* option, char *file)
773 {
774   //
775   // Creates Tree branches for the ITS.
776   //
777   //
778   Int_t buffersize = 4000;
779   char branchname[30];
780   sprintf(branchname,"%s",GetName());
781
782   AliDetector::MakeBranch(option,file);
783
784   const char *cD = strstr(option,"D");
785   const char *cR = strstr(option,"R");
786
787   if (cD) {
788   //
789   // one branch for digits per type of detector
790   //
791    const char *det[3] = {"SPD","SDD","SSD"};
792
793    char digclass[40];
794    char clclass[40];
795
796    Int_t i;
797    for (i=0; i<kNTYPES ;i++) {
798        AliITSDetType *iDetType=DetType(i); 
799        iDetType->GetClassNames(digclass,clclass);
800        // digits
801        if(!((*fDtype)[i])) (*fDtype)[i] = new TClonesArray(digclass,10000);
802        else ResetDigits(i);
803    }
804
805    for (i=0; i<kNTYPES ;i++) {
806       if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
807       else  sprintf(branchname,"%sDigits%d",GetName(),i+1);      
808       if (fDtype && gAlice->TreeD()) {
809         gAlice->MakeBranchInTree(gAlice->TreeD(), 
810                          branchname, &((*fDtype)[i]), buffersize, file);
811             cout << "Making Branch " << branchname;
812             cout << " for digits of type "<< i+1 << endl;
813       } 
814     }
815   }
816
817   if (cR) {
818   //
819   // only one branch for rec points for all detector types
820   //
821     sprintf(branchname,"%sRecPoints",GetName());
822  
823     if(!fRecPoints) fRecPoints=new TClonesArray("AliITSRecPoint",10000);
824
825     if (fRecPoints && gAlice->TreeR()) {
826       gAlice->MakeBranchInTree(gAlice->TreeR(), 
827                                branchname, &fRecPoints, buffersize, file) ;
828       cout << "Making Branch " << branchname;
829       cout << " for reconstructed space points" << endl;
830     }
831   }     
832 }
833
834 //___________________________________________
835 void AliITS::SetTreeAddress()
836 {
837
838   // Set branch address for the Trees.
839
840   char branchname[30];
841   AliDetector::SetTreeAddress();
842
843   const char *det[3] = {"SPD","SDD","SSD"};
844
845   TBranch *branch;
846   TTree *treeD = gAlice->TreeD();
847   TTree *treeR = gAlice->TreeR();
848
849   Int_t i;
850   if (treeD) {
851       for (i=0; i<kNTYPES; i++) {
852           if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
853           else  sprintf(branchname,"%sDigits%d",GetName(),i+1);
854           if (fDtype) {
855               branch = treeD->GetBranch(branchname);
856               if (branch) branch->SetAddress(&((*fDtype)[i]));
857           }
858       }
859   }
860
861  
862   if (treeR) {
863     sprintf(branchname,"%sRecPoints",GetName());
864       branch = treeR->GetBranch(branchname);
865       if (branch) branch->SetAddress(&fRecPoints);
866   }
867   
868
869 }
870
871 //____________________________________________________________________________
872 void AliITS::InitModules(Int_t size,Int_t &nmodules){
873
874   //initialize the modules array
875
876   if(fITSmodules){ 
877       fITSmodules->Delete();
878       delete fITSmodules;
879   }
880
881     Int_t nl,indexMAX,index;
882
883     if(size<=0){ // default to using data stored in AliITSgeom
884         if(fITSgeom==0) {
885             Error("AliITS::InitModules",
886                 "in AliITS::InitModule fITSgeom not defined\n");
887             return;
888         } // end if fITSgeom==0
889         nl = fITSgeom->GetNlayers();
890         indexMAX = fITSgeom->GetModuleIndex(nl,fITSgeom->GetNladders(nl),
891                                             fITSgeom->GetNdetectors(nl))+1;
892         nmodules = indexMAX;
893         fITSmodules = new TObjArray(indexMAX);
894         for(index=0;index<indexMAX;index++){
895                 fITSmodules->AddAt( new AliITSmodule(index),index);
896         } // end for index
897     }else{
898         fITSmodules = new TObjArray(size);
899         for(index=0;index<size;index++) {
900             fITSmodules->AddAt( new AliITSmodule(index),index);
901         }
902
903         nmodules = size;
904     } // end i size<=0
905 }
906
907 //____________________________________________________________________________
908 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t nmodules,Option_t *option,Text_t *filename){
909
910   // fill the modules with the sorted by module hits; add hits from background
911   // if option=Add
912
913
914     static TTree *trH1;                 //Tree with background hits
915     static TClonesArray *fHits2;        //List of hits for one track only
916
917     static Bool_t first=kTRUE;
918     static TFile *file;
919     const char *addBgr = strstr(option,"Add");
920
921
922     if (addBgr ) {
923         if(first) {
924             cout<<"filename "<<filename<<endl;
925             file=new TFile(filename);
926             cout<<"I have opened "<<filename<<" file "<<endl;
927             fHits2     = new TClonesArray("AliITShit",1000  );
928         }           
929         first=kFALSE;
930         file->cd();
931         file->ls();
932         // Get Hits Tree header from file
933         if(fHits2) fHits2->Clear();
934         if(trH1) delete trH1;
935         trH1=0;
936         
937         char treeName[20];
938         sprintf(treeName,"TreeH%d",bgrev);
939         trH1 = (TTree*)gDirectory->Get(treeName);
940         //printf("TrH1 %p of treename %s for event %d \n",trH1,treeName,bgrev);
941         
942         if (!trH1) {
943             Error("AliITS::FillModules",
944             "cannot find Hits Tree for event:%d\n",bgrev);
945         }
946         // Set branch addresses
947         TBranch *branch;
948         char branchname[20];
949         sprintf(branchname,"%s",GetName());
950         if (trH1 && fHits2) {
951             branch = trH1->GetBranch(branchname);
952             if (branch) branch->SetAddress(&fHits2);
953         }
954
955         // test
956         //Int_t ntracks1 =(Int_t)TrH1->GetEntries();
957         //printf("background - ntracks1 - %d\n",ntracks1);
958    }
959
960     //Int_t npart = gAlice->GetEvent(evnt);
961     //if(npart<=0) return;
962     TClonesArray *itsHits = this->Hits();
963     Int_t lay,lad,det,index;
964     AliITShit *itsHit=0;
965     AliITSmodule *mod=0;
966
967     TTree *iTH = gAlice->TreeH();
968     Int_t ntracks =(Int_t) iTH->GetEntries();
969
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         //printf("nhits %d\n",nhits);
976         if (!nhits) continue;
977         for(h=0; h<nhits; h++){
978             itsHit = (AliITShit *)itsHits->UncheckedAt(h);
979             itsHit->GetDetectorID(lay,lad,det);
980             // temporarily index=det-1 !!!
981             if(fITSgeom) index = fITSgeom->GetModuleIndex(lay,lad,det);
982             else index=det-1;
983             //
984             mod = this->GetModule(index);
985             mod->AddHit(itsHit,t,h);
986         } // end loop over hits 
987     } // end loop over tracks
988
989     // open the file with background
990     
991     if (addBgr ) {
992           Int_t track,i;
993           ntracks =(Int_t)trH1->GetEntries();
994             //printf("background - ntracks1 %d\n",ntracks);
995             //printf("background - Start loop over tracks \n");     
996             //   Loop over tracks
997
998             for (track=0; track<ntracks; track++) {
999
1000                 if (fHits2)       fHits2->Clear();
1001                 trH1->GetEvent(track);
1002                 //   Loop over hits
1003                 for(i=0;i<fHits2->GetEntriesFast();++i) {
1004
1005                     itsHit=(AliITShit*) (*fHits2)[i];
1006                     itsHit->GetDetectorID(lay,lad,det);
1007                     // temporarily index=det-1 !!!
1008                     if(fITSgeom) index = fITSgeom->GetModuleIndex(lay,lad,det);
1009                     else index=det-1;
1010                     //
1011                     mod = this->GetModule(index);
1012                     mod->AddHit(itsHit,track,i);
1013                }  // end loop over hits
1014             } // end loop over tracks
1015
1016             TTree *fAli=gAlice->TreeK();
1017             TFile *fileAli=0;
1018             
1019             if (fAli) fileAli =fAli->GetCurrentFile();
1020             fileAli->cd();
1021
1022     } // end if add
1023
1024     //gObjectTable->Print();
1025
1026 }
1027
1028 //____________________________________________________________________________
1029
1030 void AliITS::SDigits2Digits()
1031 {
1032   
1033   AliITSgeom *geom = GetITSgeom();
1034   
1035   // SPD
1036   AliITSDetType *iDetType;
1037   iDetType=DetType(0);
1038   AliITSsegmentationSPD *seg0=(AliITSsegmentationSPD*)iDetType->GetSegmentationModel();
1039   AliITSresponseSPD *res0 = (AliITSresponseSPD*)iDetType->GetResponseModel();
1040   AliITSsimulationSPD *sim0=new AliITSsimulationSPD(seg0,res0);
1041   SetSimulationModel(0,sim0);
1042   // test
1043   // printf("SPD dimensions %f %f \n",seg0->Dx(),seg0->Dz());
1044   // printf("SPD npixels %d %d \n",seg0->Npz(),seg0->Npx());
1045   // printf("SPD pitches %d %d \n",seg0->Dpz(0),seg0->Dpx(0));
1046   // end test
1047   // 
1048   // SDD
1049   //Set response functions
1050   // SDD compression param: 2 fDecrease, 2fTmin, 2fTmax or disable, 2 fTolerance
1051
1052   iDetType=DetType(1);
1053   AliITSresponseSDD *res1 = (AliITSresponseSDD*)iDetType->GetResponseModel();
1054   if (!res1) {
1055     res1=new AliITSresponseSDD();
1056     SetResponseModel(1,res1);
1057   }
1058   Float_t noise, baseline;
1059   res1->GetNoiseParam(noise,baseline);
1060   Float_t noise_after_el =  res1->GetNoiseAfterElectronics();
1061   Float_t fCutAmp = baseline + 2.*noise_after_el;
1062   Int_t cp[8]={0,0,(int)fCutAmp,(int)fCutAmp,0,0,0,0}; //1D
1063   res1->SetCompressParam(cp);
1064   AliITSsegmentationSDD *seg1=(AliITSsegmentationSDD*)iDetType->GetSegmentationModel();
1065   if (!seg1) {
1066     seg1 = new AliITSsegmentationSDD(geom,res1);
1067     SetSegmentationModel(1,seg1);
1068   }
1069   AliITSsimulationSDD *sim1=new AliITSsimulationSDD(seg1,res1);
1070   SetSimulationModel(1,sim1);
1071
1072   // SSD
1073   iDetType=DetType(2);
1074   AliITSsegmentationSSD *seg2=(AliITSsegmentationSSD*)iDetType->GetSegmentationModel();
1075   AliITSresponseSSD *res2 = (AliITSresponseSSD*)iDetType->GetResponseModel();
1076   res2->SetSigmaSpread(3.,2.);
1077   AliITSsimulationSSD *sim2=new AliITSsimulationSSD(seg2,res2);
1078   SetSimulationModel(2,sim2);
1079  
1080   cerr<<"Digitizing ITS...\n";
1081  
1082   TStopwatch timer;
1083   timer.Start();
1084   HitsToDigits(0,0,-1," ","All"," ");
1085   timer.Stop(); timer.Print();
1086
1087   delete sim0;
1088   delete sim1;
1089   delete sim2;
1090 }
1091
1092
1093 //____________________________________________________________________________
1094 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size, Option_t *option, Option_t *opt,Text_t *filename)
1095 {
1096     // keep galice.root for signal and name differently the file for 
1097     // background when add! otherwise the track info for signal will be lost !
1098   
1099    // the condition below will disappear when the geom class will be
1100    // initialised for all versions - for the moment it is only for v5 !
1101    // 7 is the SDD beam test version  
1102    Int_t ver = this->IsVersion(); 
1103    if(ver!=5 && ver!=7) return; 
1104
1105    const char *all = strstr(opt,"All");
1106    const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1107
1108    Int_t nmodules;
1109    InitModules(size,nmodules); 
1110    FillModules(evNumber,bgrev,nmodules,option,filename);
1111
1112    //TBranch *branch;
1113    AliITSsimulation* sim;
1114    //TObjArray *branches=gAlice->TreeD()->GetListOfBranches();
1115    AliITSgeom *geom = GetITSgeom();
1116
1117    Int_t id,module;
1118    Int_t first,last;
1119    for (id=0;id<kNTYPES;id++) {
1120         if (!all && !det[id]) continue;
1121         //branch = (TBranch*)branches->UncheckedAt(id);
1122         AliITSDetType *iDetType=DetType(id); 
1123         sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1124         if (!sim) {
1125            Error("HitsToDigits","The simulation class was not instantiated!");
1126            exit(1);
1127            // or SetDefaultSimulation();
1128         }
1129         if(geom) {
1130           first = geom->GetStartDet(id);
1131           last = geom->GetLastDet(id);
1132         } else first=last=0;
1133         cout << "det type " << id << " first, last "<< first << last << endl;
1134         for(module=first;module<=last;module++) {
1135             AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
1136             sim->DigitiseModule(mod,module,evNumber);
1137             // fills all branches - wasted disk space
1138             gAlice->TreeD()->Fill(); 
1139             ResetDigits();
1140             // try and fill only the branch 
1141             //branch->Fill();
1142             //ResetDigits(id);
1143         } // loop over modules
1144    } // loop over detector types
1145
1146    ClearModules();
1147
1148    Int_t nentries=(Int_t)gAlice->TreeD()->GetEntries();
1149    cout << "nentries in TreeD" << nentries << endl;
1150
1151    char hname[30];
1152    sprintf(hname,"TreeD%d",evNumber);
1153    gAlice->TreeD()->Write(hname,TObject::kOverwrite);
1154    // reset tree
1155    gAlice->TreeD()->Reset();
1156
1157 }
1158
1159
1160 //____________________________________________________________________________
1161 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt)
1162 {
1163   // cluster finding and reconstruction of space points
1164   
1165    // the condition below will disappear when the geom class will be
1166    // initialised for all versions - for the moment it is only for v5 !
1167    // 7 is the SDD beam test version  
1168    Int_t ver = this->IsVersion(); 
1169    if(ver!=5) return; 
1170
1171    const char *all = strstr(opt,"All");
1172    const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1173
1174    static Bool_t first=kTRUE;
1175    if (!TreeC() && first) {
1176        MakeTreeC("C");
1177        first=kFALSE;
1178    }
1179
1180    TTree *treeC=TreeC();
1181  
1182
1183    //TBranch *branch;
1184    AliITSClusterFinder* rec;
1185
1186    //TObjArray *branches=gAlice->TreeR()->GetListOfBranches();
1187    AliITSgeom *geom = GetITSgeom();
1188
1189    Int_t id,module;
1190    for (id=0;id<kNTYPES;id++) {
1191         if (!all && !det[id]) continue;
1192         //branch = (TBranch*)branches->UncheckedAt(id);
1193         AliITSDetType *iDetType=DetType(id); 
1194         rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1195         if (!rec) {
1196            Error("DigitsToRecPoints","The cluster finder class was not instantiated!");
1197            exit(1);
1198            // or SetDefaultClusterFinders();
1199         }
1200         TClonesArray *itsDigits  = this->DigitsAddress(id);
1201
1202         Int_t first,last;
1203         if(geom) {
1204           first = geom->GetStartDet(id);
1205           last = geom->GetLastDet(id);
1206         } else first=last=0;
1207         //printf("first last %d %d\n",first,last);
1208         for(module=first;module<=last;module++) {
1209               this->ResetDigits();
1210               if (all) gAlice->TreeD()->GetEvent(lastentry+module);
1211               else gAlice->TreeD()->GetEvent(lastentry+(module-first));
1212               Int_t ndigits = itsDigits->GetEntriesFast();
1213               if (ndigits) rec->FindRawClusters();
1214               gAlice->TreeR()->Fill(); 
1215               ResetRecPoints();
1216               treeC->Fill();
1217               ResetClusters();
1218               // try and fill only the branch 
1219               //branch->Fill();
1220               //ResetRecPoints(id);
1221         } // loop over modules
1222    } // loop over detector types
1223
1224
1225    Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1226    Int_t ncentries=(Int_t)treeC->GetEntries();
1227    cout << " nentries ncentries " << nentries << ncentries <<  endl;
1228
1229    char hname[30];
1230    sprintf(hname,"TreeR%d",evNumber);
1231    gAlice->TreeR()->Write(hname,TObject::kOverwrite);
1232    // reset tree
1233    gAlice->TreeR()->Reset();
1234
1235    sprintf(hname,"TreeC%d",evNumber);
1236    treeC->Write(hname,TObject::kOverwrite);
1237    treeC->Reset();
1238 }
1239
1240
1241 //____________________________________________________________________________
1242 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
1243 Option_t *option,Option_t *opt,Text_t *filename)
1244 {
1245     // keep galice.root for signal and name differently the file for 
1246     // background when add! otherwise the track info for signal will be lost !
1247   
1248
1249    // the condition below will disappear when the geom class will be
1250    // initialised for all versions - for the moment it is only for v5 !  
1251    Int_t ver = this->IsVersion(); 
1252    if(ver!=5) return; 
1253
1254    const char *all = strstr(opt,"All");
1255    const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1256
1257    Int_t nmodules;
1258    InitModules(size,nmodules);
1259    FillModules(evNumber,bgrev,nmodules,option,filename);
1260
1261
1262    AliITSsimulation* sim;
1263    AliITSgeom *geom = GetITSgeom();
1264
1265    TRandom *random=new TRandom[9];
1266    random[0].SetSeed(111);
1267    random[1].SetSeed(222);
1268    random[2].SetSeed(333);              
1269    random[3].SetSeed(444);
1270    random[4].SetSeed(555);
1271    random[5].SetSeed(666);              
1272    random[6].SetSeed(777);
1273    random[7].SetSeed(888);
1274    random[8].SetSeed(999);              
1275
1276
1277    Int_t id,module;
1278    for (id=0;id<kNTYPES;id++) {
1279         if (!all && !det[id]) continue;
1280         AliITSDetType *iDetType=DetType(id); 
1281         sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1282         if (!sim) {
1283            Error("HitsToFastPoints","The simulation class was not instantiated!");
1284            exit(1);
1285            // or SetDefaultSimulation();
1286         }
1287         Int_t first = geom->GetStartDet(id);
1288         Int_t last = geom->GetLastDet(id);
1289         for(module=first;module<=last;module++) {
1290             AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
1291             sim->CreateFastRecPoints(mod,module,random);
1292             gAlice->TreeR()->Fill(); 
1293             ResetRecPoints();
1294         } // loop over modules
1295    } // loop over detector types
1296
1297
1298    ClearModules();
1299
1300    //Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1301
1302    char hname[30];
1303    sprintf(hname,"TreeR%d",evNumber);
1304    gAlice->TreeR()->Write(hname,TObject::kOverwrite);
1305    // reset tree
1306    gAlice->TreeR()->Reset();
1307
1308    delete [] random;
1309
1310 }
1311
1312 //________________________________________________________________
1313
1314 AliITStrack  AliITS::Tracking(AliITStrack &track, AliITStrack *reference,TObjArray *fastpoints, Int_t
1315 **vettid, Bool_t flagvert,  AliITSRad *rl  ) { 
1316                                                                                 
1317 //Origin  A. Badala' and G.S. Pappalardo:  e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it 
1318
1319                                                                                     
1320   TList *list= new TList();   
1321
1322   AliITStrack tr(track);
1323   
1324   list->AddLast(&tr);
1325   
1326   Double_t Pt=(tr).GetPt();
1327   //cout << "\n Pt = " << Pt <<"\n";
1328
1329   AliITStracking obj(list, reference, this, fastpoints,TMath::Abs(Pt),vettid, flagvert, rl);
1330   list->Delete();
1331   delete list;
1332
1333   Int_t itot=-1;
1334   TVector VecTotLabref(18);
1335   Int_t lay, k;
1336   for(lay=5; lay>=0; lay--) {
1337     TVector VecLabref(3); 
1338     VecLabref=(*reference).GetLabTrack(lay);
1339     Float_t ClustZ=(*reference).GetZclusterTrack( lay);   //modified il 5-3-2001         
1340     for(k=0; k<3; k++) { // {itot++; VecTotLabref(itot)=VecLabref(k);}  // modified 5-3-2002
1341  
1342         Int_t lpp=(Int_t)VecLabref(k);
1343         if(lpp>=0) {
1344           TParticle *p=(TParticle*) gAlice->Particle(lpp);
1345           Int_t pcode=p->GetPdgCode();
1346           if(pcode==11) VecLabref(k)=p->GetFirstMother();
1347         }    
1348     
1349         itot++; VecTotLabref(itot)=VecLabref(k);
1350         if(VecLabref(k)==0. && ClustZ == 0.) VecTotLabref(itot) =-3.;    
1351     }    
1352   }
1353   Long_t labref;
1354   Int_t freq;  
1355   (*reference).Search(VecTotLabref, labref, freq);
1356     
1357   if(freq < 5) labref=-labref;          
1358   (*reference).SetLabel(labref);
1359
1360   return *reference; 
1361
1362 }
1363
1364
1365
1366 //________________________________________________________________
1367
1368
1369
1370 void AliITS::DoTracking(Int_t evNumber, Int_t min_t, Int_t max_t, TFile *file, Bool_t flagvert) {
1371
1372 //   ex macro for tracking ITS
1373
1374   printf("begin DoTracking - file %p\n",file);
1375
1376   const char *pname="75x40_100x60";
1377   
1378  Int_t imax=200,jmax=450;
1379   AliITSRad *rl = new AliITSRad(imax,jmax);
1380   //cout<<" dopo costruttore AliITSRad\n"; getchar();
1381     
1382   struct GoodTrack {
1383     Int_t lab,code;
1384     Float_t px,py,pz,x,y,z,pxg,pyg,pzg,ptg;
1385     Bool_t flag;
1386   };
1387   
1388
1389   gAlice->GetEvent(0);
1390
1391   AliTPC *TPC=(AliTPC*)gAlice->GetDetector("TPC");
1392   AliTPCParam *digp = (AliTPCParam*)file->Get(pname);
1393   if (digp!=0) TPC->SetParam(digp);
1394   
1395   GoodTrack gt[15000];
1396   Int_t ngood=0;
1397   ifstream in("itsgood_tracks");
1398
1399   cerr<<"Reading itsgood tracks...\n";
1400   while (in>>gt[ngood].lab>>gt[ngood].code
1401           >>gt[ngood].px >>gt[ngood].py>>gt[ngood].pz
1402           >>gt[ngood].x  >>gt[ngood].y >>gt[ngood].z
1403           >>gt[ngood].pxg  >>gt[ngood].pyg >>gt[ngood].pzg
1404           >>gt[ngood].ptg >>gt[ngood].flag) {
1405     ngood++;
1406     cerr<<ngood<<'\r';
1407     if (ngood==15000) {
1408       cerr<<"Too many good tracks !\n";
1409       break;
1410     }
1411   }
1412   if (!in.eof()) cerr<<"Read error (itsgood_tracks) !\n";
1413   
1414   
1415 // Load tracks
1416   TFile *tf=TFile::Open("tpctracks.root");
1417   if (!tf->IsOpen()) {cerr<<"Can't open tpctracks.root !\n"; return ;}
1418   TObjArray tracks(200000);
1419   TTree *tracktree=(TTree*)tf->Get("TreeT");
1420   TBranch *tbranch=tracktree->GetBranch("tracks");
1421   Int_t nentr=(Int_t)tracktree->GetEntries();
1422   Int_t kk;
1423   for (kk=0; kk<nentr; kk++) {
1424     AliTPCtrack *iotrack=new AliTPCtrack;
1425     tbranch->SetAddress(&iotrack);
1426     tracktree->GetEvent(kk);
1427     tracks.AddLast(iotrack);
1428   }   
1429   tf->Close();
1430
1431
1432   Int_t nt = tracks.GetEntriesFast();
1433   cerr<<"Number of found tracks "<<nt<<endl;
1434   
1435   TVector DataOut(9);
1436   Int_t kkk=0;
1437   
1438   Double_t ptg=0.,pxg=0.,pyg=0.,pzg=0.;
1439
1440   //////////////////////////////  good tracks definition in TPC  ////////////////////////////////
1441       
1442   ofstream out1 ("AliITSTrag.out");
1443   Int_t i;
1444   for (i=0; i<ngood; i++) out1 << gt[i].ptg << "\n";
1445   out1.close();
1446
1447
1448   TVector vec(5);
1449   TTree *TR=gAlice->TreeR();
1450   Int_t nent=(Int_t)TR->GetEntries();
1451   TClonesArray  *recPoints = RecPoints();
1452   Int_t numbpoints;
1453   Int_t totalpoints=0;
1454   Int_t *np = new Int_t[nent];
1455   Int_t **vettid = new Int_t* [nent];
1456   Int_t mod;
1457   for (mod=0; mod<nent; mod++) {
1458     vettid[mod]=0;
1459     this->ResetRecPoints();
1460     //gAlice->TreeR()->GetEvent(mod+1); //first entry in TreeR is empty
1461     gAlice->TreeR()->GetEvent(mod); //first entry in TreeR is empty
1462     numbpoints = recPoints->GetEntries();
1463     totalpoints+=numbpoints;
1464     np[mod] = numbpoints;
1465   //cout<<" mod = "<<mod<<"   numbpoints = "<<numbpoints<<"\n"; getchar();
1466     vettid[mod] = new Int_t[numbpoints];
1467     Int_t ii;
1468     for (ii=0;ii<numbpoints; ii++) *(vettid[mod]+ii)=0;
1469   }
1470
1471   AliTPCtrack *track;
1472
1473      
1474   if(min_t < 0) {min_t = 0; max_t = nt-1;}   
1475
1476 /*
1477   ///////////////////////////////// Definition of vertex end its error ////////////////////////////
1478   ////////////////////////// In the future it will be given by a method ///////////////////////////
1479   Double_t Vx=0.;
1480   Double_t Vy=0.;
1481   Double_t Vz=0.;
1482   
1483   Float_t sigmavx=0.0050;      // 50  microns
1484   Float_t sigmavy=0.0050;      // 50  microns
1485   Float_t sigmavz=0.010;       // 100 microns
1486
1487   //Vx+=gRandom->Gaus(0,sigmavx);  Vy+=gRandom->Gaus(0,sigmavy);  Vz+=gRandom->Gaus(0,sigmavz);
1488   TVector vertex(3), ervertex(3)
1489   vertex(0)=Vx; vertex(1)=Vy; vertex(2)=Vz;
1490   ervertex(0)=sigmavx;  ervertex(1)=sigmavy;  ervertex(2)=sigmavz;
1491   /////////////////////////////////////////////////////////////////////////////////////////////////
1492 */      
1493
1494   //TDirectory *savedir=gDirectory; 
1495
1496   TTree tracktree1("TreeT","Tree with ITS tracks");
1497   AliITSiotrack *iotrack=0;
1498   tracktree1.Branch("ITStracks","AliITSiotrack",&iotrack,32000,0);
1499
1500   ofstream out ("AliITSTra.out");
1501   
1502   Int_t j;       
1503   for (j=min_t; j<=max_t; j++) {     
1504     track=(AliTPCtrack*)tracks.UncheckedAt(j);
1505     Int_t flaglab=0;
1506     if (!track) continue;
1507     ////// elimination of not good tracks ////////////   
1508     Int_t ilab=TMath::Abs(track->GetLabel());
1509     Int_t iii;
1510     for (iii=0;iii<ngood;iii++) {
1511          //cout<<" ilab, gt[iii].lab = "<<ilab<<" "<<gt[iii].lab<<"\n"; getchar();
1512       if (ilab==gt[iii].lab) { 
1513         flaglab=1;
1514         ptg=gt[iii].ptg; 
1515         pxg=gt[iii].pxg;
1516         pyg=gt[iii].pyg;
1517         pzg=gt[iii].pzg;        
1518         break;
1519       }
1520     }
1521          //cout<<" j flaglab =  " <<j<<" "<<flaglab<<"\n";  getchar();
1522     if (!flaglab) continue;
1523          //cout<<" j =  " <<j<<"\n";  getchar();
1524   /*            
1525     ////// old propagation to the end of TPC //////////////       
1526     Double_t xk=76.;
1527     track->PropagateTo(xk);
1528     xk-=0.11;
1529     track->PropagateTo(xk,42.7,2.27); //C
1530     xk-=2.6;
1531     track->PropagateTo(xk,36.2,1.98e-3); //C02
1532     xk-=0.051;
1533     track->PropagateTo(xk,42.7,2.27); //C 
1534     /////////////////////////////////////////////////// 
1535          */
1536                  
1537          ////// new propagation to the end of TPC //////////////
1538     Double_t xk=77.415;
1539     track->PropagateTo(xk, 28.94, 1.204e-3);     //Ne
1540          xk -=0.01;
1541     track->PropagateTo(xk, 44.77, 1.71);         //Tedlar
1542          xk -=0.04;
1543     track->PropagateTo(xk, 44.86, 1.45);         //Kevlar
1544          xk -=2.0;
1545     track->PropagateTo(xk, 41.28, 0.029);        //Nomex         
1546     xk-=16;
1547     track->PropagateTo(xk,36.2,1.98e-3); //C02
1548          xk -=0.01;
1549     track->PropagateTo(xk, 24.01, 2.7);  //Al    
1550          xk -=0.01;
1551     track->PropagateTo(xk, 44.77, 1.71);         //Tedlar
1552          xk -=0.04;
1553     track->PropagateTo(xk, 44.86, 1.45);         //Kevlar
1554          xk -=0.5;
1555     track->PropagateTo(xk, 41.28, 0.029);        //Nomex                                                 
1556                      
1557        ///////////////////////////////////////////////////////////////                   
1558   
1559    ///////////////////////////////////////////////////////////////
1560     AliITStrack trackITS(*track);
1561     AliITStrack result(*track);
1562     AliITStrack primarytrack(*track); 
1563     
1564 ///////////////////////////////////////////////////////////////////////////////////////////////
1565          TVector Vgeant(3);
1566          Vgeant=result.GetVertex(); 
1567                           
1568   // Definition of Dv and Zv for vertex constraint      
1569      Double_t sigmaDv=0.0050;  Double_t sigmaZv=0.010;  
1570     //Double_t sigmaDv=0.0015;  Double_t sigmaZv=0.0015;                                  
1571         Double_t uniform= gRandom->Uniform();
1572         Double_t signdv;
1573         if(uniform<=0.5) signdv=-1.;
1574            else
1575                  signdv=1.;
1576          
1577         Double_t Vr=TMath::Sqrt(Vgeant(0)*Vgeant(0)+ Vgeant(1)*Vgeant(1));
1578           Double_t Dv=gRandom->Gaus(signdv*Vr,(Float_t)sigmaDv); 
1579     Double_t Zv=gRandom->Gaus(Vgeant(2),(Float_t)sigmaZv);
1580                                 
1581   //cout<<" Dv e Zv = "<<Dv<<" "<<Zv<<"\n";                             
1582     trackITS.SetDv(Dv);  trackITS.SetZv(Zv);
1583     trackITS.SetsigmaDv(sigmaDv); trackITS.SetsigmaZv(sigmaZv); 
1584     result.SetDv(Dv);  result.SetZv(Zv);
1585     result.SetsigmaDv(sigmaDv); result.SetsigmaZv(sigmaZv);
1586     primarytrack.SetDv(Dv);  primarytrack.SetZv(Zv);
1587     primarytrack.SetsigmaDv(sigmaDv); primarytrack.SetsigmaZv(sigmaZv);                                                                 
1588
1589 /////////////////////////////////////////////////////////////////////////////////////////////////        
1590                 
1591     primarytrack.PrimaryTrack(rl);
1592     TVector  d2=primarytrack.Getd2();
1593     TVector  tgl2=primarytrack.Gettgl2();
1594     TVector  dtgl=primarytrack.Getdtgl();
1595     trackITS.Setd2(d2); trackITS.Settgl2(tgl2);  trackITS.Setdtgl(dtgl); 
1596     result.Setd2(d2); result.Settgl2(tgl2);  result.Setdtgl(dtgl);         
1597          /*                      
1598     trackITS.SetVertex(vertex); trackITS.SetErrorVertex(ervertex);
1599     result.SetVertex(vertex);   result.SetErrorVertex(ervertex);   
1600     */                         
1601     Tracking(trackITS,&result,recPoints,vettid, flagvert,rl);  
1602              
1603     // cout<<" progressive track number = "<<j<<"\r";
1604    // cout<<j<<"\r";
1605    // cout<<" progressive track number = "<<j<<"\n";
1606     Long_t labITS=result.GetLabel();
1607    // cout << " ITS track label = " << labITS << "\n";              
1608     int lab=track->GetLabel();              
1609    // cout << " TPC track label = " << lab <<"\n";
1610     //result.GetClusters(); getchar();  //to print the cluster matrix
1611          
1612 //propagation to vertex
1613         
1614     Double_t rbeam=3.;
1615      
1616     result.Propagation(rbeam);
1617        
1618     TMatrix *cov;
1619          cov=&result.GetCMatrix();       
1620     Double_t pt=TMath::Abs(result.GetPt());
1621     Double_t Dr=result.GetD();
1622     Double_t Z=result.GetZ();
1623     Double_t tgl=result.GetTgl();
1624     Double_t C=result.GetC();
1625     Double_t Cy=C/2.;
1626     Double_t Dz=Z-(tgl/Cy)*TMath::ASin(result.arga(rbeam));
1627           Dz-=Vgeant(2);
1628           
1629          // cout<<" Dr e dz alla fine = "<<Dr<<" "<<Dz<<"\n"; getchar();
1630     Double_t phi=result.Getphi();
1631     Double_t phivertex = phi - TMath::ASin(result.argA(rbeam));
1632     Double_t duepi=2.*TMath::Pi();       
1633     if(phivertex>duepi) phivertex-=duepi;
1634     if(phivertex<0.) phivertex+=duepi;
1635     Double_t Dtot=TMath::Sqrt(Dr*Dr+Dz*Dz);
1636          
1637 //////////////////////////////////////////////////////////////////////////////////////////      
1638     Int_t NumofCluster, idmodule,idpoint;
1639     NumofCluster=result.GetNumClust();
1640     if(NumofCluster >=5)  {      
1641
1642
1643       AliITSiotrack outtrack;
1644
1645       iotrack=&outtrack;
1646
1647       iotrack->SetStatePhi(phi);
1648       iotrack->SetStateZ(Z);
1649       iotrack->SetStateD(Dr);
1650       iotrack->SetStateTgl(tgl);
1651       iotrack->SetStateC(C);
1652                 Double_t radius=result.Getrtrack();
1653                 iotrack->SetRadius(radius);
1654                 Int_t charge;
1655                 if(C>0.) charge=-1;  else charge=1;
1656                 iotrack->SetCharge(charge);
1657                 
1658                 
1659
1660       iotrack->SetCovMatrix(cov);         
1661
1662       Double_t px=pt*TMath::Cos(phi);
1663       Double_t py=pt*TMath::Sin(phi);
1664       Double_t pz=pt*tgl;
1665                 
1666       Double_t xtrack=Dr*TMath::Sin(phi);
1667       Double_t ytrack=Dr*TMath::Cos(phi);
1668       Double_t ztrack=Dz+Vgeant(2);
1669
1670
1671       iotrack->SetPx(px);
1672       iotrack->SetPy(py);
1673       iotrack->SetPz(pz);
1674       iotrack->SetX(xtrack);
1675       iotrack->SetY(ytrack);
1676       iotrack->SetZ(ztrack);
1677       iotrack->SetLabel(labITS);
1678                 
1679       Int_t il;         
1680                 for(il=0;il<6; il++){
1681                   iotrack->SetIdPoint(il,result.GetIdPoint(il));
1682                   iotrack->SetIdModule(il,result.GetIdModule(il));              
1683                 }
1684       tracktree1.Fill();
1685
1686    //cout<<" labITS = "<<labITS<<"\n";
1687         //cout<<" phi z Dr tgl C = "<<phi<<" "<<Z<<" "<<Dr<<" "<<tgl<<" "<<C<<"\n";  getchar();    
1688
1689      DataOut(kkk) = ptg; kkk++; DataOut(kkk)=labITS; kkk++; DataOut(kkk)=lab; kkk++;            
1690
1691       for (il=0;il<6;il++) {
1692         idpoint=result.GetIdPoint(il);
1693         idmodule=result.GetIdModule(il);
1694         *(vettid[idmodule]+idpoint)=1; 
1695         iotrack->SetIdPoint(il,idpoint);
1696         iotrack->SetIdModule(il,idmodule);
1697       }
1698       
1699     //  cout<<"  +++++++++++++  pt e ptg = "<<pt<<" "<<ptg<<"  ++++++++++\n";
1700       Double_t difpt= (pt-ptg)/ptg*100.;                                        
1701       DataOut(kkk)=difpt; kkk++;                                             
1702       Double_t lambdag=TMath::ATan(pzg/ptg);
1703       Double_t   lam=TMath::ATan(tgl);      
1704       Double_t diflam = (lam - lambdag)*1000.;
1705       DataOut(kkk) = diflam; kkk++;                                         
1706       Double_t phig=TMath::ATan(pyg/pxg);
1707       Double_t phi=phivertex;  
1708       Double_t difphi = (phi - phig)*1000.;
1709       DataOut(kkk)=difphi; kkk++;
1710       DataOut(kkk)=Dtot*1.e4; kkk++;
1711       DataOut(kkk)=Dr*1.e4; kkk++;
1712       DataOut(kkk)=Dz*1.e4; kkk++;
1713       Int_t r;
1714       for (r=0; r<9; r++) { out<<DataOut(r)<<" ";}
1715       out<<"\n";
1716       kkk=0;  
1717                 
1718             
1719     } // end if on NumofCluster
1720   //gObjectTable->Print();    // stampa memoria     
1721   }  //  end for (int j=min_t; j<=max_t; j++)
1722   
1723   out.close();  
1724   
1725  
1726   static Bool_t first=kTRUE;
1727   static TFile *tfile;
1728
1729         if(first) {
1730             tfile=new TFile("itstracks.root","RECREATE");
1731             //cout<<"I have opened itstracks.root file "<<endl;
1732         }           
1733         first=kFALSE;
1734         tfile->cd();
1735         tfile->ls();
1736
1737    char hname[30];
1738    sprintf(hname,"TreeT%d",evNumber);
1739
1740   tracktree1.Write(hname);
1741
1742
1743   
1744             TTree *fAli=gAlice->TreeK();
1745             TFile *fileAli=0;
1746             
1747             if (fAli) fileAli =fAli->GetCurrentFile();
1748             fileAli->cd();
1749      
1750   ////////////////////////////////////////////////////////////////////////////////////////////////
1751
1752   printf("delete vectors\n");
1753   if(np) delete [] np;
1754   if(vettid) delete [] vettid;
1755   
1756 }