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