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.13 2000/06/13 15:32:44 nilsen
19 fix compilation error on HP and DEC unix.
21 Revision 1.12 2000/06/12 23:43:16 nilsen
22 New ITS code replacing the old structure and simulations code.
24 Revision 1.9.2.8 2000/06/12 18:05:59 barbera
25 fixed posible compilation errors on HP unix
27 Revision 1.9.2.7 2000/06/11 20:20:18 barbera
28 New AliITS base clase for the new structure.
30 Revision 1.9.2.3 2000/02/02 13:42:09 barbera
31 fixed AliITS.cxx for new AliRun structure. Added ITS hits list to list of hits which will have their track numbers updated
33 Revision 1.9.2.2 2000/01/23 03:03:13 nilsen
34 //fixed FillModule. Removed fi(fabs(xl)<dx....
36 Revision 1.9.2.1 2000/01/12 19:03:32 nilsen
37 This is the version of the files after the merging done in December 1999.
38 See the ReadMe110100.txt file for details
40 Revision 1.9 1999/11/14 14:33:25 fca
41 Correct problems with distructors and pointers, thanks to I.Hrivnacova
43 Revision 1.8 1999/09/29 09:24:19 fca
44 Introduction of the Copyright and cvs Log
48 ///////////////////////////////////////////////////////////////////////////////
50 // An overview of the basic philosophy of the ITS code development
51 // and analysis is show in the figure below.
54 <img src="picts/ITS/ITS_Analysis_schema.gif">
57 <font size=+2 color=red>
58 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
59 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
65 // AliITS. Inner Traking System base class.
66 // This class contains the base procedures for the Inner Tracking System
70 <img src="picts/ITS/AliITS_Class_Diagram.gif">
73 <font size=+2 color=red>
74 <p>This show the class diagram of the different elements that are part of
82 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
85 // Modified and documented by Bjorn S. Nilsen
89 // Modified and documented by A. Bologna
92 // AliITS is the general base class for the ITS. Also see AliDetector for
93 // futher information.
95 ///////////////////////////////////////////////////////////////////////////////
100 #include <TObjArray.h>
101 #include <TClonesArray.h>
103 #include <TObjectTable.h>
109 #include "AliITSMap.h"
110 #include "AliITSClusterFinder.h"
111 #include "AliITSsimulation.h"
112 #include "AliITSsimulationFastPoints.h"
113 #include "AliITSsegmentationSPD.h"
114 #include "AliITSresponseSPD.h"
115 #include "AliITSsegmentationSDD.h"
116 #include "AliITSresponseSDD.h"
117 #include "AliITSsegmentationSSD.h"
118 #include "AliITSresponseSSD.h"
119 //#include "AliITStrack.h"
124 //_____________________________________________________________________________
125 AliITS::AliITS() : AliDetector() {
127 // Default initialiser for ITS
128 // The default constructor of the AliITS class. In addition to
129 // creating the AliITS class it zeros the variables fIshunt (a member
130 // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
131 // fITSpoints, fIdSens, and fIdName. The AliDetector default constructor
155 //_____________________________________________________________________________
156 AliITS::AliITS(const char *name, const char *title):AliDetector(name,title){
158 // Default initialiser for ITS
159 // The constructor of the AliITS class. In addition to creating the
160 // AliITS class, it allocates memory for the TClonesArrays fHits and
161 // fDigits, and for the TObjArray fITSpoints. It also zeros the variables
162 // fIshunt (a member of AliDetector class), fEuclidOut, and fIdN, and zeros
163 // the pointers fIdSens and fIdName. To help in displaying hits via the ROOT
164 // macro display.C AliITS also sets the marker color to red. The variables
165 // passes with this constructor, const char *name and *title, are used by
166 // the constructor of AliDetector class. See AliDetector class for a
167 // description of these parameters and its constructor functions.
170 fHits = new TClonesArray("AliITShit", 1560);
171 gAlice->AddHitList(fHits);
192 fDetTypes = new TObjArray(fNDetTypes);
195 for(i=0;i<fNDetTypes;i++) {
196 (*fDetTypes)[i]=new AliITSDetType();
200 SetMarkerColor(kRed);
204 //___________________________________________________________________________
205 AliITS::AliITS(AliITS &source){
206 if(this==&source) return;
207 printf("Error: You are not allowed to make a copy of the AliITS\n");
210 //____________________________________________________________________________
211 AliITS& AliITS::operator=(AliITS &source){
212 if(this==&source) return *this;
213 printf("Error: You are not allowed to make a copy of the AliITS\n");
216 //____________________________________________________________________________
217 void AliITS::ClearModules(){
218 //clear the modules TObjArray
222 Int_t indSPD = fITSgeom->GetModuleIndex(2,fITSgeom->GetNladders(2),
223 fITSgeom->GetNdetectors(2));
224 Int_t indSDD = fITSgeom->GetModuleIndex(4,fITSgeom->GetNladders(4),
225 fITSgeom->GetNdetectors(4));
226 for(i=0;i<fITSmodules->GetEntriesFast();i++){
228 delete (AliITSmodule *) fITSmodules->At(i);
230 delete (AliITSmodule *) fITSmodules->At(i);
232 delete (AliITSmodule *) fITSmodules->At(i);
234 }// end if fITSmodules!=0
237 //_____________________________________________________________________________
240 // Default distructor for ITS
241 // The default destructor of the AliITS class. In addition to deleting
242 // the AliITS class it deletes the memory pointed to by the fHits, fDigits,
243 // fIdSens, fIdName, and fITSpoints.
250 if(fIdName!=0) delete[] fIdName;
251 if(fIdSens!=0) delete[] fIdSens;
253 this->ClearModules();
255 }// end if fITSmodules!=0
260 for(i=0;i<fNDetTypes;i++) {
266 for(i=0;i<fNDetTypes;i++) {
278 if (fTreeC) delete fTreeC;
282 //___________________________________________
283 AliITSDetType* AliITS::DetType(Int_t id)
285 //return pointer to id detector type
286 return ((AliITSDetType*) (*fDetTypes)[id]);
289 //___________________________________________
290 void AliITS::SetClasses(Int_t id, char* digit, char* cluster)
292 //set the digit and cluster classes to be used for the id detector type
293 ((AliITSDetType*) (*fDetTypes)[id])->ClassNames(digit,cluster);
296 //___________________________________________
297 void AliITS::SetResponseModel(Int_t id, AliITSresponse *response)
299 //set the response model for the id detector type
301 ((AliITSDetType*) (*fDetTypes)[id])->ResponseModel(response);
305 //___________________________________________
306 void AliITS::SetSegmentationModel(Int_t id, AliITSsegmentation *seg)
308 //set the segmentation model for the id detector type
310 ((AliITSDetType*) (*fDetTypes)[id])->SegmentationModel(seg);
314 //___________________________________________
315 void AliITS::SetSimulationModel(Int_t id, AliITSsimulation *sim)
317 //set the simulation model for the id detector type
319 ((AliITSDetType*) (*fDetTypes)[id])->SimulationModel(sim);
322 //___________________________________________
323 void AliITS::SetReconstructionModel(Int_t id, AliITSClusterFinder *reconst)
325 //set the cluster finder model for the id detector type
327 ((AliITSDetType*) (*fDetTypes)[id])->ReconstructionModel(reconst);
331 //_____________________________________________________________________________
332 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
335 // The function to add information to the AliITShit class. See the
336 // AliITShit class for a full description. This function allocates the
337 // necessary new space for the hit information and passes the variable
338 // track, and the pointers *vol and *hits to the AliITShit constructor
341 TClonesArray &lhits = *fHits;
342 new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
344 //_____________________________________________________________________________
345 void AliITS::AddRealDigit(Int_t id, Int_t *digits)
347 // add a real digit - as coming from data
349 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
350 new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
353 //_____________________________________________________________________________
354 void AliITS::AddDigit(Int_t id, AliITSdigit *d)
357 // add a simulated digit
359 // should have ctors of type AliITSdigitSDD(const AliITSdigitSDD &)
361 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
366 new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
369 new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
372 new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
378 //_____________________________________________________________________________
379 void AliITS::AddDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,Float_t *charges){
381 // add a simulated digit to the list
383 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
387 new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks);
390 new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,charges);
393 new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks);
399 //_____________________________________________________________________________
400 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c)
403 // add a cluster to the list
405 // should have ctors of type AliITSRawClusterSDD(const AliITSRawClusterSDD &)
407 TClonesArray &lcl = *((TClonesArray*)(*fCtype)[id]);
412 new(lcl[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
415 new(lcl[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
418 new(lcl[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
425 //_____________________________________________________________________________
426 void AliITS::AddRecPoint(const AliITSRecPoint &r)
429 // Add a reconstructed space point to the list
431 TClonesArray &lrecp = *fRecPoints;
432 new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
436 //____________________________________________
437 void AliITS::ResetDigits()
440 // Reset number of digits and the digits array for thE ITS detector
446 for(i=0;i<fNDetTypes;i++ ) {
447 if ((*fDtype)[i]) ((TClonesArray*)(*fDtype)[i])->Clear();
448 if (fNdtype) fNdtype[i]=0;
452 //____________________________________________
453 void AliITS::ResetDigits(Int_t i)
456 // Reset number of digits and the digits array for this branch
458 if ((*fDtype)[i]) ((TClonesArray*)(*fDtype)[i])->Clear();
459 if (fNdtype) fNdtype[i]=0;
462 //____________________________________________
463 void AliITS::ResetClusters()
466 // Reset number of clusters and the clusters array for ITS
470 for(i=0;i<fNDetTypes;i++ ) {
471 if ((*fCtype)[i]) ((TClonesArray*)(*fCtype)[i])->Clear();
472 if (fNctype) fNctype[i]=0;
476 //____________________________________________
477 void AliITS::ResetClusters(Int_t i)
480 // Reset number of clusters and the clusters array for this branch
482 if ((*fCtype)[i]) ((TClonesArray*)(*fCtype)[i])->Clear();
483 if (fNctype) fNctype[i]=0;
488 //____________________________________________
489 void AliITS::ResetRecPoints()
492 // Reset number of rec points and the rec points array
495 if (fRecPoints) fRecPoints->Clear();
499 //_____________________________________________________________________________
500 Int_t AliITS::DistancetoPrimitive(Int_t , Int_t ){
502 // Distance from mouse to ITS on the screen. Dummy routine
503 // A dummy routine used by the ROOT macro display.C to allow for the
504 // use of the mouse (pointing device) in the macro. In general this should
505 // never be called. If it is it returns the number 9999 for any value of
511 //_____________________________________________________________________________
514 // Initialise ITS after it has been built
515 // This routine initializes the AliITS class. It is intended to be called
516 // from the Init function in AliITSv?. Besides displaying a banner
517 // indicating that it has been called it initializes the array fIdSens
518 // and sets the default segmentation, response, digit and raw cluster classes
519 // Therefore it should be called after a call to CreateGeometry.
528 for(i=0;i<35;i++) printf("*");
529 printf(" ITS_INIT ");
530 for(i=0;i<35;i++) printf("*");
534 for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
536 for(i=0;i<80;i++) printf("*");
540 //_____________________________________________________________________________
541 void AliITS::SetDefaults()
543 // sets the default segmentation, response, digit and raw cluster classes
546 AliITSDetType *iDetType;
550 AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(fITSgeom);
551 AliITSresponseSPD *resp0=new AliITSresponseSPD();
553 if (!iDetType->GetSegmentationModel()) SetSegmentationModel(0,seg0);
554 if (!iDetType->GetResponseModel()) SetResponseModel(0,resp0);
555 // set digit and raw cluster classes to be used
556 const char *kData=resp0->DataType();
557 if (strstr(kData,"real")) {
558 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSPD");
559 } else iDetType->ClassNames("AliITSdigitSPD","AliITSRawClusterSPD");
562 AliITSresponseSDD *resp1=new AliITSresponseSDD();
563 AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(fITSgeom,resp1);
565 if (!iDetType->GetSegmentationModel()) SetSegmentationModel(1,seg1);
566 if (!iDetType->GetResponseModel()) SetResponseModel(1,resp1);
567 kData=resp1->DataType();
568 Option_t *opt=resp1->ZeroSuppOption();
569 if ((!strstr(opt,"2D")) && (!strstr(opt,"1D")) || strstr(kData,"real") ) {
570 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
571 } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
574 AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(fITSgeom);
575 AliITSresponseSSD *resp2=new AliITSresponseSSD();
577 if (!iDetType->GetSegmentationModel()) SetSegmentationModel(2,seg2);
578 if (!iDetType->GetResponseModel()) SetResponseModel(2,resp2);
579 kData=resp2->DataType();
580 if (strstr(kData,"real")) {
581 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSSD");
582 } else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
585 Warning("SetDefaults","Only the three basic detector types are initialised!");
591 //_____________________________________________________________________________
593 void AliITS::MakeTreeC(Option_t *option)
595 // create a separate tree to store the clusters
597 char *optC = strstr(option,"C");
598 if (optC && !fTreeC) fTreeC = new TTree("TC","Clusters in ITS");
600 Int_t buffersize = 4000;
603 char *det[3] = {"SPD","SDD","SSD"};
605 // one branch for Clusters per type of detector
607 for(i=0; i<fNDetTypes ;i++) {
608 if (fNDetTypes==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
609 else sprintf(branchname,"%sClusters%d",GetName(),i+1);
610 if (fCtype && fTreeC) {
611 TreeC()->Branch(branchname,&((*fCtype)[i]), buffersize);
612 printf("Making Branch %s for Clusters of detector type %d\n",branchname,i+1);
618 //_____________________________________________________________________________
619 void AliITS::GetTreeC(Int_t event)
622 // get the clusters tree for this event and set the branch address
626 char *det[3] = {"SPD","SDD","SSD"};
633 sprintf(treeName,"TreeC%d",event);
634 fTreeC = (TTree*)gDirectory->Get(treeName);
640 for(i=0; i<fNDetTypes; i++) {
641 if (fNDetTypes==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
642 else sprintf(branchname,"%sClusters%d",GetName(),i+1);
644 branch = fTreeC->GetBranch(branchname);
645 if (branch) branch->SetAddress(&((*fCtype)[i]));
649 printf("ERROR: cannot find Clusters Tree for event:%d\n",event);
653 //_____________________________________________________________________________
654 void AliITS::MakeBranch(Option_t* option){
656 // Creates Tree branches for the ITS.
660 Int_t buffersize = 4000;
662 sprintf(branchname,"%s",GetName());
664 AliDetector::MakeBranch(option);
667 // one branch for digits per type of detector
669 char *det[3] = {"SPD","SDD","SSD"};
671 fNdtype = new Int_t[fNDetTypes];
672 fDtype = new TObjArray(fNDetTypes);
674 fNctype = new Int_t[fNDetTypes];
675 fCtype = new TObjArray(fNDetTypes);
681 for(i=0; i<fNDetTypes ;i++) {
682 AliITSDetType *iDetType=DetType(i);
683 iDetType->GetClassNames(kDigclass,kClclass);
684 //printf("i, digclass, recclass %d %s %s\n",i,kDigclass,kClclass);
686 (*fDtype)[i] = new TClonesArray(kDigclass,100);
689 (*fCtype)[i] = new TClonesArray(kClclass,100);
694 for(i=0; i<fNDetTypes ;i++) {
695 if (fNDetTypes==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
696 else sprintf(branchname,"%sDigits%d",GetName(),i+1);
698 if (fDtype && gAlice->TreeD()) {
699 gAlice->TreeD()->Branch(branchname,&((*fDtype)[i]), buffersize);
700 printf("Making Branch %s for digits of type %d\n",branchname,i+1);
704 // only one branch for rec points for all detector types
705 sprintf(branchname,"%sRecPoints",GetName());
707 fRecPoints=new TClonesArray("AliITSRecPoint",1000);
709 if (fRecPoints && gAlice->TreeR()) {
710 gAlice->TreeR()->Branch(branchname,&fRecPoints, buffersize);
711 printf("Making Branch %s for reconstructed space points\n",branchname);
717 //___________________________________________
718 void AliITS::SetTreeAddress()
721 // Set branch address for the Trees.
724 AliDetector::SetTreeAddress();
726 char *det[3] = {"SPD","SDD","SSD"};
729 TTree *treeD = gAlice->TreeD();
730 TTree *treeR = gAlice->TreeR();
734 for(i=0; i<fNDetTypes; i++) {
735 if (fNDetTypes==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
736 else sprintf(branchname,"%sDigits%d",GetName(),i+1);
738 branch = treeD->GetBranch(branchname);
739 if (branch) branch->SetAddress(&((*fDtype)[i]));
746 sprintf(branchname,"%sRecPoints",GetName());
748 branch = treeR->GetBranch(branchname);
749 if (branch) branch->SetAddress(&fRecPoints);
756 //____________________________________________________________________________
757 void AliITS::InitModules(Int_t size,Int_t &nmodules){
759 //initialize the modules array
761 Int_t nl,indexMAX,index;
764 if(size<=0){ // default to using data stored in AliITSgeom
766 printf("Error in AliITS::InitModule fITSgeom not defined\n");
768 } // end if fITSgeom==0
769 nl = fITSgeom->GetNlayers();
770 indexMAX = fITSgeom->GetModuleIndex(nl,fITSgeom->GetNladders(nl),
771 fITSgeom->GetNdetectors(nl))+1;
773 fITSmodules = new TObjArray(indexMAX);
774 indSPD = fITSgeom->GetModuleIndex(2,fITSgeom->GetNladders(2),
775 fITSgeom->GetNdetectors(2));
776 indSDD = fITSgeom->GetModuleIndex(4,fITSgeom->GetNladders(4),
777 fITSgeom->GetNdetectors(4));
778 for(index=0;index<indexMAX;index++){
780 fITSmodules->AddAt( new AliITSmodule(index),index);
781 else if(index<=indSDD)
782 fITSmodules->AddAt( new AliITSmodule(index),index);
784 fITSmodules->AddAt( new AliITSmodule(index),index);
787 fITSmodules = new TObjArray(size);
792 //____________________________________________________________________________
793 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t lastev,Int_t nmodules,Option_t *option,Text_t *filename){
795 // fill the modules with the sorted by module hits; add hits from background
798 static TTree *trH1; //Tree with background hits
799 static TClonesArray *fHits2; //List of hits for one track only
801 static Bool_t first=kTRUE;
803 char *add = strstr(option,"Add");
808 cout<<"filename "<<filename<<endl;
809 file=new TFile(filename);
810 cout<<"I have opened "<<filename<<" file "<<endl;
811 fHits2 = new TClonesArray("AliITShit",1000 );
816 // Get Hits Tree header from file
817 if(fHits2) fHits2->Clear();
818 if(trH1) delete trH1;
822 sprintf(treeName,"TreeH%d",bgrev);
823 trH1 = (TTree*)gDirectory->Get(treeName);
824 //printf("TrH1 %p of treename %s for event %d \n",trH1,treeName,bgrev);
827 printf("ERROR: cannot find Hits Tree for event:%d\n",bgrev);
829 // Set branch addresses
832 sprintf(branchname,"%s",GetName());
833 if (trH1 && fHits2) {
834 branch = trH1->GetBranch(branchname);
835 if (branch) branch->SetAddress(&fHits2);
839 //Int_t ntracks1 =(Int_t)TrH1->GetEntries();
840 //printf("background - ntracks1 - %d\n",ntracks1);
843 // store temporarily coordinates of signal particles - see where delete
844 static Int_t *signal;
846 if(!signal) signal=new Int_t[nmodules];
847 memset(signal,0,sizeof(int)*nmodules);
848 // Float_t xhit[nmodules][4];
849 Float_t **xhit = new Float_t*[nmodules];
850 for(i=0;i<nmodules;i++) xhit[i] = new Float_t[4];
851 // Float_t yhit[nmodules][4];
852 Float_t **yhit = new Float_t*[nmodules];
853 for(i=0;i<nmodules;i++) yhit[i] = new Float_t[4];
855 Int_t npart = gAlice->GetEvent(evnt);
857 TClonesArray *itsHits = this->Hits();
858 Int_t lay,lad,det,index;
862 TTree *iTH = gAlice->TreeH();
863 Int_t ntracks =(Int_t) iTH->GetEntries();
866 for(t=0; t<ntracks; t++){
869 Int_t nhits = itsHits->GetEntriesFast();
870 if (!nhits) continue;
871 // cout << nhits << " hits in track " << t << endl;
872 for(h=0; h<nhits; h++){
873 itsHit = (AliITShit *)itsHits->UncheckedAt(h);
874 itsHit->GetDetectorID(lay,lad,det);
875 index = fITSgeom->GetModuleIndex(lay,lad,det);
876 mod = this->GetModule(index);
878 xhit[index][signal[index]]=itsHit->fX;
879 yhit[index][signal[index]]=itsHit->fY;
881 if (signal[index] >4)
882 printf("index,nsignal %d %d\n",index,signal[index]);
884 if(lay == 1 || lay == 2)
885 mod->AddHit((AliITShit *) itsHit,t,h);
886 else if(lay == 3 || lay == 4)
887 mod->AddHit((AliITShit *) itsHit,t,h);
888 else if(lay == 5 || lay ==6)
889 mod->AddHit((AliITShit *)itsHit,t,h);
891 mod->AddHit(itsHit,t,h);
893 } // end loop over hits
894 } // end loop over tracks
896 // open the file with background
900 ntracks =(Int_t)trH1->GetEntries();
901 //printf("background - ntracks1 %d\n",ntracks);
902 //printf("background - Start loop over tracks \n");
905 for(track=0; track<ntracks; track++) {
907 if (fHits2) fHits2->Clear();
908 trH1->GetEvent(track);
910 for(i=0;i<fHits2->GetEntriesFast();++i) {
912 itsHit=(AliITShit*) (*fHits2)[i];
913 itsHit->GetDetectorID(lay,lad,det);
914 index = fITSgeom->GetModuleIndex(lay,lad,det);
915 mod = this->GetModule(index);
917 Float_t xbgr=itsHit->fX;
918 Float_t ybgr=itsHit->fY;
919 Float_t ebgr=itsHit->GetIonization();
922 for(isig =0; isig < signal[index]; isig++) {
924 (xbgr-xhit[index][isig])*(xbgr-xhit[index][isig])
925 +(ybgr-yhit[index][isig])*(ybgr-yhit[index][isig]);
926 if (dist<0.2&& ebgr!=0) cond=kTRUE; // check this number for ITS!
929 if(lay == 1 || lay == 2)
930 mod->AddHit((AliITShit *) itsHit,track,i);
931 else if(lay == 3 || lay == 4)
932 mod->AddHit((AliITShit *) itsHit,track,i);
933 else if(lay == 5 || lay ==6)
934 mod->AddHit((AliITShit *)itsHit,track,i);
936 mod->AddHit(itsHit,track,i);
939 } // end loop over hits
940 } // end loop over tracks
942 TTree *fAli=gAlice->TreeK();
945 if (fAli) fileAli =fAli->GetCurrentFile();
946 //printf("fAli, file %p %p\n",fAli,file);
952 for(i=0;i<nmodules;i++) delete[] xhit[i];
954 for(i=0;i<nmodules;i++) delete[] yhit[i];
956 //if (evnt==lastev) {delete [] signal; delete signal;}
958 //gObjectTable->Print();
963 //____________________________________________________________________________
964 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t lastev,Int_t size,
965 Option_t *option,Option_t *opt,Text_t *filename)
967 // keep galice.root for signal and name differently the file for
968 // background when add! otherwise the track info for signal will be lost !
970 char *all = strstr(opt,"All");
971 char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
972 //printf("Det 1 2 3 %s %s %s \n",det[0],det[1],det[2]);
975 InitModules(size,nmodules);
976 FillModules(evNumber,bgrev,lastev,nmodules,option,filename);
977 //printf("nmodules %d\n",nmodules);
980 AliITSsimulation* sim;
982 TObjArray *branches=gAlice->TreeD()->GetListOfBranches();
983 AliITSgeom *geom = GetITSgeom();
986 for(id=0;id<3;id++) {
987 //printf("id %d All %s det[id] %s \n",id,all,det[id]);
988 if (!all && !det[id]) continue;
989 branch = (TBranch*)branches->UncheckedAt(id);
990 AliITSDetType *iDetType=DetType(id);
991 sim = (AliITSsimulation*)iDetType->GetSimulationModel();
993 Error("HitsToDigits","The simulation class was not instantiated!");
995 // or SetDefaultSimulation(id,iDetType*);
997 Int_t first = geom->GetStartDet(id);
998 Int_t last = geom->GetLastDet(id);
999 //printf("det type %d first, last %d %d \n",id,first,last);
1000 for(module=first;module<=last;module++) {
1001 AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
1002 sim->DigitiseModule(mod,module,evNumber);
1003 // fills all branches - wasted disk space
1004 gAlice->TreeD()->Fill();
1006 // try and fill only the branch
1009 } // loop over modules
1010 } // loop over detector types
1015 //Int_t nentries=(Int_t)gAlice->TreeD()->GetEntries();
1018 sprintf(hname,"TreeD%d",evNumber);
1019 gAlice->TreeD()->Write(hname);
1021 gAlice->TreeD()->Reset();
1026 //____________________________________________________________________________
1027 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt)
1030 char *all = strstr(opt,"All");
1031 char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1033 static Bool_t first=kTRUE;
1041 AliITSClusterFinder* rec;
1043 TObjArray *branches=gAlice->TreeR()->GetListOfBranches();
1044 AliITSgeom *geom = GetITSgeom();
1047 for(id=0;id<fNDetTypes;id++) {
1048 if (!all && !det[id]) continue;
1049 branch = (TBranch*)branches->UncheckedAt(id);
1050 AliITSDetType *iDetType=DetType(id);
1051 rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1053 Error("DigitsToRecPoints","The cluster finder class was not instantiated!");
1055 // or SetDefaultClusterFinder(id,iDetType*);
1057 TClonesArray *itsDigits = this->DigitsAddress(id);
1059 Int_t first = geom->GetStartDet(id);
1060 Int_t last = geom->GetLastDet(id);
1061 //printf("det type %d first, last %d %d \n",id,first,last);
1062 for(module=first;module<=last;module++) {
1063 //printf("AliITS: module=%d\n",module);
1064 this->ResetDigits();
1065 if (all) gAlice->TreeD()->GetEvent(lastentry+module);
1066 else gAlice->TreeD()->GetEvent(lastentry+(module-first));
1067 Int_t ndigits = itsDigits->GetEntriesFast();
1068 if (ndigits) rec->FindRawClusters();
1069 gAlice->TreeR()->Fill();
1073 // try and fill only the branch
1075 //ResetRecPoints(id);
1076 } // loop over modules
1077 } // loop over detector types
1080 //Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1082 //Int_t ncentries=(Int_t)TC->GetEntries();
1085 sprintf(hname,"TreeR%d",evNumber);
1086 gAlice->TreeR()->Write(hname);
1088 gAlice->TreeR()->Reset();
1090 sprintf(hname,"TreeC%d",evNumber);
1096 //____________________________________________________________________________
1097 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t lastev,Int_t size,
1098 Option_t *option,Option_t *opt,Text_t *filename)
1100 // keep galice.root for signal and name differently the file for
1101 // background when add! otherwise the track info for signal will be lost !
1103 char *all = strstr(opt,"All");
1104 char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1107 InitModules(size,nmodules);
1108 FillModules(evNumber,bgrev,lastev,nmodules,option,filename);
1110 static AliITSsimulationFastPoints* sim=0;
1111 if (!sim) sim = new AliITSsimulationFastPoints();
1114 AliITSgeom *geom = GetITSgeom();
1117 for(id=0;id<3;id++) {
1118 if (!all && !det[id]) continue;
1119 Int_t first = geom->GetStartDet(id);
1120 Int_t last = geom->GetLastDet(id);
1121 for(module=first;module<=last;module++) {
1122 AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
1123 sim->CreateFastRecPoints(mod);
1124 gAlice->TreeR()->Fill();
1126 } // loop over modules
1127 } // loop over detector types
1132 //Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1135 sprintf(hname,"TreeR%d",evNumber);
1136 gAlice->TreeR()->Write(hname);
1138 gAlice->TreeR()->Reset();
1142 //____________________________________________________________________________
1143 void AliITS::Streamer(TBuffer &R__b){
1144 // Stream an object of class AliITS.
1148 AliITSresponse *response;
1149 AliITSsegmentation *seg;
1150 TClonesArray *digitsaddress;
1151 TClonesArray *clustaddr;
1154 if (R__b.IsReading()) {
1155 Version_t R__v = R__b.ReadVersion();
1157 AliDetector::Streamer(R__b);
1161 if(fIdSens!=0) delete[] fIdSens;
1162 if(fIdName!=0) delete[] fIdName;
1163 fIdSens = new Int_t[fIdN];
1164 fIdName = new char*[fIdN];
1165 for(i=0;i<fIdN;i++) R__b >> fIdSens[i];
1166 for(i=0;i<fIdN;i++){
1168 fIdName[i] = new char[l+1]; // add room for null character.
1169 for(j=0;j<l;j++) R__b >> fIdName[i][j];
1170 fIdName[i][l] = '\0'; // Null terminate this string.
1172 R__b >> fMajorVersion;
1173 R__b >> fMinorVersion;
1175 R__b.ReadArray(fNdtype);
1177 R__b.ReadArray(fNctype);
1181 R__b >> fNRecPoints;
1182 // stream out only response and segmentation
1183 for(i=0;i<fNDetTypes;i++) {
1184 det=(AliITSDetType*)(*fDetTypes)[i];
1185 det->Streamer(R__b);
1186 response=det->GetResponseModel();
1187 if (response) response->Streamer(R__b);
1188 seg=det->GetSegmentationModel();
1189 if (seg) seg->Streamer(R__b);
1190 digitsaddress=(TClonesArray*) (*fDtype)[i];
1191 digitsaddress->Streamer(R__b);
1192 clustaddr=(TClonesArray*) (*fCtype)[i];
1193 clustaddr->Streamer(R__b);
1198 R__b.WriteVersion(AliITS::IsA());
1199 AliDetector::Streamer(R__b);
1203 for(i=0;i<fIdN;i++) R__b <<fIdSens[i];
1204 for(i=0;i<fIdN;i++){
1205 l = strlen(fIdName[i]);
1207 for(j=0;j<l;j++) R__b << fIdName[i][j];
1209 R__b << fMajorVersion;
1210 R__b << fMinorVersion;
1212 R__b.WriteArray(fNdtype, fNDetTypes);
1214 R__b.WriteArray(fNctype, fNDetTypes);
1218 R__b << fNRecPoints;
1219 for(i=0;i<fNDetTypes;i++) {
1220 det=(AliITSDetType*)(*fDetTypes)[i];
1221 det->Streamer(R__b);
1222 response=det->GetResponseModel();
1223 if(response) response->Streamer(R__b);
1224 seg=det->GetSegmentationModel();
1225 if (seg) seg->Streamer(R__b);
1226 digitsaddress=(TClonesArray*) (*fDtype)[i];
1227 digitsaddress->Streamer(R__b);
1228 clustaddr=(TClonesArray*) (*fCtype)[i];
1229 clustaddr->Streamer(R__b);
1236 ClassImp(AliITSRecPoint)
1237 ClassImp(AliITSsegmentation)
1238 ClassImp(AliITSresponse)
1239 //ClassImp(AliITStrack)