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