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.9.2.3 2000/02/02 13:42:09 barbera
19 fixed AliITS.cxx for new AliRun structure. Added ITS hits list to list of hits which will have their track numbers updated
21 Revision 1.9.2.2 2000/01/23 03:03:13 nilsen
22 //fixed FillModule. Removed fi(fabs(xl)<dx....
24 Revision 1.9.2.1 2000/01/12 19:03:32 nilsen
25 This is the version of the files after the merging done in December 1999.
26 See the ReadMe110100.txt file for details
28 Revision 1.9 1999/11/14 14:33:25 fca
29 Correct problems with distructors and pointers, thanks to I.Hrivnacova
31 Revision 1.8 1999/09/29 09:24:19 fca
32 Introduction of the Copyright and cvs Log
36 ///////////////////////////////////////////////////////////////////////////////
38 // An overview of the basic philosophy of the ITS code development
39 // and analysis is show in the figure below.
42 <img src="picts/ITS/ITS_Analysis_schema.gif">
45 <font size=+2 color=red>
46 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
47 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
53 // AliITS. Inner Traking System base class.
54 // This class contains the base procedures for the Inner Tracking System
58 <img src="picts/ITS/AliITS_Class_Diagram.gif">
61 <font size=+2 color=red>
62 <p>This show the class diagram of the different elements that are part of
70 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
73 // Modified and documented by Bjorn S. Nilsen
77 // Modified and documented by A. Bologna
80 // AliITS is the general base class for the ITS. Also see AliDetector for
81 // futher information.
83 ///////////////////////////////////////////////////////////////////////////////
88 #include <TObjArray.h>
90 #include <TObjectTable.h>
96 #include "AliITSMap.h"
97 #include "AliITSDetType.h"
98 #include "AliITSClusterFinder.h"
99 #include "AliITSsimulation.h"
100 #include "AliITSsegmentationSPD.h"
101 #include "AliITSresponseSPD.h"
102 #include "AliITSsegmentationSDD.h"
103 #include "AliITSresponseSDD.h"
104 #include "AliITSsegmentationSSD.h"
105 #include "AliITSresponseSSD.h"
111 //_____________________________________________________________________________
112 AliITS::AliITS() : AliDetector() {
114 // Default initialiser for ITS
115 // The default constructor of the AliITS class. In addition to
116 // creating the AliITS class it zeros the variables fIshunt (a member
117 // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
118 // fITSpoints, fIdSens, and fIdName. The AliDetector default constructor
125 //fNDetTypes = fgkNTYPES;
144 //_____________________________________________________________________________
145 AliITS::AliITS(const char *name, const char *title):AliDetector(name,title){
147 // Default initialiser for ITS
148 // The constructor of the AliITS class. In addition to creating the
149 // AliITS class, it allocates memory for the TClonesArrays fHits and
150 // fDigits, and for the TObjArray fITSpoints. It also zeros the variables
151 // fIshunt (a member of AliDetector class), fEuclidOut, and fIdN, and zeros
152 // the pointers fIdSens and fIdName. To help in displaying hits via the ROOT
153 // macro display.C AliITS also sets the marker color to red. The variables
154 // passes with this constructor, const char *name and *title, are used by
155 // the constructor of AliDetector class. See AliDetector class for a
156 // description of these parameters and its constructor functions.
159 fHits = new TClonesArray("AliITShit", 1560);
160 gAlice->AddHitList(fHits);
162 //fNDetTypes = fgkNTYPES;
164 fNdtype = new Int_t[fgkNTYPES];
165 fDtype = new TObjArray(fgkNTYPES);
167 fNctype = new Int_t[fgkNTYPES];
168 fCtype = new TObjArray(fgkNTYPES);
184 fDetTypes = new TObjArray(fgkNTYPES);
187 for(i=0;i<fgkNTYPES;i++) {
188 (*fDetTypes)[i]=new AliITSDetType();
194 SetMarkerColor(kRed);
198 //___________________________________________________________________________
199 AliITS::AliITS(AliITS &source){
201 if(this==&source) return;
202 printf("Error: You are not allowed to make a copy of the AliITS\n");
205 //____________________________________________________________________________
206 AliITS& AliITS::operator=(AliITS &source){
207 // assignment operator
208 if(this==&source) return *this;
209 printf("Error: You are not allowed to make a copy of the AliITS\n");
211 return *this; //fake return
213 //____________________________________________________________________________
214 void AliITS::ClearModules(){
216 //clear the modules TObjArray
220 Int_t indSPD = fITSgeom->GetModuleIndex(2,fITSgeom->GetNladders(2),
221 fITSgeom->GetNdetectors(2));
222 Int_t indSDD = fITSgeom->GetModuleIndex(4,fITSgeom->GetNladders(4),
223 fITSgeom->GetNdetectors(4));
224 for(i=0;i<fITSmodules->GetEntriesFast();i++){
225 if(!fITSmodules->At(i)) continue;
227 delete (AliITSmodule *) fITSmodules->At(i);
229 delete (AliITSmodule *) fITSmodules->At(i);
231 delete (AliITSmodule *) fITSmodules->At(i);
233 }// end if fITSmodules!=0
236 //_____________________________________________________________________________
239 // Default distructor for ITS
240 // The default destructor of the AliITS class. In addition to deleting
241 // the AliITS class it deletes the memory pointed to by the fHits, fDigits,
242 // fIdSens, fIdName, and fITSpoints.
249 if(fIdName!=0) delete[] fIdName;
250 if(fIdSens!=0) delete[] fIdSens;
252 this->ClearModules();
254 }// end if fITSmodules!=0
259 for (i=0;i<fgkNTYPES;i++) {
265 for (i=0;i<fgkNTYPES;i++) {
277 if (fTreeC) delete fTreeC;
281 //___________________________________________
282 AliITSDetType* AliITS::DetType(Int_t id)
284 //return pointer to id detector type
285 return ((AliITSDetType*) (*fDetTypes)[id]);
288 //___________________________________________
289 void AliITS::SetClasses(Int_t id, const char *digit, const char *cluster)
291 //set the digit and cluster classes to be used for the id detector type
292 ((AliITSDetType*) (*fDetTypes)[id])->ClassNames(digit,cluster);
295 //___________________________________________
296 void AliITS::SetResponseModel(Int_t id, AliITSresponse *response)
298 //set the response model for the id detector type
300 ((AliITSDetType*) (*fDetTypes)[id])->ResponseModel(response);
304 //___________________________________________
305 void AliITS::SetSegmentationModel(Int_t id, AliITSsegmentation *seg)
307 //set the segmentation model for the id detector type
309 ((AliITSDetType*) (*fDetTypes)[id])->SegmentationModel(seg);
313 //___________________________________________
314 void AliITS::SetSimulationModel(Int_t id, AliITSsimulation *sim)
316 //set the simulation model for the id detector type
318 ((AliITSDetType*) (*fDetTypes)[id])->SimulationModel(sim);
321 //___________________________________________
322 void AliITS::SetReconstructionModel(Int_t id, AliITSClusterFinder *reconst)
324 //set the cluster finder model for the id detector type
326 ((AliITSDetType*) (*fDetTypes)[id])->ReconstructionModel(reconst);
330 //_____________________________________________________________________________
331 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
334 // The function to add information to the AliITShit class. See the
335 // AliITShit class for a full description. This function allocates the
336 // necessary new space for the hit information and passes the variable
337 // track, and the pointers *vol and *hits to the AliITShit constructor
340 TClonesArray &lhits = *fHits;
341 new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
343 //_____________________________________________________________________________
344 void AliITS::AddRealDigit(Int_t id, Int_t *digits)
346 // add a real digit - as coming from data
348 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
349 new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
352 //_____________________________________________________________________________
353 void AliITS::AddSimDigit(Int_t id, AliITSdigit *d)
356 // add a simulated digit
358 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
363 new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
366 new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
369 new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
375 //_____________________________________________________________________________
376 void AliITS::AddSimDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,Int_t *hits,Float_t *charges){
378 // add a simulated digit to the list
380 TClonesArray &ldigits = *((TClonesArray*)(*fDtype)[id]);
384 new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks,hits);
387 new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,hits,charges);
390 new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks,hits);
396 //_____________________________________________________________________________
397 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c)
400 // add a cluster to the list
402 TClonesArray &lcl = *((TClonesArray*)(*fCtype)[id]);
407 new(lcl[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
410 new(lcl[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
413 new(lcl[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
420 //_____________________________________________________________________________
421 void AliITS::AddRecPoint(const AliITSRecPoint &r)
424 // Add a reconstructed space point to the list
426 TClonesArray &lrecp = *fRecPoints;
427 new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
431 //____________________________________________
432 void AliITS::ResetDigits()
435 // Reset number of digits and the digits array for thE ITS detector
441 for (i=0;i<fgkNTYPES;i++ ) {
442 if ((*fDtype)[i]) ((TClonesArray*)(*fDtype)[i])->Clear();
443 if (fNdtype) fNdtype[i]=0;
447 //____________________________________________
448 void AliITS::ResetDigits(Int_t i)
451 // Reset number of digits and the digits array for this branch
453 if ((*fDtype)[i]) ((TClonesArray*)(*fDtype)[i])->Clear();
454 if (fNdtype) fNdtype[i]=0;
458 //____________________________________________
459 void AliITS::ResetClusters()
462 // Reset number of clusters and the clusters array for ITS
466 for (i=0;i<fgkNTYPES;i++ ) {
467 if ((*fCtype)[i]) ((TClonesArray*)(*fCtype)[i])->Clear();
468 if (fNctype) fNctype[i]=0;
473 //____________________________________________
474 void AliITS::ResetClusters(Int_t i)
477 // Reset number of clusters and the clusters array for this branch
479 if ((*fCtype)[i]) ((TClonesArray*)(*fCtype)[i])->Clear();
480 if (fNctype) fNctype[i]=0;
485 //____________________________________________
486 void AliITS::ResetRecPoints()
489 // Reset number of rec points and the rec points array
492 if (fRecPoints) fRecPoints->Clear();
496 //_____________________________________________________________________________
497 Int_t AliITS::DistancetoPrimitive(Int_t , Int_t ){
499 // Distance from mouse to ITS on the screen. Dummy routine
500 // A dummy routine used by the ROOT macro display.C to allow for the
501 // use of the mouse (pointing device) in the macro. In general this should
502 // never be called. If it is it returns the number 9999 for any value of
508 //_____________________________________________________________________________
511 // Initialise ITS after it has been built
512 // This routine initializes the AliITS class. It is intended to be called
513 // from the Init function in AliITSv?. Besides displaying a banner
514 // indicating that it has been called it initializes the array fIdSens
515 // and sets the default segmentation, response, digit and raw cluster classes
516 // Therefore it should be called after a call to CreateGeometry.
524 for(i=0;i<35;i++) printf("*");
525 printf(" ITS_INIT ");
526 for(i=0;i<35;i++) printf("*");
530 for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
532 for(i=0;i<80;i++) printf("*");
536 //_____________________________________________________________________________
537 void AliITS::SetDefaults()
539 // sets the default segmentation, response, digit and raw cluster classes
542 AliITSDetType *iDetType;
546 AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(fITSgeom);
547 AliITSresponseSPD *resp0=new AliITSresponseSPD();
549 if (!iDetType->GetSegmentationModel()) SetSegmentationModel(0,seg0);
550 if (!iDetType->GetResponseModel()) SetResponseModel(0,resp0);
551 // set digit and raw cluster classes to be used
552 const char *kData0=resp0->DataType();
553 if (strstr(kData0,"real")) {
554 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSPD");
555 } else iDetType->ClassNames("AliITSdigitSPD","AliITSRawClusterSPD");
558 AliITSresponseSDD *resp1=new AliITSresponseSDD();
559 AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(fITSgeom,resp1);
561 if (!iDetType->GetSegmentationModel()) SetSegmentationModel(1,seg1);
562 if (!iDetType->GetResponseModel()) SetResponseModel(1,resp1);
563 const char *kData1=resp1->DataType();
564 const char *kopt=resp1->ZeroSuppOption();
565 if ((!strstr(kopt,"2D")) && (!strstr(kopt,"1D")) || strstr(kData1,"real") ) {
566 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
567 } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
570 AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(fITSgeom);
571 AliITSresponseSSD *resp2=new AliITSresponseSSD();
573 if (!iDetType->GetSegmentationModel()) SetSegmentationModel(2,seg2);
574 if (!iDetType->GetResponseModel()) SetResponseModel(2,resp2);
575 const char *kData2=resp2->DataType();
576 if (strstr(kData2,"real")) {
577 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSSD");
578 } else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
581 Warning("SetDefaults","Only the three basic detector types are initialised!");
587 //_____________________________________________________________________________
588 void AliITS::SetDefaultSimulation()
593 //_____________________________________________________________________________
594 void AliITS::SetDefaultClusterFinders()
599 //_____________________________________________________________________________
601 void AliITS::MakeTreeC(Option_t *option)
603 // create a separate tree to store the clusters
605 char *optC = strstr(option,"C");
606 if (optC && !fTreeC) fTreeC = new TTree("TC","Clusters in ITS");
608 Int_t buffersize = 4000;
611 char *det[3] = {"SPD","SDD","SSD"};
613 // one branch for Clusters per type of detector
615 for (i=0; i<fgkNTYPES ;i++) {
616 if (fgkNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
617 else sprintf(branchname,"%sClusters%d",GetName(),i+1);
618 if (fCtype && fTreeC) {
619 TreeC()->Branch(branchname,&((*fCtype)[i]), buffersize);
620 printf("Making Branch %s for Clusters of detector type %d\n",branchname,i+1);
626 //_____________________________________________________________________________
627 void AliITS::GetTreeC(Int_t event)
630 // get the clusters tree for this event and set the branch address
634 char *det[3] = {"SPD","SDD","SSD"};
641 sprintf(treeName,"TreeC%d",event);
642 fTreeC = (TTree*)gDirectory->Get(treeName);
648 for (i=0; i<fgkNTYPES; i++) {
649 if (fgkNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
650 else sprintf(branchname,"%sClusters%d",GetName(),i+1);
652 branch = fTreeC->GetBranch(branchname);
653 if (branch) branch->SetAddress(&((*fCtype)[i]));
657 printf("ERROR: cannot find Clusters Tree for event:%d\n",event);
661 //_____________________________________________________________________________
662 void AliITS::MakeBranch(Option_t* option){
664 // Creates Tree branches for the ITS.
668 Int_t buffersize = 4000;
670 sprintf(branchname,"%s",GetName());
672 AliDetector::MakeBranch(option);
675 // one branch for digits per type of detector
677 char *det[3] = {"SPD","SDD","SSD"};
683 for (i=0; i<fgkNTYPES ;i++) {
684 AliITSDetType *iDetType=DetType(i);
685 iDetType->GetClassNames(digclass,clclass);
686 //printf("i, digclass, recclass %d %s %s\n",i,digclass,clclass);
688 (*fDtype)[i] = new TClonesArray(digclass,10000);
690 (*fCtype)[i] = new TClonesArray(clclass,10000);
694 for (i=0; i<fgkNTYPES ;i++) {
695 if (fgkNTYPES==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",10000);
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<fgkNTYPES; i++) {
735 if (fgkNTYPES==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
762 //this->ClearModules();
766 Int_t nl,indexMAX,index;
769 if(size<=0){ // default to using data stored in AliITSgeom
771 printf("Error in AliITS::InitModule fITSgeom not defined\n");
773 } // end if fITSgeom==0
774 nl = fITSgeom->GetNlayers();
775 indexMAX = fITSgeom->GetModuleIndex(nl,fITSgeom->GetNladders(nl),
776 fITSgeom->GetNdetectors(nl))+1;
778 fITSmodules = new TObjArray(indexMAX);
779 indSPD = fITSgeom->GetModuleIndex(2,fITSgeom->GetNladders(2),
780 fITSgeom->GetNdetectors(2));
781 indSDD = fITSgeom->GetModuleIndex(4,fITSgeom->GetNladders(4),
782 fITSgeom->GetNdetectors(4));
783 for(index=0;index<indexMAX;index++){
785 fITSmodules->AddAt( new AliITSmodule(index),index);
786 else if(index<=indSDD)
787 fITSmodules->AddAt( new AliITSmodule(index),index);
789 fITSmodules->AddAt( new AliITSmodule(index),index);
792 fITSmodules = new TObjArray(size);
797 //____________________________________________________________________________
798 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t nmodules,Option_t *option,Text_t *filename){
800 // fill the modules with the sorted by module hits; add hits from background
803 static TTree *trH1; //Tree with background hits
804 static TClonesArray *fHits2; //List of hits for one track only
806 static Bool_t first=kTRUE;
808 char *addBgr = strstr(option,"Add");
813 cout<<"filename "<<filename<<endl;
814 file=new TFile(filename);
815 cout<<"I have opened "<<filename<<" file "<<endl;
816 fHits2 = new TClonesArray("AliITShit",1000 );
821 // Get Hits Tree header from file
822 if(fHits2) fHits2->Clear();
823 if(trH1) delete trH1;
827 sprintf(treeName,"TreeH%d",bgrev);
828 trH1 = (TTree*)gDirectory->Get(treeName);
829 //printf("TrH1 %p of treename %s for event %d \n",trH1,treeName,bgrev);
832 printf("ERROR: cannot find Hits Tree for event:%d\n",bgrev);
834 // Set branch addresses
837 sprintf(branchname,"%s",GetName());
838 if (trH1 && fHits2) {
839 branch = trH1->GetBranch(branchname);
840 if (branch) branch->SetAddress(&fHits2);
844 //Int_t ntracks1 =(Int_t)TrH1->GetEntries();
845 //printf("background - ntracks1 - %d\n",ntracks1);
848 Int_t npart = gAlice->GetEvent(evnt);
850 TClonesArray *itsHits = this->Hits();
851 Int_t lay,lad,det,index;
855 TTree *iTH = gAlice->TreeH();
856 Int_t ntracks =(Int_t) iTH->GetEntries();
859 for(t=0; t<ntracks; t++){
862 Int_t nhits = itsHits->GetEntriesFast();
863 if (!nhits) continue;
864 for(h=0; h<nhits; h++){
865 itsHit = (AliITShit *)itsHits->UncheckedAt(h);
866 itsHit->GetDetectorID(lay,lad,det);
867 index = fITSgeom->GetModuleIndex(lay,lad,det);
868 mod = this->GetModule(index);
869 if(lay == 1 || lay == 2)
870 mod->AddHit((AliITShit *) itsHit,t,h);
871 else if(lay == 3 || lay == 4)
872 mod->AddHit((AliITShit *) itsHit,t,h);
873 else if(lay == 5 || lay ==6)
874 mod->AddHit((AliITShit *)itsHit,t,h);
876 mod->AddHit(itsHit,t,h);
878 } // end loop over hits
879 } // end loop over tracks
881 // open the file with background
885 ntracks =(Int_t)trH1->GetEntries();
886 //printf("background - ntracks1 %d\n",ntracks);
887 //printf("background - Start loop over tracks \n");
890 for (track=0; track<ntracks; track++) {
892 if (fHits2) fHits2->Clear();
893 trH1->GetEvent(track);
895 for(i=0;i<fHits2->GetEntriesFast();++i) {
897 itsHit=(AliITShit*) (*fHits2)[i];
898 itsHit->GetDetectorID(lay,lad,det);
899 index = fITSgeom->GetModuleIndex(lay,lad,det);
900 mod = this->GetModule(index);
902 if(lay == 1 || lay == 2)
903 mod->AddHit((AliITShit *) itsHit,track,i);
904 else if(lay == 3 || lay == 4)
905 mod->AddHit((AliITShit *) itsHit,track,i);
906 else if(lay == 5 || lay ==6)
907 mod->AddHit((AliITShit *)itsHit,track,i);
909 mod->AddHit(itsHit,track,i);
912 } // end loop over hits
913 } // end loop over tracks
915 TTree *fAli=gAlice->TreeK();
918 if (fAli) fileAli =fAli->GetCurrentFile();
923 //gObjectTable->Print();
928 //____________________________________________________________________________
929 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size, Option_t *option, Option_t *opt,Text_t *filename)
931 // keep galice.root for signal and name differently the file for
932 // background when add! otherwise the track info for signal will be lost !
934 char *all = strstr(opt,"All");
935 char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
938 InitModules(size,nmodules);
939 FillModules(evNumber,bgrev,nmodules,option,filename);
942 AliITSsimulation* sim;
944 TObjArray *branches=gAlice->TreeD()->GetListOfBranches();
945 AliITSgeom *geom = GetITSgeom();
948 for (id=0;id<fgkNTYPES;id++) {
949 if (!all && !det[id]) continue;
950 branch = (TBranch*)branches->UncheckedAt(id);
951 AliITSDetType *iDetType=DetType(id);
952 sim = (AliITSsimulation*)iDetType->GetSimulationModel();
954 Error("HitsToDigits","The simulation class was not instantiated!");
956 // or SetDefaultSimulation();
958 Int_t first = geom->GetStartDet(id);
959 Int_t last = geom->GetLastDet(id);
960 printf("det type %d first, last %d %d \n",id,first,last);
961 for(module=first;module<=last;module++) {
962 AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
963 sim->DigitiseModule(mod,module,evNumber);
964 // fills all branches - wasted disk space
965 gAlice->TreeD()->Fill();
967 // try and fill only the branch
970 } // loop over modules
971 } // loop over detector types
975 Int_t nentries=(Int_t)gAlice->TreeD()->GetEntries();
976 printf("nentries in TreeD %d\n",nentries);
979 sprintf(hname,"TreeD%d",evNumber);
980 gAlice->TreeD()->Write(hname);
982 gAlice->TreeD()->Reset();
987 //____________________________________________________________________________
988 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt)
990 // cluster finding and reconstruction of space points
992 char *all = strstr(opt,"All");
993 char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
995 static Bool_t first=kTRUE;
1003 AliITSClusterFinder* rec;
1005 TObjArray *branches=gAlice->TreeR()->GetListOfBranches();
1006 AliITSgeom *geom = GetITSgeom();
1009 for (id=0;id<fgkNTYPES;id++) {
1010 if (!all && !det[id]) continue;
1011 branch = (TBranch*)branches->UncheckedAt(id);
1012 AliITSDetType *iDetType=DetType(id);
1013 rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1015 Error("DigitsToRecPoints","The cluster finder class was not instantiated!");
1017 // or SetDefaultClusterFinders();
1019 TClonesArray *itsDigits = this->DigitsAddress(id);
1021 Int_t first = geom->GetStartDet(id);
1022 Int_t last = geom->GetLastDet(id);
1023 for(module=first;module<=last;module++) {
1024 this->ResetDigits();
1025 if (all) gAlice->TreeD()->GetEvent(lastentry+module);
1026 else gAlice->TreeD()->GetEvent(lastentry+(module-first));
1027 Int_t ndigits = itsDigits->GetEntriesFast();
1028 if (ndigits) rec->FindRawClusters();
1029 gAlice->TreeR()->Fill();
1033 // try and fill only the branch
1035 //ResetRecPoints(id);
1036 } // loop over modules
1037 } // loop over detector types
1040 Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1041 Int_t ncentries=(Int_t)iTC->GetEntries();
1042 printf(" nentries ncentries %d %d\n", nentries, ncentries);
1045 sprintf(hname,"TreeR%d",evNumber);
1046 gAlice->TreeR()->Write(hname);
1048 gAlice->TreeR()->Reset();
1050 sprintf(hname,"TreeC%d",evNumber);
1056 //____________________________________________________________________________
1057 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
1058 Option_t *option,Option_t *opt,Text_t *filename)
1060 // keep galice.root for signal and name differently the file for
1061 // background when add! otherwise the track info for signal will be lost !
1063 char *all = strstr(opt,"All");
1064 char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),strstr(opt,"SSD")};
1067 InitModules(size,nmodules);
1068 FillModules(evNumber,bgrev,nmodules,option,filename);
1071 AliITSsimulation* sim;
1072 AliITSgeom *geom = GetITSgeom();
1074 TRandom *random=new TRandom[9];
1075 random[0].SetSeed(111);
1076 random[1].SetSeed(222);
1077 random[2].SetSeed(333);
1078 random[3].SetSeed(444);
1079 random[4].SetSeed(555);
1080 random[5].SetSeed(666);
1081 random[6].SetSeed(777);
1082 random[7].SetSeed(888);
1083 random[8].SetSeed(999);
1087 for (id=0;id<fgkNTYPES;id++) {
1088 if (!all && !det[id]) continue;
1089 AliITSDetType *iDetType=DetType(id);
1090 sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1092 Error("HitsToFastPoints","The simulation class was not instantiated!");
1094 // or SetDefaultSimulation();
1096 Int_t first = geom->GetStartDet(id);
1097 Int_t last = geom->GetLastDet(id);
1098 for(module=first;module<=last;module++) {
1099 AliITSmodule *mod = (AliITSmodule *)fITSmodules->At(module);
1100 sim->CreateFastRecPoints(mod,module,random);
1101 gAlice->TreeR()->Fill();
1103 } // loop over modules
1104 } // loop over detector types
1109 //Int_t nentries=(Int_t)gAlice->TreeR()->GetEntries();
1112 sprintf(hname,"TreeR%d",evNumber);
1113 gAlice->TreeR()->Write(hname);
1115 gAlice->TreeR()->Reset();
1121 //____________________________________________________________________________
1122 void AliITS::Streamer(TBuffer &R__b){
1123 // Stream an object of class AliITS.
1127 if (R__b.IsReading()) {
1128 Version_t R__v = R__b.ReadVersion();
1130 AliDetector::Streamer(R__b);
1132 R__b.ReadArray(fIdSens);
1133 for(i=0;i<fIdN;i++) fIdName[i].Streamer(R__b);
1135 R__b >> fITSmodules;
1137 R__b >> fMajorVersion;
1138 R__b >> fMinorVersion;
1142 fNdtype = new Int_t[fgkNTYPES];
1143 R__b.ReadFastArray(fNdtype,fgkNTYPES);
1146 fNctype = new Int_t[fgkNTYPES];
1147 R__b.ReadFastArray(fNctype,fgkNTYPES);
1149 R__b >> fNRecPoints;
1153 R__b.WriteVersion(AliITS::IsA());
1154 AliDetector::Streamer(R__b);
1156 R__b.WriteArray(fIdSens,fIdN);
1157 for(i=0;i<fIdN;i++) fIdName[i].Streamer(R__b);
1159 R__b << fITSmodules;
1161 R__b << fMajorVersion;
1162 R__b << fMinorVersion;
1165 R__b.WriteFastArray(fNdtype,fgkNTYPES);
1167 R__b.WriteFastArray(fNctype,fgkNTYPES);
1169 R__b << fNRecPoints;