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.40 2001/03/17 15:07:06 mariana
19 Update SDD response parameters
21 Revision 1.39 2001/03/12 17:45:32 hristov
22 Changes needed on Sun with CC 5.0
24 Revision 1.38 2001/03/07 14:04:51 barbera
25 Some vector dimensions increased to cope with full events
27 Revision 1.37 2001/03/07 12:36:35 barbera
28 A change added in the tracking part to manage delta rays
30 Revision 1.36 2001/03/02 19:44:11 barbera
31 modified to taking into account new version tracking v1
33 Revision 1.35 2001/02/28 18:16:46 mariana
34 Make the code compatible with the new AliRun
36 Revision 1.34 2001/02/11 15:51:39 mariana
37 Set protection in MakeBranch
39 Revision 1.33 2001/02/10 22:26:39 mariana
40 Move the initialization of the containers for raw clusters in MakeTreeC()
42 Revision 1.32 2001/02/08 23:55:31 nilsen
43 Removed fMajor/MinorVersion variables in favor of variables in derived classes.
44 Set arrays char *det[3] = {"SPD","SDD","SSD"} as const.
46 Revision 1.31 2001/02/02 23:57:28 nilsen
47 Added include file that are no londer included in AliITSgeom.h
49 Revision 1.30 2001/01/30 09:23:13 hristov
50 Streamers removed (R.Brun)
52 Revision 1.29 2001/01/26 20:01:09 hristov
53 Major upgrade of AliRoot code
55 Revision 1.28 2000/12/18 14:02:00 barbera
56 new version of the ITS tracking to take into account the new TPC track parametrization
58 Revision 1.27 2000/12/08 13:49:27 barbera
59 Hidden declaration in a for loop removed to be compliant with HP-UX compiler
61 Revision 1.26 2000/11/27 13:12:13 barbera
62 New version containing the files for tracking
64 Revision 1.25 2000/11/12 22:38:05 barbera
65 Added header file for the SPD Bari model
67 Revision 1.24 2000/10/09 22:18:12 barbera
68 Bug fixes from MAriana to le AliITStest.C run correctly
70 Revision 1.23 2000/10/05 20:47:42 nilsen
71 fixed dependencies of include files. Tryed but failed to get a root automaticly
72 generates streamer function to work. Modified SetDefaults.
74 Revision 1.9.2.15 2000/10/04 16:56:40 nilsen
75 Needed to include stdlib.h
78 Revision 1.22 2000/10/04 19:45:52 barbera
79 Corrected by F. Carminati for v3.04
81 Revision 1.21 2000/10/02 21:28:08 fca
82 Removal of useless dependecies via forward declarations
84 Revision 1.20 2000/10/02 16:31:39 barbera
87 Revision 1.9.2.14 2000/10/02 15:43:51 barbera
88 General code clean-up (e.g., printf -> cout)
90 Revision 1.19 2000/09/22 12:13:25 nilsen
91 Patches and updates for fixes to this and other routines.
93 Revision 1.18 2000/07/12 05:32:20 fca
94 Correcting several syntax problem with static members
96 Revision 1.17 2000/07/10 16:07:18 fca
97 Release version of ITS code
99 Revision 1.9.2.3 2000/02/02 13:42:09 barbera
100 fixed AliITS.cxx for new AliRun structure. Added ITS hits list to list of hits which will have their track numbers updated
102 Revision 1.9.2.2 2000/01/23 03:03:13 nilsen
103 //fixed FillModule. Removed fi(fabs(xl)<dx....
105 Revision 1.9.2.1 2000/01/12 19:03:32 nilsen
106 This is the version of the files after the merging done in December 1999.
107 See the ReadMe110100.txt file for details
109 Revision 1.9 1999/11/14 14:33:25 fca
110 Correct problems with distructors and pointers, thanks to I.Hrivnacova
112 Revision 1.8 1999/09/29 09:24:19 fca
113 Introduction of the Copyright and cvs Log
117 ///////////////////////////////////////////////////////////////////////////////
119 // An overview of the basic philosophy of the ITS code development
120 // and analysis is show in the figure below.
123 <img src="picts/ITS/ITS_Analysis_schema.gif">
126 <font size=+2 color=red>
127 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
128 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
134 // AliITS. Inner Traking System base class.
135 // This class contains the base procedures for the Inner Tracking System
139 <img src="picts/ITS/AliITS_Class_Diagram.gif">
142 <font size=+2 color=red>
143 <p>This show the class diagram of the different elements that are part of
151 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
154 // Modified and documented by Bjorn S. Nilsen
158 // Modified and documented by A. Bologna
161 // AliITS is the general base class for the ITS. Also see AliDetector for
162 // futher information.
164 ///////////////////////////////////////////////////////////////////////////////
165 #include <iostream.h>
173 #include <TClonesArray.h>
175 #include <TObjectTable.h>
179 #include <TParticle.h>
184 #include "AliITSMap.h"
185 #include "AliITSDetType.h"
186 #include "AliITSClusterFinder.h"
187 //#include "AliITSsimulation.h"
188 #include "AliITSsimulationSPD.h"
189 #include "AliITSsimulationSDD.h"
190 #include "AliITSsimulationSSD.h"
191 #include "AliITSresponse.h"
192 #include "AliITSsegmentationSPD.h"
193 #include "AliITSresponseSPD.h"
194 #include "AliITSresponseSPDbari.h"
195 #include "AliITSsegmentationSDD.h"
196 #include "AliITSresponseSDD.h"
197 #include "AliITSsegmentationSSD.h"
198 #include "AliITSresponseSSD.h"
199 #include "AliITShit.h"
200 #include "AliITSgeom.h"
201 #include "AliITSdigit.h"
202 #include "AliITSmodule.h"
203 #include "AliITSRecPoint.h"
204 #include "AliITSRawCluster.h"
208 #include "AliITStrack.h"
209 #include "AliITSiotrack.h"
210 #include "AliITStracking.h"
211 #include "AliITSRad.h"
212 #include "../TPC/AliTPC.h"
213 #include "../TPC/AliTPCParam.h"
218 //_____________________________________________________________________________
219 AliITS::AliITS() : AliDetector() {
221 // Default initialiser for ITS
222 // The default constructor of the AliITS class. In addition to
223 // creating the AliITS class it zeros the variables fIshunt (a member
224 // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
225 // fITSpoints, fIdSens, and fIdName. The AliDetector default constructor
233 fNDetTypes = kNTYPES;
252 //_____________________________________________________________________________
253 AliITS::AliITS(const char *name, const char *title):AliDetector(name,title){
255 // Default initialiser for ITS
256 // The constructor of the AliITS class. In addition to creating the
257 // AliITS class, it allocates memory for the TClonesArrays fHits and
258 // fDigits, and for the TObjArray fITSpoints. It also zeros the variables
259 // fIshunt (a member of AliDetector class), fEuclidOut, and fIdN, and zeros
260 // the pointers fIdSens and fIdName. To help in displaying hits via the ROOT
261 // macro display.C AliITS also sets the marker color to red. The variables
262 // passes with this constructor, const char *name and *title, are used by
263 // the constructor of AliDetector class. See AliDetector class for a
264 // description of these parameters and its constructor functions.
268 fHits = new TClonesArray("AliITShit", 1560);
269 gAlice->AddHitList(fHits);
271 fNDetTypes = kNTYPES;
273 fNdtype = new Int_t[kNTYPES];
274 fDtype = new TObjArray(kNTYPES);
276 fNctype = new Int_t[kNTYPES];
277 fCtype = new TObjArray(kNTYPES);
293 fDetTypes = new TObjArray(kNTYPES);
296 for(i=0;i<kNTYPES;i++) {
297 (*fDetTypes)[i]=new AliITSDetType();
303 SetMarkerColor(kRed);
307 //___________________________________________________________________________
308 AliITS::AliITS(AliITS &source){
310 if(this==&source) return;
311 Error("AliITS::Copy constructor",
312 "You are not allowed to make a copy of the AliITS");
315 //____________________________________________________________________________
316 AliITS& AliITS::operator=(AliITS &source){
317 // assignment operator
318 if(this==&source) return *this;
319 Error("AliITS::operator=",
320 "You are not allowed to make a copy of the AliITS");
322 return *this; //fake return
324 //____________________________________________________________________________
325 void AliITS::ClearModules(){
326 //clear the modules TObjArray
328 if(fITSmodules) fITSmodules->Delete();
331 //_____________________________________________________________________________
334 // Default distructor for ITS
335 // The default destructor of the AliITS class. In addition to deleting
336 // the AliITS class it deletes the memory pointed to by the fHits, fDigits,
337 // fIdSens, fIdName, and fITSpoints.
344 // delete fIdName; // TObjArray of TObjStrings
345 if(fIdName!=0) delete[] fIdName; // Array of TStrings
346 if(fIdSens!=0) delete[] fIdSens;
348 this->ClearModules();
350 }// end if fITSmodules!=0
370 if (fTreeC) delete fTreeC;
372 if (fITSgeom) delete fITSgeom;
376 //___________________________________________
377 AliITSDetType* AliITS::DetType(Int_t id)
379 //return pointer to id detector type
380 return ((AliITSDetType*) (*fDetTypes)[id]);
383 //___________________________________________
384 void AliITS::SetClasses(Int_t id, const char *digit, const char *cluster)
386 //set the digit and cluster classes to be used for the id detector type
387 ((AliITSDetType*) (*fDetTypes)[id])->ClassNames(digit,cluster);
390 //___________________________________________
391 void AliITS::SetResponseModel(Int_t id, AliITSresponse *response)
393 //set the response model for the id detector type
395 ((AliITSDetType*) (*fDetTypes)[id])->ResponseModel(response);
399 //___________________________________________
400 void AliITS::SetSegmentationModel(Int_t id, AliITSsegmentation *seg)
402 //set the segmentation model for the id detector type
404 ((AliITSDetType*) (*fDetTypes)[id])->SegmentationModel(seg);
408 //___________________________________________
409 void AliITS::SetSimulationModel(Int_t id, AliITSsimulation *sim)
411 //set the simulation model for the id detector type
413 ((AliITSDetType*) (*fDetTypes)[id])->SimulationModel(sim);
416 //___________________________________________
417 void AliITS::SetReconstructionModel(Int_t id, AliITSClusterFinder *reconst)
419 //set the cluster finder model for the id detector type
421 ((AliITSDetType*) (*fDetTypes)[id])->ReconstructionModel(reconst);
425 //_____________________________________________________________________________
426 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
429 // The function to add information to the AliITShit class. See the
430 // AliITShit class for a full description. This function allocates the
431 // necessary new space for the hit information and passes the variable
432 // track, and the pointers *vol and *hits to the AliITShit constructor
435 TClonesArray &lhits = *fHits;
436 new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
438 //_____________________________________________________________________________
439 void AliITS::AddRealDigit(Int_t id, Int_t *digits)
441 // add a real digit - as coming from data
443 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
444 new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
447 //_____________________________________________________________________________
448 void AliITS::AddSimDigit(Int_t id, AliITSdigit *d)
451 // add a simulated digit
453 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
458 new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
461 new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
464 new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
470 //_____________________________________________________________________________
471 void AliITS::AddSimDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,Int_t *hits,Float_t *charges){
473 // add a simulated digit to the list
475 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
479 new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks,hits);
482 new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,hits,charges);
485 new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks,hits);
491 //_____________________________________________________________________________
492 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c)
495 // add a cluster to the list
497 TClonesArray &lcl = *((TClonesArray*)(*fCtype)[id]);
502 new(lcl[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
505 new(lcl[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
508 new(lcl[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
515 //_____________________________________________________________________________
516 void AliITS::AddRecPoint(const AliITSRecPoint &r)
519 // Add a reconstructed space point to the list
521 TClonesArray &lrecp = *fRecPoints;
522 new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
526 //____________________________________________
527 void AliITS::ResetDigits()
530 // Reset number of digits and the digits array for the ITS detector
536 for (i=0;i<kNTYPES;i++ ) {
537 if ((*fDtype)[i]) ((TClonesArray*)(*fDtype)[i])->Clear();
538 if (fNdtype) fNdtype[i]=0;
542 //____________________________________________
543 void AliITS::ResetDigits(Int_t i)
546 // Reset number of digits and the digits array for this branch
548 if ((*fDtype)[i]) ((TClonesArray*)(*fDtype)[i])->Clear();
549 if (fNdtype) fNdtype[i]=0;
553 //____________________________________________
554 void AliITS::ResetClusters()
557 // Reset number of clusters and the clusters array for ITS
561 for (i=0;i<kNTYPES;i++ ) {
562 if ((*fCtype)[i]) ((TClonesArray*)(*fCtype)[i])->Clear();
563 if (fNctype) fNctype[i]=0;
568 //____________________________________________
569 void AliITS::ResetClusters(Int_t i)
572 // Reset number of clusters and the clusters array for this branch
574 if ((*fCtype)[i]) ((TClonesArray*)(*fCtype)[i])->Clear();
575 if (fNctype) fNctype[i]=0;
580 //____________________________________________
581 void AliITS::ResetRecPoints()
584 // Reset number of rec points and the rec points array
586 if (fRecPoints) fRecPoints->Clear();
591 //_____________________________________________________________________________
592 Int_t AliITS::DistancetoPrimitive(Int_t , Int_t ){
594 // Distance from mouse to ITS on the screen. Dummy routine
595 // A dummy routine used by the ROOT macro display.C to allow for the
596 // use of the mouse (pointing device) in the macro. In general this should
597 // never be called. If it is it returns the number 9999 for any value of
603 //_____________________________________________________________________________
606 // Initialise ITS after it has been built
607 // This routine initializes the AliITS class. It is intended to be called
608 // from the Init function in AliITSv?. Besides displaying a banner
609 // indicating that it has been called it initializes the array fIdSens
610 // and sets the default segmentation, response, digit and raw cluster classes
611 // Therefore it should be called after a call to CreateGeometry.
618 for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
622 //_____________________________________________________________________________
623 void AliITS::SetDefaults()
625 // sets the default segmentation, response, digit and raw cluster classes
627 printf("SetDefaults\n");
629 AliITSDetType *iDetType;
635 if (!iDetType->GetSegmentationModel()) {
636 AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(fITSgeom);
637 SetSegmentationModel(0,seg0);
639 if (!iDetType->GetResponseModel()) {
640 SetResponseModel(0,new AliITSresponseSPD());
642 // set digit and raw cluster classes to be used
644 const char *kData0=(iDetType->GetResponseModel())->DataType();
645 if (strstr(kData0,"real")) {
646 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSPD");
647 } else iDetType->ClassNames("AliITSdigitSPD","AliITSRawClusterSPD");
651 if (!iDetType->GetResponseModel()) {
652 SetResponseModel(1,new AliITSresponseSDD());
654 AliITSresponse *resp1=iDetType->GetResponseModel();
655 if (!iDetType->GetSegmentationModel()) {
656 AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(fITSgeom,resp1);
657 SetSegmentationModel(1,seg1);
659 const char *kData1=(iDetType->GetResponseModel())->DataType();
660 const char *kopt=iDetType->GetResponseModel()->ZeroSuppOption();
661 if ((!strstr(kopt,"2D")) && (!strstr(kopt,"1D")) || strstr(kData1,"real") ) {
662 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
663 } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
667 if (!iDetType->GetSegmentationModel()) {
668 AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(fITSgeom);
669 SetSegmentationModel(2,seg2);
671 if (!iDetType->GetResponseModel()) {
672 SetResponseModel(2,new AliITSresponseSSD());
674 const char *kData2=(iDetType->GetResponseModel())->DataType();
675 if (strstr(kData2,"real")) {
676 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSSD");
677 } else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
680 Warning("SetDefaults","Only the three basic detector types are initialised!");
687 //_____________________________________________________________________________
688 void AliITS::SetDefaultSimulation()
693 //_____________________________________________________________________________
694 void AliITS::SetDefaultClusterFinders()
699 //_____________________________________________________________________________
701 void AliITS::MakeTreeC(Option_t *option)
703 // create a separate tree to store the clusters
705 cout << "AliITS::MakeTreeC" << endl;
707 const char *optC = strstr(option,"C");
708 if (optC && !fTreeC) fTreeC = new TTree("TC","Clusters in ITS");
711 Int_t buffersize = 4000;
714 const char *det[3] = {"SPD","SDD","SSD"};
719 // one branch for Clusters per type of detector
721 for (i=0; i<kNTYPES ;i++) {
722 AliITSDetType *iDetType=DetType(i);
723 iDetType->GetClassNames(digclass,clclass);
725 (*fCtype)[i] = new TClonesArray(clclass,10000);
726 if (kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
727 else sprintf(branchname,"%sClusters%d",GetName(),i+1);
728 if (fCtype && fTreeC) {
729 TreeC()->Branch(branchname,&((*fCtype)[i]), buffersize);
730 cout << "Making Branch " << branchname;
731 cout << " for Clusters of detector type " << i+1 << endl;
737 //_____________________________________________________________________________
738 void AliITS::GetTreeC(Int_t event)
741 cout << "AliITS::GetTreeC" << endl;
743 // get the clusters tree for this event and set the branch address
747 const char *det[3] = {"SPD","SDD","SSD"};
754 sprintf(treeName,"TreeC%d",event);
755 fTreeC = (TTree*)gDirectory->Get(treeName);
760 for (i=0; i<kNTYPES; i++) {
761 if (kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
762 else sprintf(branchname,"%sClusters%d",GetName(),i+1);
764 branch = fTreeC->GetBranch(branchname);
765 if (branch) branch->SetAddress(&((*fCtype)[i]));
769 Error("AliITS::GetTreeC",
770 "cannot find Clusters Tree for event:%d\n",event);
774 //_____________________________________________________________________________
775 void AliITS::MakeBranch(Option_t* option, char *file)
778 // Creates Tree branches for the ITS.
781 Int_t buffersize = 4000;
783 sprintf(branchname,"%s",GetName());
785 AliDetector::MakeBranch(option,file);
787 const char *cD = strstr(option,"D");
788 const char *cR = strstr(option,"R");
792 // one branch for digits per type of detector
794 const char *det[3] = {"SPD","SDD","SSD"};
800 for (i=0; i<kNTYPES ;i++) {
801 AliITSDetType *iDetType=DetType(i);
802 iDetType->GetClassNames(digclass,clclass);
804 if(!((*fDtype)[i])) (*fDtype)[i] = new TClonesArray(digclass,10000);
808 for (i=0; i<kNTYPES ;i++) {
809 if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
810 else sprintf(branchname,"%sDigits%d",GetName(),i+1);
811 if (fDtype && gAlice->TreeD()) {
812 gAlice->MakeBranchInTree(gAlice->TreeD(),
813 branchname, &((*fDtype)[i]), buffersize, file);
814 cout << "Making Branch " << branchname;
815 cout << " for digits of type "<< i+1 << endl;
822 // only one branch for rec points for all detector types
824 sprintf(branchname,"%sRecPoints",GetName());
826 if(!fRecPoints) fRecPoints=new TClonesArray("AliITSRecPoint",10000);
828 if (fRecPoints && gAlice->TreeR()) {
829 gAlice->MakeBranchInTree(gAlice->TreeR(),
830 branchname, &fRecPoints, buffersize, file) ;
831 cout << "Making Branch " << branchname;
832 cout << " for reconstructed space points" << endl;
837 //___________________________________________
838 void AliITS::SetTreeAddress()
841 // Set branch address for the Trees.
844 AliDetector::SetTreeAddress();
846 const char *det[3] = {"SPD","SDD","SSD"};
849 TTree *treeD = gAlice->TreeD();
850 TTree *treeR = gAlice->TreeR();
854 for (i=0; i<kNTYPES; i++) {
855 if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
856 else sprintf(branchname,"%sDigits%d",GetName(),i+1);
858 branch = treeD->GetBranch(branchname);
859 if (branch) branch->SetAddress(&((*fDtype)[i]));
866 sprintf(branchname,"%sRecPoints",GetName());
867 branch = treeR->GetBranch(branchname);
868 if (branch) branch->SetAddress(&fRecPoints);
874 //____________________________________________________________________________
875 void AliITS::InitModules(Int_t size,Int_t &nmodules){
877 //initialize the modules array
880 fITSmodules->Delete();
884 Int_t nl,indexMAX,index;
886 if(size<=0){ // default to using data stored in AliITSgeom
888 Error("AliITS::InitModules",
889 "in AliITS::InitModule fITSgeom not defined\n");
891 } // end if fITSgeom==0
892 nl = fITSgeom->GetNlayers();
893 indexMAX = fITSgeom->GetModuleIndex(nl,fITSgeom->GetNladders(nl),
894 fITSgeom->GetNdetectors(nl))+1;
896 fITSmodules = new TObjArray(indexMAX);
897 for(index=0;index<indexMAX;index++){
898 fITSmodules->AddAt( new AliITSmodule(index),index);
901 fITSmodules = new TObjArray(size);
902 for(index=0;index<size;index++) {
903 fITSmodules->AddAt( new AliITSmodule(index),index);
910 //____________________________________________________________________________
911 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t nmodules,Option_t *option,Text_t *filename){
913 // fill the modules with the sorted by module hits; add hits from background
917 static TTree *trH1; //Tree with background hits
918 static TClonesArray *fHits2; //List of hits for one track only
920 static Bool_t first=kTRUE;
922 const char *addBgr = strstr(option,"Add");
927 cout<<"filename "<<filename<<endl;
928 file=new TFile(filename);
929 cout<<"I have opened "<<filename<<" file "<<endl;
930 fHits2 = new TClonesArray("AliITShit",1000 );
935 // Get Hits Tree header from file
936 if(fHits2) fHits2->Clear();
937 if(trH1) delete trH1;
941 sprintf(treeName,"TreeH%d",bgrev);
942 trH1 = (TTree*)gDirectory->Get(treeName);
943 //printf("TrH1 %p of treename %s for event %d \n",trH1,treeName,bgrev);
946 Error("AliITS::FillModules",
947 "cannot find Hits Tree for event:%d\n",bgrev);
949 // Set branch addresses
952 sprintf(branchname,"%s",GetName());
953 if (trH1 && fHits2) {
954 branch = trH1->GetBranch(branchname);
955 if (branch) branch->SetAddress(&fHits2);
959 //Int_t ntracks1 =(Int_t)TrH1->GetEntries();
960 //printf("background - ntracks1 - %d\n",ntracks1);
963 //Int_t npart = gAlice->GetEvent(evnt);
964 //if(npart<=0) return;
965 TClonesArray *itsHits = this->Hits();
966 Int_t lay,lad,det,index;
970 TTree *iTH = gAlice->TreeH();
971 Int_t ntracks =(Int_t) iTH->GetEntries();
974 for(t=0; t<ntracks; t++){
977 Int_t nhits = itsHits->GetEntriesFast();
978 //printf("nhits %d\n",nhits);
979 if (!nhits) continue;
980 for(h=0; h<nhits; h++){
981 itsHit = (AliITShit *)itsHits->UncheckedAt(h);
982 itsHit->GetDetectorID(lay,lad,det);
983 // temporarily index=det-1 !!!
984 if(fITSgeom) index = fITSgeom->GetModuleIndex(lay,lad,det);
987 mod = this->GetModule(index);
988 mod->AddHit(itsHit,t,h);
989 } // end loop over hits
990 } // end loop over tracks
992 // open the file with background
996 ntracks =(Int_t)trH1->GetEntries();
997 //printf("background - ntracks1 %d\n",ntracks);
998 //printf("background - Start loop over tracks \n");
1001 for (track=0; track<ntracks; track++) {
1003 if (fHits2) fHits2->Clear();
1004 trH1->GetEvent(track);
1006 for(i=0;i<fHits2->GetEntriesFast();++i) {
1008 itsHit=(AliITShit*) (*fHits2)[i];
1009 itsHit->GetDetectorID(lay,lad,det);
1010 // temporarily index=det-1 !!!
1011 if(fITSgeom) index = fITSgeom->GetModuleIndex(lay,lad,det);
1014 mod = this->GetModule(index);
1015 mod->AddHit(itsHit,track,i);
1016 } // end loop over hits
1017 } // end loop over tracks
1019 TTree *fAli=gAlice->TreeK();
1022 if (fAli) fileAli =fAli->GetCurrentFile();
1027 //gObjectTable->Print();
1031 //____________________________________________________________________________
1033 void AliITS::SDigits2Digits()
1036 AliITSgeom *geom = GetITSgeom();
1039 AliITSDetType *iDetType;
1040 iDetType=DetType(0);
1041 AliITSsegmentationSPD *seg0=(AliITSsegmentationSPD*)iDetType->GetSegmentationModel();
1042 AliITSresponseSPD *res0 = (AliITSresponseSPD*)iDetType->GetResponseModel();
1043 AliITSsimulationSPD *sim0=new AliITSsimulationSPD(seg0,res0);
1044 SetSimulationModel(0,sim0);
1046 // printf("SPD dimensions %f %f \n",seg0->Dx(),seg0->Dz());
1047 // printf("SPD npixels %d %d \n",seg0->Npz(),seg0->Npx());
1048 // printf("SPD pitches %d %d \n",seg0->Dpz(0),seg0->Dpx(0));
1052 //Set response functions
1053 // SDD compression param: 2 fDecrease, 2fTmin, 2fTmax or disable, 2 fTolerance
1055 iDetType=DetType(1);
1056 AliITSresponseSDD *res1 = (AliITSresponseSDD*)iDetType->GetResponseModel();
1058 res1=new AliITSresponseSDD();
1059 SetResponseModel(1,res1);
1061 Float_t noise, baseline;
1062 res1->GetNoiseParam(noise,baseline);
1063 Float_t noise_after_el = res1->GetNoiseAfterElectronics();
1064 Float_t fCutAmp = baseline + 2.*noise_after_el;
1065 Int_t cp[8]={0,0,(int)fCutAmp,(int)fCutAmp,0,0,0,0}; //1D
1066 res1->SetCompressParam(cp);
1067 AliITSsegmentationSDD *seg1=(AliITSsegmentationSDD*)iDetType->GetSegmentationModel();
1069 seg1 = new AliITSsegmentationSDD(geom,res1);
1070 SetSegmentationModel(1,seg1);
1072 AliITSsimulationSDD *sim1=new AliITSsimulationSDD(seg1,res1);
1073 SetSimulationModel(1,sim1);
1076 iDetType=DetType(2);
1077 AliITSsegmentationSSD *seg2=(AliITSsegmentationSSD*)iDetType->GetSegmentationModel();
1078 AliITSresponseSSD *res2 = (AliITSresponseSSD*)iDetType->GetResponseModel();
1079 res2->SetSigmaSpread(3.,2.);
1080 AliITSsimulationSSD *sim2=new AliITSsimulationSSD(seg2,res2);
1081 SetSimulationModel(2,sim2);
1083 cerr<<"Digitizing ITS...\n";
1087 HitsToDigits(0,0,-1," ","All"," ");
1088 timer.Stop(); timer.Print();
1096 //____________________________________________________________________________
1097 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size, Option_t *option, Option_t *opt,Text_t *filename)
1099 // keep galice.root for signal and name differently the file for
1100 // background when add! otherwise the track info for signal will be lost !
1102 // the condition below will disappear when the geom class will be
1103 // initialised for all versions - for the moment it is only for v5 !
1104 // 7 is the SDD beam test version
1105 Int_t ver = this->IsVersion();
1106 if(ver!=5 && ver!=7) return;
1108 const char *all = strstr(opt,"All");
1109 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1110 cout<<" 1 AliITS "<<endl;
1112 InitModules(size,nmodules);
1113 cout<<" 2 AliITS "<<endl;
1114 FillModules(evNumber,bgrev,nmodules,option,filename);
1115 cout<<" 3 AliITS "<<endl;
1118 AliITSsimulation* sim;
1119 //TObjArray *branches=gAlice->TreeD()->GetListOfBranches();
1120 AliITSgeom *geom = GetITSgeom();
1123 Int_t lay, lad, detect;
1125 for (id=0;id<kNTYPES;id++) {
1126 if (!all && !det[id]) continue;
1127 //branch = (TBranch*)branches->UncheckedAt(id);
1128 AliITSDetType *iDetType=DetType(id);
1129 sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1132 Error("HitsToDigits","The simulation class was not instantiated!");
1134 // or SetDefaultSimulation();
1138 first = geom->GetStartDet(id);
1139 last = geom->GetLastDet(id);
1140 } else first=last=0;
1141 cout << "det type " << id << " first, last "<< first << last << endl;
1142 for(module=first;module<=last;module++) {
1143 AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
1145 geom->GetModuleId(module,lay, lad, detect);
1147 ((AliITSsegmentationSSD*)(((AliITSsimulationSSD*)sim)->GetSegmentation()))->SetLayer(6);
1149 ((AliITSsegmentationSSD*)(((AliITSsimulationSSD*)sim)->GetSegmentation()))->SetLayer(5);
1151 sim->DigitiseModule(mod,module,evNumber);
1152 // fills all branches - wasted disk space
1153 gAlice->TreeD()->Fill();
1155 // try and fill only the branch
1158 } // loop over modules
1159 } // loop over detector types
1163 Int_t nentries=(Int_t)gAlice->TreeD()->GetEntries();
1164 cout << "nentries in TreeD" << nentries << endl;
1167 sprintf(hname,"TreeD%d",evNumber);
1168 gAlice->TreeD()->Write(hname,TObject::kOverwrite);
1170 gAlice->TreeD()->Reset();
1175 //____________________________________________________________________________
1176 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt)
1178 // cluster finding and reconstruction of space points
1180 // the condition below will disappear when the geom class will be
1181 // initialised for all versions - for the moment it is only for v5 !
1182 // 7 is the SDD beam test version
1183 Int_t ver = this->IsVersion();
1186 const char *all = strstr(opt,"All");
1187 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1189 static Bool_t first=kTRUE;
1190 if (!TreeC() && first) {
1195 TTree *treeC=TreeC();
1199 AliITSClusterFinder* rec;
1201 //TObjArray *branches=gAlice->TreeR()->GetListOfBranches();
1202 AliITSgeom *geom = GetITSgeom();
1205 for (id=0;id<kNTYPES;id++) {
1206 if (!all && !det[id]) continue;
1207 //branch = (TBranch*)branches->UncheckedAt(id);
1208 AliITSDetType *iDetType=DetType(id);
1209 rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1211 Error("DigitsToRecPoints","The cluster finder class was not instantiated!");
1213 // or SetDefaultClusterFinders();
1215 TClonesArray *itsDigits = this->DigitsAddress(id);
1219 first = geom->GetStartDet(id);
1220 last = geom->GetLastDet(id);
1221 } else first=last=0;
1222 //printf("first last %d %d\n",first,last);
1223 for(module=first;module<=last;module++) {
1224 this->ResetDigits();
1225 if (all) gAlice->TreeD()->GetEvent(lastentry+module);
1226 else gAlice->TreeD()->GetEvent(lastentry+(module-first));
1227 Int_t ndigits = itsDigits->GetEntriesFast();
1228 if (ndigits) rec->FindRawClusters();
1229 gAlice->TreeR()->Fill();
1233 // try and fill only the branch
1235 //ResetRecPoints(id);
1236 } // loop over modules
1237 } // loop over detector types
1240 Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1241 Int_t ncentries=(Int_t)treeC->GetEntries();
1242 cout << " nentries ncentries " << nentries << ncentries << endl;
1245 sprintf(hname,"TreeR%d",evNumber);
1246 gAlice->TreeR()->Write(hname,TObject::kOverwrite);
1248 gAlice->TreeR()->Reset();
1250 sprintf(hname,"TreeC%d",evNumber);
1251 treeC->Write(hname,TObject::kOverwrite);
1256 //____________________________________________________________________________
1257 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
1258 Option_t *option,Option_t *opt,Text_t *filename)
1260 // keep galice.root for signal and name differently the file for
1261 // background when add! otherwise the track info for signal will be lost !
1264 // the condition below will disappear when the geom class will be
1265 // initialised for all versions - for the moment it is only for v5 !
1266 Int_t ver = this->IsVersion();
1269 const char *all = strstr(opt,"All");
1270 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1273 InitModules(size,nmodules);
1274 FillModules(evNumber,bgrev,nmodules,option,filename);
1277 AliITSsimulation* sim;
1278 AliITSgeom *geom = GetITSgeom();
1280 TRandom *random=new TRandom[9];
1281 random[0].SetSeed(111);
1282 random[1].SetSeed(222);
1283 random[2].SetSeed(333);
1284 random[3].SetSeed(444);
1285 random[4].SetSeed(555);
1286 random[5].SetSeed(666);
1287 random[6].SetSeed(777);
1288 random[7].SetSeed(888);
1289 random[8].SetSeed(999);
1293 for (id=0;id<kNTYPES;id++) {
1294 if (!all && !det[id]) continue;
1295 AliITSDetType *iDetType=DetType(id);
1296 sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1298 Error("HitsToFastPoints","The simulation class was not instantiated!");
1300 // or SetDefaultSimulation();
1302 Int_t first = geom->GetStartDet(id);
1303 Int_t last = geom->GetLastDet(id);
1304 for(module=first;module<=last;module++) {
1305 AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
1306 sim->CreateFastRecPoints(mod,module,random);
1307 gAlice->TreeR()->Fill();
1309 } // loop over modules
1310 } // loop over detector types
1315 //Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1318 sprintf(hname,"TreeR%d",evNumber);
1319 gAlice->TreeR()->Write(hname,TObject::kOverwrite);
1321 gAlice->TreeR()->Reset();
1327 //________________________________________________________________
1329 AliITStrack AliITS::Tracking(AliITStrack &track, AliITStrack *reference,TObjArray *fastpoints, Int_t
1330 **vettid, Bool_t flagvert, AliITSRad *rl ) {
1332 //Origin A. Badala' and G.S. Pappalardo: e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1335 TList *list= new TList();
1337 AliITStrack tr(track);
1341 Double_t Pt=(tr).GetPt();
1342 //cout << "\n Pt = " << Pt <<"\n";
1344 AliITStracking obj(list, reference, this, fastpoints,TMath::Abs(Pt),vettid, flagvert, rl);
1349 TVector VecTotLabref(18);
1351 for(lay=5; lay>=0; lay--) {
1352 TVector VecLabref(3);
1353 VecLabref=(*reference).GetLabTrack(lay);
1354 Float_t ClustZ=(*reference).GetZclusterTrack( lay); //modified il 5-3-2001
1355 for(k=0; k<3; k++) { // {itot++; VecTotLabref(itot)=VecLabref(k);} // modified 5-3-2002
1357 Int_t lpp=(Int_t)VecLabref(k);
1359 TParticle *p=(TParticle*) gAlice->Particle(lpp);
1360 Int_t pcode=p->GetPdgCode();
1361 if(pcode==11) VecLabref(k)=p->GetFirstMother();
1364 itot++; VecTotLabref(itot)=VecLabref(k);
1365 if(VecLabref(k)==0. && ClustZ == 0.) VecTotLabref(itot) =-3.;
1370 (*reference).Search(VecTotLabref, labref, freq);
1372 if(freq < 5) labref=-labref;
1373 (*reference).SetLabel(labref);
1381 //________________________________________________________________
1385 void AliITS::DoTracking(Int_t evNumber, Int_t min_t, Int_t max_t, TFile *file, Bool_t flagvert) {
1387 // ex macro for tracking ITS
1389 printf("begin DoTracking - file %p\n",file);
1391 const char *pname="75x40_100x60";
1393 Int_t imax=200,jmax=450;
1394 AliITSRad *rl = new AliITSRad(imax,jmax);
1395 //cout<<" dopo costruttore AliITSRad\n"; getchar();
1399 Float_t px,py,pz,x,y,z,pxg,pyg,pzg,ptg;
1404 gAlice->GetEvent(0);
1406 AliTPC *TPC=(AliTPC*)gAlice->GetDetector("TPC");
1407 AliTPCParam *digp = (AliTPCParam*)file->Get(pname);
1408 if (digp!=0) TPC->SetParam(digp);
1410 GoodTrack gt[15000];
1412 ifstream in("itsgood_tracks");
1414 cerr<<"Reading itsgood tracks...\n";
1415 while (in>>gt[ngood].lab>>gt[ngood].code
1416 >>gt[ngood].px >>gt[ngood].py>>gt[ngood].pz
1417 >>gt[ngood].x >>gt[ngood].y >>gt[ngood].z
1418 >>gt[ngood].pxg >>gt[ngood].pyg >>gt[ngood].pzg
1419 >>gt[ngood].ptg >>gt[ngood].flag) {
1423 cerr<<"Too many good tracks !\n";
1427 if (!in.eof()) cerr<<"Read error (itsgood_tracks) !\n";
1431 TFile *tf=TFile::Open("tpctracks.root");
1432 if (!tf->IsOpen()) {cerr<<"Can't open tpctracks.root !\n"; return ;}
1433 TObjArray tracks(200000);
1434 TTree *tracktree=(TTree*)tf->Get("TreeT");
1435 TBranch *tbranch=tracktree->GetBranch("tracks");
1436 Int_t nentr=(Int_t)tracktree->GetEntries();
1438 for (kk=0; kk<nentr; kk++) {
1439 AliTPCtrack *iotrack=new AliTPCtrack;
1440 tbranch->SetAddress(&iotrack);
1441 tracktree->GetEvent(kk);
1442 tracks.AddLast(iotrack);
1447 Int_t nt = tracks.GetEntriesFast();
1448 cerr<<"Number of found tracks "<<nt<<endl;
1453 Double_t ptg=0.,pxg=0.,pyg=0.,pzg=0.;
1455 ////////////////////////////// good tracks definition in TPC ////////////////////////////////
1457 ofstream out1 ("AliITSTrag.out");
1459 for (i=0; i<ngood; i++) out1 << gt[i].ptg << "\n";
1464 TTree *TR=gAlice->TreeR();
1465 Int_t nent=(Int_t)TR->GetEntries();
1466 TClonesArray *recPoints = RecPoints();
1468 Int_t totalpoints=0;
1469 Int_t *np = new Int_t[nent];
1470 Int_t **vettid = new Int_t* [nent];
1472 for (mod=0; mod<nent; mod++) {
1474 this->ResetRecPoints();
1475 //gAlice->TreeR()->GetEvent(mod+1); //first entry in TreeR is empty
1476 gAlice->TreeR()->GetEvent(mod); //first entry in TreeR is empty
1477 numbpoints = recPoints->GetEntries();
1478 totalpoints+=numbpoints;
1479 np[mod] = numbpoints;
1480 //cout<<" mod = "<<mod<<" numbpoints = "<<numbpoints<<"\n"; getchar();
1481 vettid[mod] = new Int_t[numbpoints];
1483 for (ii=0;ii<numbpoints; ii++) *(vettid[mod]+ii)=0;
1489 if(min_t < 0) {min_t = 0; max_t = nt-1;}
1492 ///////////////////////////////// Definition of vertex end its error ////////////////////////////
1493 ////////////////////////// In the future it will be given by a method ///////////////////////////
1498 Float_t sigmavx=0.0050; // 50 microns
1499 Float_t sigmavy=0.0050; // 50 microns
1500 Float_t sigmavz=0.010; // 100 microns
1502 //Vx+=gRandom->Gaus(0,sigmavx); Vy+=gRandom->Gaus(0,sigmavy); Vz+=gRandom->Gaus(0,sigmavz);
1503 TVector vertex(3), ervertex(3)
1504 vertex(0)=Vx; vertex(1)=Vy; vertex(2)=Vz;
1505 ervertex(0)=sigmavx; ervertex(1)=sigmavy; ervertex(2)=sigmavz;
1506 /////////////////////////////////////////////////////////////////////////////////////////////////
1509 //TDirectory *savedir=gDirectory;
1511 TTree tracktree1("TreeT","Tree with ITS tracks");
1512 AliITSiotrack *iotrack=0;
1513 tracktree1.Branch("ITStracks","AliITSiotrack",&iotrack,32000,0);
1515 ofstream out ("AliITSTra.out");
1518 for (j=min_t; j<=max_t; j++) {
1519 track=(AliTPCtrack*)tracks.UncheckedAt(j);
1521 if (!track) continue;
1522 ////// elimination of not good tracks ////////////
1523 Int_t ilab=TMath::Abs(track->GetLabel());
1525 for (iii=0;iii<ngood;iii++) {
1526 //cout<<" ilab, gt[iii].lab = "<<ilab<<" "<<gt[iii].lab<<"\n"; getchar();
1527 if (ilab==gt[iii].lab) {
1536 //cout<<" j flaglab = " <<j<<" "<<flaglab<<"\n"; getchar();
1537 if (!flaglab) continue;
1538 //cout<<" j = " <<j<<"\n"; getchar();
1540 ////// old propagation to the end of TPC //////////////
1542 track->PropagateTo(xk);
1544 track->PropagateTo(xk,42.7,2.27); //C
1546 track->PropagateTo(xk,36.2,1.98e-3); //C02
1548 track->PropagateTo(xk,42.7,2.27); //C
1549 ///////////////////////////////////////////////////
1552 ////// new propagation to the end of TPC //////////////
1554 track->PropagateTo(xk, 28.94, 1.204e-3); //Ne
1556 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
1558 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
1560 track->PropagateTo(xk, 41.28, 0.029); //Nomex
1562 track->PropagateTo(xk,36.2,1.98e-3); //C02
1564 track->PropagateTo(xk, 24.01, 2.7); //Al
1566 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
1568 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
1570 track->PropagateTo(xk, 41.28, 0.029); //Nomex
1572 ///////////////////////////////////////////////////////////////
1574 ///////////////////////////////////////////////////////////////
1575 AliITStrack trackITS(*track);
1576 AliITStrack result(*track);
1577 AliITStrack primarytrack(*track);
1579 ///////////////////////////////////////////////////////////////////////////////////////////////
1581 Vgeant=result.GetVertex();
1583 // Definition of Dv and Zv for vertex constraint
1584 Double_t sigmaDv=0.0050; Double_t sigmaZv=0.010;
1585 //Double_t sigmaDv=0.0015; Double_t sigmaZv=0.0015;
1586 Double_t uniform= gRandom->Uniform();
1588 if(uniform<=0.5) signdv=-1.;
1592 Double_t Vr=TMath::Sqrt(Vgeant(0)*Vgeant(0)+ Vgeant(1)*Vgeant(1));
1593 Double_t Dv=gRandom->Gaus(signdv*Vr,(Float_t)sigmaDv);
1594 Double_t Zv=gRandom->Gaus(Vgeant(2),(Float_t)sigmaZv);
1596 //cout<<" Dv e Zv = "<<Dv<<" "<<Zv<<"\n";
1597 trackITS.SetDv(Dv); trackITS.SetZv(Zv);
1598 trackITS.SetsigmaDv(sigmaDv); trackITS.SetsigmaZv(sigmaZv);
1599 result.SetDv(Dv); result.SetZv(Zv);
1600 result.SetsigmaDv(sigmaDv); result.SetsigmaZv(sigmaZv);
1601 primarytrack.SetDv(Dv); primarytrack.SetZv(Zv);
1602 primarytrack.SetsigmaDv(sigmaDv); primarytrack.SetsigmaZv(sigmaZv);
1604 /////////////////////////////////////////////////////////////////////////////////////////////////
1606 primarytrack.PrimaryTrack(rl);
1607 TVector d2=primarytrack.Getd2();
1608 TVector tgl2=primarytrack.Gettgl2();
1609 TVector dtgl=primarytrack.Getdtgl();
1610 trackITS.Setd2(d2); trackITS.Settgl2(tgl2); trackITS.Setdtgl(dtgl);
1611 result.Setd2(d2); result.Settgl2(tgl2); result.Setdtgl(dtgl);
1613 trackITS.SetVertex(vertex); trackITS.SetErrorVertex(ervertex);
1614 result.SetVertex(vertex); result.SetErrorVertex(ervertex);
1616 Tracking(trackITS,&result,recPoints,vettid, flagvert,rl);
1618 // cout<<" progressive track number = "<<j<<"\r";
1620 // cout<<" progressive track number = "<<j<<"\n";
1621 Long_t labITS=result.GetLabel();
1622 // cout << " ITS track label = " << labITS << "\n";
1623 int lab=track->GetLabel();
1624 // cout << " TPC track label = " << lab <<"\n";
1625 //result.GetClusters(); getchar(); //to print the cluster matrix
1627 //propagation to vertex
1631 result.Propagation(rbeam);
1634 cov=&result.GetCMatrix();
1635 Double_t pt=TMath::Abs(result.GetPt());
1636 Double_t Dr=result.GetD();
1637 Double_t Z=result.GetZ();
1638 Double_t tgl=result.GetTgl();
1639 Double_t C=result.GetC();
1641 Double_t Dz=Z-(tgl/Cy)*TMath::ASin(result.arga(rbeam));
1644 // cout<<" Dr e dz alla fine = "<<Dr<<" "<<Dz<<"\n"; getchar();
1645 Double_t phi=result.Getphi();
1646 Double_t phivertex = phi - TMath::ASin(result.argA(rbeam));
1647 Double_t duepi=2.*TMath::Pi();
1648 if(phivertex>duepi) phivertex-=duepi;
1649 if(phivertex<0.) phivertex+=duepi;
1650 Double_t Dtot=TMath::Sqrt(Dr*Dr+Dz*Dz);
1652 //////////////////////////////////////////////////////////////////////////////////////////
1653 Int_t NumofCluster, idmodule,idpoint;
1654 NumofCluster=result.GetNumClust();
1655 if(NumofCluster >=5) {
1658 AliITSiotrack outtrack;
1662 iotrack->SetStatePhi(phi);
1663 iotrack->SetStateZ(Z);
1664 iotrack->SetStateD(Dr);
1665 iotrack->SetStateTgl(tgl);
1666 iotrack->SetStateC(C);
1667 Double_t radius=result.Getrtrack();
1668 iotrack->SetRadius(radius);
1670 if(C>0.) charge=-1; else charge=1;
1671 iotrack->SetCharge(charge);
1675 iotrack->SetCovMatrix(cov);
1677 Double_t px=pt*TMath::Cos(phi);
1678 Double_t py=pt*TMath::Sin(phi);
1681 Double_t xtrack=Dr*TMath::Sin(phi);
1682 Double_t ytrack=Dr*TMath::Cos(phi);
1683 Double_t ztrack=Dz+Vgeant(2);
1689 iotrack->SetX(xtrack);
1690 iotrack->SetY(ytrack);
1691 iotrack->SetZ(ztrack);
1692 iotrack->SetLabel(labITS);
1695 for(il=0;il<6; il++){
1696 iotrack->SetIdPoint(il,result.GetIdPoint(il));
1697 iotrack->SetIdModule(il,result.GetIdModule(il));
1701 //cout<<" labITS = "<<labITS<<"\n";
1702 //cout<<" phi z Dr tgl C = "<<phi<<" "<<Z<<" "<<Dr<<" "<<tgl<<" "<<C<<"\n"; getchar();
1704 DataOut(kkk) = ptg; kkk++; DataOut(kkk)=labITS; kkk++; DataOut(kkk)=lab; kkk++;
1706 for (il=0;il<6;il++) {
1707 idpoint=result.GetIdPoint(il);
1708 idmodule=result.GetIdModule(il);
1709 *(vettid[idmodule]+idpoint)=1;
1710 iotrack->SetIdPoint(il,idpoint);
1711 iotrack->SetIdModule(il,idmodule);
1714 // cout<<" +++++++++++++ pt e ptg = "<<pt<<" "<<ptg<<" ++++++++++\n";
1715 Double_t difpt= (pt-ptg)/ptg*100.;
1716 DataOut(kkk)=difpt; kkk++;
1717 Double_t lambdag=TMath::ATan(pzg/ptg);
1718 Double_t lam=TMath::ATan(tgl);
1719 Double_t diflam = (lam - lambdag)*1000.;
1720 DataOut(kkk) = diflam; kkk++;
1721 Double_t phig=TMath::ATan(pyg/pxg);
1722 Double_t phi=phivertex;
1723 Double_t difphi = (phi - phig)*1000.;
1724 DataOut(kkk)=difphi; kkk++;
1725 DataOut(kkk)=Dtot*1.e4; kkk++;
1726 DataOut(kkk)=Dr*1.e4; kkk++;
1727 DataOut(kkk)=Dz*1.e4; kkk++;
1729 for (r=0; r<9; r++) { out<<DataOut(r)<<" ";}
1734 } // end if on NumofCluster
1735 //gObjectTable->Print(); // stampa memoria
1736 } // end for (int j=min_t; j<=max_t; j++)
1741 static Bool_t first=kTRUE;
1742 static TFile *tfile;
1745 tfile=new TFile("itstracks.root","RECREATE");
1746 //cout<<"I have opened itstracks.root file "<<endl;
1753 sprintf(hname,"TreeT%d",evNumber);
1755 tracktree1.Write(hname);
1759 TTree *fAli=gAlice->TreeK();
1762 if (fAli) fileAli =fAli->GetCurrentFile();
1765 ////////////////////////////////////////////////////////////////////////////////////////////////
1767 printf("delete vectors\n");
1768 if(np) delete [] np;
1769 if(vettid) delete [] vettid;