1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.39 2001/03/12 17:45:32 hristov
19 Changes needed on Sun with CC 5.0
21 Revision 1.38 2001/03/07 14:04:51 barbera
22 Some vector dimensions increased to cope with full events
24 Revision 1.37 2001/03/07 12:36:35 barbera
25 A change added in the tracking part to manage delta rays
27 Revision 1.36 2001/03/02 19:44:11 barbera
28 modified to taking into account new version tracking v1
30 Revision 1.35 2001/02/28 18:16:46 mariana
31 Make the code compatible with the new AliRun
33 Revision 1.34 2001/02/11 15:51:39 mariana
34 Set protection in MakeBranch
36 Revision 1.33 2001/02/10 22:26:39 mariana
37 Move the initialization of the containers for raw clusters in MakeTreeC()
39 Revision 1.32 2001/02/08 23:55:31 nilsen
40 Removed fMajor/MinorVersion variables in favor of variables in derived classes.
41 Set arrays char *det[3] = {"SPD","SDD","SSD"} as const.
43 Revision 1.31 2001/02/02 23:57:28 nilsen
44 Added include file that are no londer included in AliITSgeom.h
46 Revision 1.30 2001/01/30 09:23:13 hristov
47 Streamers removed (R.Brun)
49 Revision 1.29 2001/01/26 20:01:09 hristov
50 Major upgrade of AliRoot code
52 Revision 1.28 2000/12/18 14:02:00 barbera
53 new version of the ITS tracking to take into account the new TPC track parametrization
55 Revision 1.27 2000/12/08 13:49:27 barbera
56 Hidden declaration in a for loop removed to be compliant with HP-UX compiler
58 Revision 1.26 2000/11/27 13:12:13 barbera
59 New version containing the files for tracking
61 Revision 1.25 2000/11/12 22:38:05 barbera
62 Added header file for the SPD Bari model
64 Revision 1.24 2000/10/09 22:18:12 barbera
65 Bug fixes from MAriana to le AliITStest.C run correctly
67 Revision 1.23 2000/10/05 20:47:42 nilsen
68 fixed dependencies of include files. Tryed but failed to get a root automaticly
69 generates streamer function to work. Modified SetDefaults.
71 Revision 1.9.2.15 2000/10/04 16:56:40 nilsen
72 Needed to include stdlib.h
75 Revision 1.22 2000/10/04 19:45:52 barbera
76 Corrected by F. Carminati for v3.04
78 Revision 1.21 2000/10/02 21:28:08 fca
79 Removal of useless dependecies via forward declarations
81 Revision 1.20 2000/10/02 16:31:39 barbera
84 Revision 1.9.2.14 2000/10/02 15:43:51 barbera
85 General code clean-up (e.g., printf -> cout)
87 Revision 1.19 2000/09/22 12:13:25 nilsen
88 Patches and updates for fixes to this and other routines.
90 Revision 1.18 2000/07/12 05:32:20 fca
91 Correcting several syntax problem with static members
93 Revision 1.17 2000/07/10 16:07:18 fca
94 Release version of ITS code
96 Revision 1.9.2.3 2000/02/02 13:42:09 barbera
97 fixed AliITS.cxx for new AliRun structure. Added ITS hits list to list of hits which will have their track numbers updated
99 Revision 1.9.2.2 2000/01/23 03:03:13 nilsen
100 //fixed FillModule. Removed fi(fabs(xl)<dx....
102 Revision 1.9.2.1 2000/01/12 19:03:32 nilsen
103 This is the version of the files after the merging done in December 1999.
104 See the ReadMe110100.txt file for details
106 Revision 1.9 1999/11/14 14:33:25 fca
107 Correct problems with distructors and pointers, thanks to I.Hrivnacova
109 Revision 1.8 1999/09/29 09:24:19 fca
110 Introduction of the Copyright and cvs Log
114 ///////////////////////////////////////////////////////////////////////////////
116 // An overview of the basic philosophy of the ITS code development
117 // and analysis is show in the figure below.
120 <img src="picts/ITS/ITS_Analysis_schema.gif">
123 <font size=+2 color=red>
124 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
125 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
131 // AliITS. Inner Traking System base class.
132 // This class contains the base procedures for the Inner Tracking System
136 <img src="picts/ITS/AliITS_Class_Diagram.gif">
139 <font size=+2 color=red>
140 <p>This show the class diagram of the different elements that are part of
148 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
151 // Modified and documented by Bjorn S. Nilsen
155 // Modified and documented by A. Bologna
158 // AliITS is the general base class for the ITS. Also see AliDetector for
159 // futher information.
161 ///////////////////////////////////////////////////////////////////////////////
162 #include <iostream.h>
170 #include <TClonesArray.h>
172 #include <TObjectTable.h>
176 #include <TParticle.h>
181 #include "AliITSMap.h"
182 #include "AliITSDetType.h"
183 #include "AliITSClusterFinder.h"
184 //#include "AliITSsimulation.h"
185 #include "AliITSsimulationSPD.h"
186 #include "AliITSsimulationSDD.h"
187 #include "AliITSsimulationSSD.h"
188 #include "AliITSresponse.h"
189 #include "AliITSsegmentationSPD.h"
190 #include "AliITSresponseSPD.h"
191 #include "AliITSresponseSPDbari.h"
192 #include "AliITSsegmentationSDD.h"
193 #include "AliITSresponseSDD.h"
194 #include "AliITSsegmentationSSD.h"
195 #include "AliITSresponseSSD.h"
196 #include "AliITShit.h"
197 #include "AliITSgeom.h"
198 #include "AliITSdigit.h"
199 #include "AliITSmodule.h"
200 #include "AliITSRecPoint.h"
201 #include "AliITSRawCluster.h"
205 #include "AliITStrack.h"
206 #include "AliITSiotrack.h"
207 #include "AliITStracking.h"
208 #include "AliITSRad.h"
209 #include "../TPC/AliTPC.h"
210 #include "../TPC/AliTPCParam.h"
215 //_____________________________________________________________________________
216 AliITS::AliITS() : AliDetector() {
218 // Default initialiser for ITS
219 // The default constructor of the AliITS class. In addition to
220 // creating the AliITS class it zeros the variables fIshunt (a member
221 // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
222 // fITSpoints, fIdSens, and fIdName. The AliDetector default constructor
230 fNDetTypes = kNTYPES;
249 //_____________________________________________________________________________
250 AliITS::AliITS(const char *name, const char *title):AliDetector(name,title){
252 // Default initialiser for ITS
253 // The constructor of the AliITS class. In addition to creating the
254 // AliITS class, it allocates memory for the TClonesArrays fHits and
255 // fDigits, and for the TObjArray fITSpoints. It also zeros the variables
256 // fIshunt (a member of AliDetector class), fEuclidOut, and fIdN, and zeros
257 // the pointers fIdSens and fIdName. To help in displaying hits via the ROOT
258 // macro display.C AliITS also sets the marker color to red. The variables
259 // passes with this constructor, const char *name and *title, are used by
260 // the constructor of AliDetector class. See AliDetector class for a
261 // description of these parameters and its constructor functions.
265 fHits = new TClonesArray("AliITShit", 1560);
266 gAlice->AddHitList(fHits);
268 fNDetTypes = kNTYPES;
270 fNdtype = new Int_t[kNTYPES];
271 fDtype = new TObjArray(kNTYPES);
273 fNctype = new Int_t[kNTYPES];
274 fCtype = new TObjArray(kNTYPES);
290 fDetTypes = new TObjArray(kNTYPES);
293 for(i=0;i<kNTYPES;i++) {
294 (*fDetTypes)[i]=new AliITSDetType();
300 SetMarkerColor(kRed);
304 //___________________________________________________________________________
305 AliITS::AliITS(AliITS &source){
307 if(this==&source) return;
308 Error("AliITS::Copy constructor",
309 "You are not allowed to make a copy of the AliITS");
312 //____________________________________________________________________________
313 AliITS& AliITS::operator=(AliITS &source){
314 // assignment operator
315 if(this==&source) return *this;
316 Error("AliITS::operator=",
317 "You are not allowed to make a copy of the AliITS");
319 return *this; //fake return
321 //____________________________________________________________________________
322 void AliITS::ClearModules(){
323 //clear the modules TObjArray
325 if(fITSmodules) fITSmodules->Delete();
328 //_____________________________________________________________________________
331 // Default distructor for ITS
332 // The default destructor of the AliITS class. In addition to deleting
333 // the AliITS class it deletes the memory pointed to by the fHits, fDigits,
334 // fIdSens, fIdName, and fITSpoints.
341 // delete fIdName; // TObjArray of TObjStrings
342 if(fIdName!=0) delete[] fIdName; // Array of TStrings
343 if(fIdSens!=0) delete[] fIdSens;
345 this->ClearModules();
347 }// end if fITSmodules!=0
367 if (fTreeC) delete fTreeC;
369 if (fITSgeom) delete fITSgeom;
373 //___________________________________________
374 AliITSDetType* AliITS::DetType(Int_t id)
376 //return pointer to id detector type
377 return ((AliITSDetType*) (*fDetTypes)[id]);
380 //___________________________________________
381 void AliITS::SetClasses(Int_t id, const char *digit, const char *cluster)
383 //set the digit and cluster classes to be used for the id detector type
384 ((AliITSDetType*) (*fDetTypes)[id])->ClassNames(digit,cluster);
387 //___________________________________________
388 void AliITS::SetResponseModel(Int_t id, AliITSresponse *response)
390 //set the response model for the id detector type
392 ((AliITSDetType*) (*fDetTypes)[id])->ResponseModel(response);
396 //___________________________________________
397 void AliITS::SetSegmentationModel(Int_t id, AliITSsegmentation *seg)
399 //set the segmentation model for the id detector type
401 ((AliITSDetType*) (*fDetTypes)[id])->SegmentationModel(seg);
405 //___________________________________________
406 void AliITS::SetSimulationModel(Int_t id, AliITSsimulation *sim)
408 //set the simulation model for the id detector type
410 ((AliITSDetType*) (*fDetTypes)[id])->SimulationModel(sim);
413 //___________________________________________
414 void AliITS::SetReconstructionModel(Int_t id, AliITSClusterFinder *reconst)
416 //set the cluster finder model for the id detector type
418 ((AliITSDetType*) (*fDetTypes)[id])->ReconstructionModel(reconst);
422 //_____________________________________________________________________________
423 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
426 // The function to add information to the AliITShit class. See the
427 // AliITShit class for a full description. This function allocates the
428 // necessary new space for the hit information and passes the variable
429 // track, and the pointers *vol and *hits to the AliITShit constructor
432 TClonesArray &lhits = *fHits;
433 new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
435 //_____________________________________________________________________________
436 void AliITS::AddRealDigit(Int_t id, Int_t *digits)
438 // add a real digit - as coming from data
440 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
441 new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
444 //_____________________________________________________________________________
445 void AliITS::AddSimDigit(Int_t id, AliITSdigit *d)
448 // add a simulated digit
450 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
455 new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
458 new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
461 new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
467 //_____________________________________________________________________________
468 void AliITS::AddSimDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,Int_t *hits,Float_t *charges){
470 // add a simulated digit to the list
472 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
476 new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks,hits);
479 new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,hits,charges);
482 new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks,hits);
488 //_____________________________________________________________________________
489 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c)
492 // add a cluster to the list
494 TClonesArray &lcl = *((TClonesArray*)(*fCtype)[id]);
499 new(lcl[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
502 new(lcl[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
505 new(lcl[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
512 //_____________________________________________________________________________
513 void AliITS::AddRecPoint(const AliITSRecPoint &r)
516 // Add a reconstructed space point to the list
518 TClonesArray &lrecp = *fRecPoints;
519 new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
523 //____________________________________________
524 void AliITS::ResetDigits()
527 // Reset number of digits and the digits array for the ITS detector
533 for (i=0;i<kNTYPES;i++ ) {
534 if ((*fDtype)[i]) ((TClonesArray*)(*fDtype)[i])->Clear();
535 if (fNdtype) fNdtype[i]=0;
539 //____________________________________________
540 void AliITS::ResetDigits(Int_t i)
543 // Reset number of digits and the digits array for this branch
545 if ((*fDtype)[i]) ((TClonesArray*)(*fDtype)[i])->Clear();
546 if (fNdtype) fNdtype[i]=0;
550 //____________________________________________
551 void AliITS::ResetClusters()
554 // Reset number of clusters and the clusters array for ITS
558 for (i=0;i<kNTYPES;i++ ) {
559 if ((*fCtype)[i]) ((TClonesArray*)(*fCtype)[i])->Clear();
560 if (fNctype) fNctype[i]=0;
565 //____________________________________________
566 void AliITS::ResetClusters(Int_t i)
569 // Reset number of clusters and the clusters array for this branch
571 if ((*fCtype)[i]) ((TClonesArray*)(*fCtype)[i])->Clear();
572 if (fNctype) fNctype[i]=0;
577 //____________________________________________
578 void AliITS::ResetRecPoints()
581 // Reset number of rec points and the rec points array
583 if (fRecPoints) fRecPoints->Clear();
588 //_____________________________________________________________________________
589 Int_t AliITS::DistancetoPrimitive(Int_t , Int_t ){
591 // Distance from mouse to ITS on the screen. Dummy routine
592 // A dummy routine used by the ROOT macro display.C to allow for the
593 // use of the mouse (pointing device) in the macro. In general this should
594 // never be called. If it is it returns the number 9999 for any value of
600 //_____________________________________________________________________________
603 // Initialise ITS after it has been built
604 // This routine initializes the AliITS class. It is intended to be called
605 // from the Init function in AliITSv?. Besides displaying a banner
606 // indicating that it has been called it initializes the array fIdSens
607 // and sets the default segmentation, response, digit and raw cluster classes
608 // Therefore it should be called after a call to CreateGeometry.
615 for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
619 //_____________________________________________________________________________
620 void AliITS::SetDefaults()
622 // sets the default segmentation, response, digit and raw cluster classes
624 printf("SetDefaults\n");
626 AliITSDetType *iDetType;
632 if (!iDetType->GetSegmentationModel()) {
633 AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(fITSgeom);
634 SetSegmentationModel(0,seg0);
636 if (!iDetType->GetResponseModel()) {
637 SetResponseModel(0,new AliITSresponseSPD());
639 // set digit and raw cluster classes to be used
641 const char *kData0=(iDetType->GetResponseModel())->DataType();
642 if (strstr(kData0,"real")) {
643 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSPD");
644 } else iDetType->ClassNames("AliITSdigitSPD","AliITSRawClusterSPD");
648 if (!iDetType->GetResponseModel()) {
649 SetResponseModel(1,new AliITSresponseSDD());
651 AliITSresponse *resp1=iDetType->GetResponseModel();
652 if (!iDetType->GetSegmentationModel()) {
653 AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(fITSgeom,resp1);
654 SetSegmentationModel(1,seg1);
656 const char *kData1=(iDetType->GetResponseModel())->DataType();
657 const char *kopt=iDetType->GetResponseModel()->ZeroSuppOption();
658 if ((!strstr(kopt,"2D")) && (!strstr(kopt,"1D")) || strstr(kData1,"real") ) {
659 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
660 } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
664 if (!iDetType->GetSegmentationModel()) {
665 AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(fITSgeom);
666 SetSegmentationModel(2,seg2);
668 if (!iDetType->GetResponseModel()) {
669 SetResponseModel(2,new AliITSresponseSSD());
671 const char *kData2=(iDetType->GetResponseModel())->DataType();
672 if (strstr(kData2,"real")) {
673 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSSD");
674 } else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
677 Warning("SetDefaults","Only the three basic detector types are initialised!");
684 //_____________________________________________________________________________
685 void AliITS::SetDefaultSimulation()
690 //_____________________________________________________________________________
691 void AliITS::SetDefaultClusterFinders()
696 //_____________________________________________________________________________
698 void AliITS::MakeTreeC(Option_t *option)
700 // create a separate tree to store the clusters
702 cout << "AliITS::MakeTreeC" << endl;
704 const char *optC = strstr(option,"C");
705 if (optC && !fTreeC) fTreeC = new TTree("TC","Clusters in ITS");
708 Int_t buffersize = 4000;
711 const char *det[3] = {"SPD","SDD","SSD"};
716 // one branch for Clusters per type of detector
718 for (i=0; i<kNTYPES ;i++) {
719 AliITSDetType *iDetType=DetType(i);
720 iDetType->GetClassNames(digclass,clclass);
722 (*fCtype)[i] = new TClonesArray(clclass,10000);
723 if (kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
724 else sprintf(branchname,"%sClusters%d",GetName(),i+1);
725 if (fCtype && fTreeC) {
726 TreeC()->Branch(branchname,&((*fCtype)[i]), buffersize);
727 cout << "Making Branch " << branchname;
728 cout << " for Clusters of detector type " << i+1 << endl;
734 //_____________________________________________________________________________
735 void AliITS::GetTreeC(Int_t event)
738 cout << "AliITS::GetTreeC" << endl;
740 // get the clusters tree for this event and set the branch address
744 const char *det[3] = {"SPD","SDD","SSD"};
751 sprintf(treeName,"TreeC%d",event);
752 fTreeC = (TTree*)gDirectory->Get(treeName);
757 for (i=0; i<kNTYPES; i++) {
758 if (kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
759 else sprintf(branchname,"%sClusters%d",GetName(),i+1);
761 branch = fTreeC->GetBranch(branchname);
762 if (branch) branch->SetAddress(&((*fCtype)[i]));
766 Error("AliITS::GetTreeC",
767 "cannot find Clusters Tree for event:%d\n",event);
771 //_____________________________________________________________________________
772 void AliITS::MakeBranch(Option_t* option, char *file)
775 // Creates Tree branches for the ITS.
778 Int_t buffersize = 4000;
780 sprintf(branchname,"%s",GetName());
782 AliDetector::MakeBranch(option,file);
784 const char *cD = strstr(option,"D");
785 const char *cR = strstr(option,"R");
789 // one branch for digits per type of detector
791 const char *det[3] = {"SPD","SDD","SSD"};
797 for (i=0; i<kNTYPES ;i++) {
798 AliITSDetType *iDetType=DetType(i);
799 iDetType->GetClassNames(digclass,clclass);
801 if(!((*fDtype)[i])) (*fDtype)[i] = new TClonesArray(digclass,10000);
805 for (i=0; i<kNTYPES ;i++) {
806 if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
807 else sprintf(branchname,"%sDigits%d",GetName(),i+1);
808 if (fDtype && gAlice->TreeD()) {
809 gAlice->MakeBranchInTree(gAlice->TreeD(),
810 branchname, &((*fDtype)[i]), buffersize, file);
811 cout << "Making Branch " << branchname;
812 cout << " for digits of type "<< i+1 << endl;
819 // only one branch for rec points for all detector types
821 sprintf(branchname,"%sRecPoints",GetName());
823 if(!fRecPoints) fRecPoints=new TClonesArray("AliITSRecPoint",10000);
825 if (fRecPoints && gAlice->TreeR()) {
826 gAlice->MakeBranchInTree(gAlice->TreeR(),
827 branchname, &fRecPoints, buffersize, file) ;
828 cout << "Making Branch " << branchname;
829 cout << " for reconstructed space points" << endl;
834 //___________________________________________
835 void AliITS::SetTreeAddress()
838 // Set branch address for the Trees.
841 AliDetector::SetTreeAddress();
843 const char *det[3] = {"SPD","SDD","SSD"};
846 TTree *treeD = gAlice->TreeD();
847 TTree *treeR = gAlice->TreeR();
851 for (i=0; i<kNTYPES; i++) {
852 if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
853 else sprintf(branchname,"%sDigits%d",GetName(),i+1);
855 branch = treeD->GetBranch(branchname);
856 if (branch) branch->SetAddress(&((*fDtype)[i]));
863 sprintf(branchname,"%sRecPoints",GetName());
864 branch = treeR->GetBranch(branchname);
865 if (branch) branch->SetAddress(&fRecPoints);
871 //____________________________________________________________________________
872 void AliITS::InitModules(Int_t size,Int_t &nmodules){
874 //initialize the modules array
877 fITSmodules->Delete();
881 Int_t nl,indexMAX,index;
883 if(size<=0){ // default to using data stored in AliITSgeom
885 Error("AliITS::InitModules",
886 "in AliITS::InitModule fITSgeom not defined\n");
888 } // end if fITSgeom==0
889 nl = fITSgeom->GetNlayers();
890 indexMAX = fITSgeom->GetModuleIndex(nl,fITSgeom->GetNladders(nl),
891 fITSgeom->GetNdetectors(nl))+1;
893 fITSmodules = new TObjArray(indexMAX);
894 for(index=0;index<indexMAX;index++){
895 fITSmodules->AddAt( new AliITSmodule(index),index);
898 fITSmodules = new TObjArray(size);
899 for(index=0;index<size;index++) {
900 fITSmodules->AddAt( new AliITSmodule(index),index);
907 //____________________________________________________________________________
908 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t nmodules,Option_t *option,Text_t *filename){
910 // fill the modules with the sorted by module hits; add hits from background
914 static TTree *trH1; //Tree with background hits
915 static TClonesArray *fHits2; //List of hits for one track only
917 static Bool_t first=kTRUE;
919 const char *addBgr = strstr(option,"Add");
924 cout<<"filename "<<filename<<endl;
925 file=new TFile(filename);
926 cout<<"I have opened "<<filename<<" file "<<endl;
927 fHits2 = new TClonesArray("AliITShit",1000 );
932 // Get Hits Tree header from file
933 if(fHits2) fHits2->Clear();
934 if(trH1) delete trH1;
938 sprintf(treeName,"TreeH%d",bgrev);
939 trH1 = (TTree*)gDirectory->Get(treeName);
940 //printf("TrH1 %p of treename %s for event %d \n",trH1,treeName,bgrev);
943 Error("AliITS::FillModules",
944 "cannot find Hits Tree for event:%d\n",bgrev);
946 // Set branch addresses
949 sprintf(branchname,"%s",GetName());
950 if (trH1 && fHits2) {
951 branch = trH1->GetBranch(branchname);
952 if (branch) branch->SetAddress(&fHits2);
956 //Int_t ntracks1 =(Int_t)TrH1->GetEntries();
957 //printf("background - ntracks1 - %d\n",ntracks1);
960 //Int_t npart = gAlice->GetEvent(evnt);
961 //if(npart<=0) return;
962 TClonesArray *itsHits = this->Hits();
963 Int_t lay,lad,det,index;
967 TTree *iTH = gAlice->TreeH();
968 Int_t ntracks =(Int_t) iTH->GetEntries();
971 for(t=0; t<ntracks; t++){
974 Int_t nhits = itsHits->GetEntriesFast();
975 //printf("nhits %d\n",nhits);
976 if (!nhits) continue;
977 for(h=0; h<nhits; h++){
978 itsHit = (AliITShit *)itsHits->UncheckedAt(h);
979 itsHit->GetDetectorID(lay,lad,det);
980 // temporarily index=det-1 !!!
981 if(fITSgeom) index = fITSgeom->GetModuleIndex(lay,lad,det);
984 mod = this->GetModule(index);
985 mod->AddHit(itsHit,t,h);
986 } // end loop over hits
987 } // end loop over tracks
989 // open the file with background
993 ntracks =(Int_t)trH1->GetEntries();
994 //printf("background - ntracks1 %d\n",ntracks);
995 //printf("background - Start loop over tracks \n");
998 for (track=0; track<ntracks; track++) {
1000 if (fHits2) fHits2->Clear();
1001 trH1->GetEvent(track);
1003 for(i=0;i<fHits2->GetEntriesFast();++i) {
1005 itsHit=(AliITShit*) (*fHits2)[i];
1006 itsHit->GetDetectorID(lay,lad,det);
1007 // temporarily index=det-1 !!!
1008 if(fITSgeom) index = fITSgeom->GetModuleIndex(lay,lad,det);
1011 mod = this->GetModule(index);
1012 mod->AddHit(itsHit,track,i);
1013 } // end loop over hits
1014 } // end loop over tracks
1016 TTree *fAli=gAlice->TreeK();
1019 if (fAli) fileAli =fAli->GetCurrentFile();
1024 //gObjectTable->Print();
1028 //____________________________________________________________________________
1030 void AliITS::SDigits2Digits()
1033 AliITSgeom *geom = GetITSgeom();
1036 AliITSDetType *iDetType;
1037 iDetType=DetType(0);
1038 AliITSsegmentationSPD *seg0=(AliITSsegmentationSPD*)iDetType->GetSegmentationModel();
1039 AliITSresponseSPD *res0 = (AliITSresponseSPD*)iDetType->GetResponseModel();
1040 AliITSsimulationSPD *sim0=new AliITSsimulationSPD(seg0,res0);
1041 SetSimulationModel(0,sim0);
1043 // printf("SPD dimensions %f %f \n",seg0->Dx(),seg0->Dz());
1044 // printf("SPD npixels %d %d \n",seg0->Npz(),seg0->Npx());
1045 // printf("SPD pitches %d %d \n",seg0->Dpz(0),seg0->Dpx(0));
1049 //Set response functions
1050 // SDD compression param: 2 fDecrease, 2fTmin, 2fTmax or disable, 2 fTolerance
1052 iDetType=DetType(1);
1053 AliITSresponseSDD *res1 = (AliITSresponseSDD*)iDetType->GetResponseModel();
1055 res1=new AliITSresponseSDD();
1056 SetResponseModel(1,res1);
1058 Float_t noise, baseline;
1059 res1->GetNoiseParam(noise,baseline);
1060 Float_t noise_after_el = res1->GetNoiseAfterElectronics();
1061 Float_t fCutAmp = baseline + 2.*noise_after_el;
1062 Int_t cp[8]={0,0,(int)fCutAmp,(int)fCutAmp,0,0,0,0}; //1D
1063 res1->SetCompressParam(cp);
1064 AliITSsegmentationSDD *seg1=(AliITSsegmentationSDD*)iDetType->GetSegmentationModel();
1066 seg1 = new AliITSsegmentationSDD(geom,res1);
1067 SetSegmentationModel(1,seg1);
1069 AliITSsimulationSDD *sim1=new AliITSsimulationSDD(seg1,res1);
1070 SetSimulationModel(1,sim1);
1073 iDetType=DetType(2);
1074 AliITSsegmentationSSD *seg2=(AliITSsegmentationSSD*)iDetType->GetSegmentationModel();
1075 AliITSresponseSSD *res2 = (AliITSresponseSSD*)iDetType->GetResponseModel();
1076 res2->SetSigmaSpread(3.,2.);
1077 AliITSsimulationSSD *sim2=new AliITSsimulationSSD(seg2,res2);
1078 SetSimulationModel(2,sim2);
1080 cerr<<"Digitizing ITS...\n";
1084 HitsToDigits(0,0,-1," ","All"," ");
1085 timer.Stop(); timer.Print();
1093 //____________________________________________________________________________
1094 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size, Option_t *option, Option_t *opt,Text_t *filename)
1096 // keep galice.root for signal and name differently the file for
1097 // background when add! otherwise the track info for signal will be lost !
1099 // the condition below will disappear when the geom class will be
1100 // initialised for all versions - for the moment it is only for v5 !
1101 // 7 is the SDD beam test version
1102 Int_t ver = this->IsVersion();
1103 if(ver!=5 && ver!=7) return;
1105 const char *all = strstr(opt,"All");
1106 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1109 InitModules(size,nmodules);
1110 FillModules(evNumber,bgrev,nmodules,option,filename);
1113 AliITSsimulation* sim;
1114 //TObjArray *branches=gAlice->TreeD()->GetListOfBranches();
1115 AliITSgeom *geom = GetITSgeom();
1119 for (id=0;id<kNTYPES;id++) {
1120 if (!all && !det[id]) continue;
1121 //branch = (TBranch*)branches->UncheckedAt(id);
1122 AliITSDetType *iDetType=DetType(id);
1123 sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1125 Error("HitsToDigits","The simulation class was not instantiated!");
1127 // or SetDefaultSimulation();
1130 first = geom->GetStartDet(id);
1131 last = geom->GetLastDet(id);
1132 } else first=last=0;
1133 cout << "det type " << id << " first, last "<< first << last << endl;
1134 for(module=first;module<=last;module++) {
1135 AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
1136 sim->DigitiseModule(mod,module,evNumber);
1137 // fills all branches - wasted disk space
1138 gAlice->TreeD()->Fill();
1140 // try and fill only the branch
1143 } // loop over modules
1144 } // loop over detector types
1148 Int_t nentries=(Int_t)gAlice->TreeD()->GetEntries();
1149 cout << "nentries in TreeD" << nentries << endl;
1152 sprintf(hname,"TreeD%d",evNumber);
1153 gAlice->TreeD()->Write(hname,TObject::kOverwrite);
1155 gAlice->TreeD()->Reset();
1160 //____________________________________________________________________________
1161 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt)
1163 // cluster finding and reconstruction of space points
1165 // the condition below will disappear when the geom class will be
1166 // initialised for all versions - for the moment it is only for v5 !
1167 // 7 is the SDD beam test version
1168 Int_t ver = this->IsVersion();
1171 const char *all = strstr(opt,"All");
1172 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1174 static Bool_t first=kTRUE;
1175 if (!TreeC() && first) {
1180 TTree *treeC=TreeC();
1184 AliITSClusterFinder* rec;
1186 //TObjArray *branches=gAlice->TreeR()->GetListOfBranches();
1187 AliITSgeom *geom = GetITSgeom();
1190 for (id=0;id<kNTYPES;id++) {
1191 if (!all && !det[id]) continue;
1192 //branch = (TBranch*)branches->UncheckedAt(id);
1193 AliITSDetType *iDetType=DetType(id);
1194 rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1196 Error("DigitsToRecPoints","The cluster finder class was not instantiated!");
1198 // or SetDefaultClusterFinders();
1200 TClonesArray *itsDigits = this->DigitsAddress(id);
1204 first = geom->GetStartDet(id);
1205 last = geom->GetLastDet(id);
1206 } else first=last=0;
1207 //printf("first last %d %d\n",first,last);
1208 for(module=first;module<=last;module++) {
1209 this->ResetDigits();
1210 if (all) gAlice->TreeD()->GetEvent(lastentry+module);
1211 else gAlice->TreeD()->GetEvent(lastentry+(module-first));
1212 Int_t ndigits = itsDigits->GetEntriesFast();
1213 if (ndigits) rec->FindRawClusters();
1214 gAlice->TreeR()->Fill();
1218 // try and fill only the branch
1220 //ResetRecPoints(id);
1221 } // loop over modules
1222 } // loop over detector types
1225 Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1226 Int_t ncentries=(Int_t)treeC->GetEntries();
1227 cout << " nentries ncentries " << nentries << ncentries << endl;
1230 sprintf(hname,"TreeR%d",evNumber);
1231 gAlice->TreeR()->Write(hname,TObject::kOverwrite);
1233 gAlice->TreeR()->Reset();
1235 sprintf(hname,"TreeC%d",evNumber);
1236 treeC->Write(hname,TObject::kOverwrite);
1241 //____________________________________________________________________________
1242 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
1243 Option_t *option,Option_t *opt,Text_t *filename)
1245 // keep galice.root for signal and name differently the file for
1246 // background when add! otherwise the track info for signal will be lost !
1249 // the condition below will disappear when the geom class will be
1250 // initialised for all versions - for the moment it is only for v5 !
1251 Int_t ver = this->IsVersion();
1254 const char *all = strstr(opt,"All");
1255 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1258 InitModules(size,nmodules);
1259 FillModules(evNumber,bgrev,nmodules,option,filename);
1262 AliITSsimulation* sim;
1263 AliITSgeom *geom = GetITSgeom();
1265 TRandom *random=new TRandom[9];
1266 random[0].SetSeed(111);
1267 random[1].SetSeed(222);
1268 random[2].SetSeed(333);
1269 random[3].SetSeed(444);
1270 random[4].SetSeed(555);
1271 random[5].SetSeed(666);
1272 random[6].SetSeed(777);
1273 random[7].SetSeed(888);
1274 random[8].SetSeed(999);
1278 for (id=0;id<kNTYPES;id++) {
1279 if (!all && !det[id]) continue;
1280 AliITSDetType *iDetType=DetType(id);
1281 sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1283 Error("HitsToFastPoints","The simulation class was not instantiated!");
1285 // or SetDefaultSimulation();
1287 Int_t first = geom->GetStartDet(id);
1288 Int_t last = geom->GetLastDet(id);
1289 for(module=first;module<=last;module++) {
1290 AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
1291 sim->CreateFastRecPoints(mod,module,random);
1292 gAlice->TreeR()->Fill();
1294 } // loop over modules
1295 } // loop over detector types
1300 //Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1303 sprintf(hname,"TreeR%d",evNumber);
1304 gAlice->TreeR()->Write(hname,TObject::kOverwrite);
1306 gAlice->TreeR()->Reset();
1312 //________________________________________________________________
1314 AliITStrack AliITS::Tracking(AliITStrack &track, AliITStrack *reference,TObjArray *fastpoints, Int_t
1315 **vettid, Bool_t flagvert, AliITSRad *rl ) {
1317 //Origin A. Badala' and G.S. Pappalardo: e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1320 TList *list= new TList();
1322 AliITStrack tr(track);
1326 Double_t Pt=(tr).GetPt();
1327 //cout << "\n Pt = " << Pt <<"\n";
1329 AliITStracking obj(list, reference, this, fastpoints,TMath::Abs(Pt),vettid, flagvert, rl);
1334 TVector VecTotLabref(18);
1336 for(lay=5; lay>=0; lay--) {
1337 TVector VecLabref(3);
1338 VecLabref=(*reference).GetLabTrack(lay);
1339 Float_t ClustZ=(*reference).GetZclusterTrack( lay); //modified il 5-3-2001
1340 for(k=0; k<3; k++) { // {itot++; VecTotLabref(itot)=VecLabref(k);} // modified 5-3-2002
1342 Int_t lpp=(Int_t)VecLabref(k);
1344 TParticle *p=(TParticle*) gAlice->Particle(lpp);
1345 Int_t pcode=p->GetPdgCode();
1346 if(pcode==11) VecLabref(k)=p->GetFirstMother();
1349 itot++; VecTotLabref(itot)=VecLabref(k);
1350 if(VecLabref(k)==0. && ClustZ == 0.) VecTotLabref(itot) =-3.;
1355 (*reference).Search(VecTotLabref, labref, freq);
1357 if(freq < 5) labref=-labref;
1358 (*reference).SetLabel(labref);
1366 //________________________________________________________________
1370 void AliITS::DoTracking(Int_t evNumber, Int_t min_t, Int_t max_t, TFile *file, Bool_t flagvert) {
1372 // ex macro for tracking ITS
1374 printf("begin DoTracking - file %p\n",file);
1376 const char *pname="75x40_100x60";
1378 Int_t imax=200,jmax=450;
1379 AliITSRad *rl = new AliITSRad(imax,jmax);
1380 //cout<<" dopo costruttore AliITSRad\n"; getchar();
1384 Float_t px,py,pz,x,y,z,pxg,pyg,pzg,ptg;
1389 gAlice->GetEvent(0);
1391 AliTPC *TPC=(AliTPC*)gAlice->GetDetector("TPC");
1392 AliTPCParam *digp = (AliTPCParam*)file->Get(pname);
1393 if (digp!=0) TPC->SetParam(digp);
1395 GoodTrack gt[15000];
1397 ifstream in("itsgood_tracks");
1399 cerr<<"Reading itsgood tracks...\n";
1400 while (in>>gt[ngood].lab>>gt[ngood].code
1401 >>gt[ngood].px >>gt[ngood].py>>gt[ngood].pz
1402 >>gt[ngood].x >>gt[ngood].y >>gt[ngood].z
1403 >>gt[ngood].pxg >>gt[ngood].pyg >>gt[ngood].pzg
1404 >>gt[ngood].ptg >>gt[ngood].flag) {
1408 cerr<<"Too many good tracks !\n";
1412 if (!in.eof()) cerr<<"Read error (itsgood_tracks) !\n";
1416 TFile *tf=TFile::Open("tpctracks.root");
1417 if (!tf->IsOpen()) {cerr<<"Can't open tpctracks.root !\n"; return ;}
1418 TObjArray tracks(200000);
1419 TTree *tracktree=(TTree*)tf->Get("TreeT");
1420 TBranch *tbranch=tracktree->GetBranch("tracks");
1421 Int_t nentr=(Int_t)tracktree->GetEntries();
1423 for (kk=0; kk<nentr; kk++) {
1424 AliTPCtrack *iotrack=new AliTPCtrack;
1425 tbranch->SetAddress(&iotrack);
1426 tracktree->GetEvent(kk);
1427 tracks.AddLast(iotrack);
1432 Int_t nt = tracks.GetEntriesFast();
1433 cerr<<"Number of found tracks "<<nt<<endl;
1438 Double_t ptg=0.,pxg=0.,pyg=0.,pzg=0.;
1440 ////////////////////////////// good tracks definition in TPC ////////////////////////////////
1442 ofstream out1 ("AliITSTrag.out");
1444 for (i=0; i<ngood; i++) out1 << gt[i].ptg << "\n";
1449 TTree *TR=gAlice->TreeR();
1450 Int_t nent=(Int_t)TR->GetEntries();
1451 TClonesArray *recPoints = RecPoints();
1453 Int_t totalpoints=0;
1454 Int_t *np = new Int_t[nent];
1455 Int_t **vettid = new Int_t* [nent];
1457 for (mod=0; mod<nent; mod++) {
1459 this->ResetRecPoints();
1460 //gAlice->TreeR()->GetEvent(mod+1); //first entry in TreeR is empty
1461 gAlice->TreeR()->GetEvent(mod); //first entry in TreeR is empty
1462 numbpoints = recPoints->GetEntries();
1463 totalpoints+=numbpoints;
1464 np[mod] = numbpoints;
1465 //cout<<" mod = "<<mod<<" numbpoints = "<<numbpoints<<"\n"; getchar();
1466 vettid[mod] = new Int_t[numbpoints];
1468 for (ii=0;ii<numbpoints; ii++) *(vettid[mod]+ii)=0;
1474 if(min_t < 0) {min_t = 0; max_t = nt-1;}
1477 ///////////////////////////////// Definition of vertex end its error ////////////////////////////
1478 ////////////////////////// In the future it will be given by a method ///////////////////////////
1483 Float_t sigmavx=0.0050; // 50 microns
1484 Float_t sigmavy=0.0050; // 50 microns
1485 Float_t sigmavz=0.010; // 100 microns
1487 //Vx+=gRandom->Gaus(0,sigmavx); Vy+=gRandom->Gaus(0,sigmavy); Vz+=gRandom->Gaus(0,sigmavz);
1488 TVector vertex(3), ervertex(3)
1489 vertex(0)=Vx; vertex(1)=Vy; vertex(2)=Vz;
1490 ervertex(0)=sigmavx; ervertex(1)=sigmavy; ervertex(2)=sigmavz;
1491 /////////////////////////////////////////////////////////////////////////////////////////////////
1494 //TDirectory *savedir=gDirectory;
1496 TTree tracktree1("TreeT","Tree with ITS tracks");
1497 AliITSiotrack *iotrack=0;
1498 tracktree1.Branch("ITStracks","AliITSiotrack",&iotrack,32000,0);
1500 ofstream out ("AliITSTra.out");
1503 for (j=min_t; j<=max_t; j++) {
1504 track=(AliTPCtrack*)tracks.UncheckedAt(j);
1506 if (!track) continue;
1507 ////// elimination of not good tracks ////////////
1508 Int_t ilab=TMath::Abs(track->GetLabel());
1510 for (iii=0;iii<ngood;iii++) {
1511 //cout<<" ilab, gt[iii].lab = "<<ilab<<" "<<gt[iii].lab<<"\n"; getchar();
1512 if (ilab==gt[iii].lab) {
1521 //cout<<" j flaglab = " <<j<<" "<<flaglab<<"\n"; getchar();
1522 if (!flaglab) continue;
1523 //cout<<" j = " <<j<<"\n"; getchar();
1525 ////// old propagation to the end of TPC //////////////
1527 track->PropagateTo(xk);
1529 track->PropagateTo(xk,42.7,2.27); //C
1531 track->PropagateTo(xk,36.2,1.98e-3); //C02
1533 track->PropagateTo(xk,42.7,2.27); //C
1534 ///////////////////////////////////////////////////
1537 ////// new propagation to the end of TPC //////////////
1539 track->PropagateTo(xk, 28.94, 1.204e-3); //Ne
1541 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
1543 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
1545 track->PropagateTo(xk, 41.28, 0.029); //Nomex
1547 track->PropagateTo(xk,36.2,1.98e-3); //C02
1549 track->PropagateTo(xk, 24.01, 2.7); //Al
1551 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
1553 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
1555 track->PropagateTo(xk, 41.28, 0.029); //Nomex
1557 ///////////////////////////////////////////////////////////////
1559 ///////////////////////////////////////////////////////////////
1560 AliITStrack trackITS(*track);
1561 AliITStrack result(*track);
1562 AliITStrack primarytrack(*track);
1564 ///////////////////////////////////////////////////////////////////////////////////////////////
1566 Vgeant=result.GetVertex();
1568 // Definition of Dv and Zv for vertex constraint
1569 Double_t sigmaDv=0.0050; Double_t sigmaZv=0.010;
1570 //Double_t sigmaDv=0.0015; Double_t sigmaZv=0.0015;
1571 Double_t uniform= gRandom->Uniform();
1573 if(uniform<=0.5) signdv=-1.;
1577 Double_t Vr=TMath::Sqrt(Vgeant(0)*Vgeant(0)+ Vgeant(1)*Vgeant(1));
1578 Double_t Dv=gRandom->Gaus(signdv*Vr,(Float_t)sigmaDv);
1579 Double_t Zv=gRandom->Gaus(Vgeant(2),(Float_t)sigmaZv);
1581 //cout<<" Dv e Zv = "<<Dv<<" "<<Zv<<"\n";
1582 trackITS.SetDv(Dv); trackITS.SetZv(Zv);
1583 trackITS.SetsigmaDv(sigmaDv); trackITS.SetsigmaZv(sigmaZv);
1584 result.SetDv(Dv); result.SetZv(Zv);
1585 result.SetsigmaDv(sigmaDv); result.SetsigmaZv(sigmaZv);
1586 primarytrack.SetDv(Dv); primarytrack.SetZv(Zv);
1587 primarytrack.SetsigmaDv(sigmaDv); primarytrack.SetsigmaZv(sigmaZv);
1589 /////////////////////////////////////////////////////////////////////////////////////////////////
1591 primarytrack.PrimaryTrack(rl);
1592 TVector d2=primarytrack.Getd2();
1593 TVector tgl2=primarytrack.Gettgl2();
1594 TVector dtgl=primarytrack.Getdtgl();
1595 trackITS.Setd2(d2); trackITS.Settgl2(tgl2); trackITS.Setdtgl(dtgl);
1596 result.Setd2(d2); result.Settgl2(tgl2); result.Setdtgl(dtgl);
1598 trackITS.SetVertex(vertex); trackITS.SetErrorVertex(ervertex);
1599 result.SetVertex(vertex); result.SetErrorVertex(ervertex);
1601 Tracking(trackITS,&result,recPoints,vettid, flagvert,rl);
1603 // cout<<" progressive track number = "<<j<<"\r";
1605 // cout<<" progressive track number = "<<j<<"\n";
1606 Long_t labITS=result.GetLabel();
1607 // cout << " ITS track label = " << labITS << "\n";
1608 int lab=track->GetLabel();
1609 // cout << " TPC track label = " << lab <<"\n";
1610 //result.GetClusters(); getchar(); //to print the cluster matrix
1612 //propagation to vertex
1616 result.Propagation(rbeam);
1619 cov=&result.GetCMatrix();
1620 Double_t pt=TMath::Abs(result.GetPt());
1621 Double_t Dr=result.GetD();
1622 Double_t Z=result.GetZ();
1623 Double_t tgl=result.GetTgl();
1624 Double_t C=result.GetC();
1626 Double_t Dz=Z-(tgl/Cy)*TMath::ASin(result.arga(rbeam));
1629 // cout<<" Dr e dz alla fine = "<<Dr<<" "<<Dz<<"\n"; getchar();
1630 Double_t phi=result.Getphi();
1631 Double_t phivertex = phi - TMath::ASin(result.argA(rbeam));
1632 Double_t duepi=2.*TMath::Pi();
1633 if(phivertex>duepi) phivertex-=duepi;
1634 if(phivertex<0.) phivertex+=duepi;
1635 Double_t Dtot=TMath::Sqrt(Dr*Dr+Dz*Dz);
1637 //////////////////////////////////////////////////////////////////////////////////////////
1638 Int_t NumofCluster, idmodule,idpoint;
1639 NumofCluster=result.GetNumClust();
1640 if(NumofCluster >=5) {
1643 AliITSiotrack outtrack;
1647 iotrack->SetStatePhi(phi);
1648 iotrack->SetStateZ(Z);
1649 iotrack->SetStateD(Dr);
1650 iotrack->SetStateTgl(tgl);
1651 iotrack->SetStateC(C);
1652 Double_t radius=result.Getrtrack();
1653 iotrack->SetRadius(radius);
1655 if(C>0.) charge=-1; else charge=1;
1656 iotrack->SetCharge(charge);
1660 iotrack->SetCovMatrix(cov);
1662 Double_t px=pt*TMath::Cos(phi);
1663 Double_t py=pt*TMath::Sin(phi);
1666 Double_t xtrack=Dr*TMath::Sin(phi);
1667 Double_t ytrack=Dr*TMath::Cos(phi);
1668 Double_t ztrack=Dz+Vgeant(2);
1674 iotrack->SetX(xtrack);
1675 iotrack->SetY(ytrack);
1676 iotrack->SetZ(ztrack);
1677 iotrack->SetLabel(labITS);
1680 for(il=0;il<6; il++){
1681 iotrack->SetIdPoint(il,result.GetIdPoint(il));
1682 iotrack->SetIdModule(il,result.GetIdModule(il));
1686 //cout<<" labITS = "<<labITS<<"\n";
1687 //cout<<" phi z Dr tgl C = "<<phi<<" "<<Z<<" "<<Dr<<" "<<tgl<<" "<<C<<"\n"; getchar();
1689 DataOut(kkk) = ptg; kkk++; DataOut(kkk)=labITS; kkk++; DataOut(kkk)=lab; kkk++;
1691 for (il=0;il<6;il++) {
1692 idpoint=result.GetIdPoint(il);
1693 idmodule=result.GetIdModule(il);
1694 *(vettid[idmodule]+idpoint)=1;
1695 iotrack->SetIdPoint(il,idpoint);
1696 iotrack->SetIdModule(il,idmodule);
1699 // cout<<" +++++++++++++ pt e ptg = "<<pt<<" "<<ptg<<" ++++++++++\n";
1700 Double_t difpt= (pt-ptg)/ptg*100.;
1701 DataOut(kkk)=difpt; kkk++;
1702 Double_t lambdag=TMath::ATan(pzg/ptg);
1703 Double_t lam=TMath::ATan(tgl);
1704 Double_t diflam = (lam - lambdag)*1000.;
1705 DataOut(kkk) = diflam; kkk++;
1706 Double_t phig=TMath::ATan(pyg/pxg);
1707 Double_t phi=phivertex;
1708 Double_t difphi = (phi - phig)*1000.;
1709 DataOut(kkk)=difphi; kkk++;
1710 DataOut(kkk)=Dtot*1.e4; kkk++;
1711 DataOut(kkk)=Dr*1.e4; kkk++;
1712 DataOut(kkk)=Dz*1.e4; kkk++;
1714 for (r=0; r<9; r++) { out<<DataOut(r)<<" ";}
1719 } // end if on NumofCluster
1720 //gObjectTable->Print(); // stampa memoria
1721 } // end for (int j=min_t; j<=max_t; j++)
1726 static Bool_t first=kTRUE;
1727 static TFile *tfile;
1730 tfile=new TFile("itstracks.root","RECREATE");
1731 //cout<<"I have opened itstracks.root file "<<endl;
1738 sprintf(hname,"TreeT%d",evNumber);
1740 tracktree1.Write(hname);
1744 TTree *fAli=gAlice->TreeK();
1747 if (fAli) fileAli =fAli->GetCurrentFile();
1750 ////////////////////////////////////////////////////////////////////////////////////////////////
1752 printf("delete vectors\n");
1753 if(np) delete [] np;
1754 if(vettid) delete [] vettid;