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 **************************************************************************/
19 ///////////////////////////////////////////////////////////////////////////////
21 // An overview of the basic philosophy of the ITS code development
22 // and analysis is show in the figure below.
25 <img src="picts/ITS/ITS_Analysis_schema.gif">
28 <font size=+2 color=red>
29 <p>Roberto Barbera is in charge of the ITS Offline code (1999).
30 <a href="mailto:roberto.barbera@ct.infn.it">Roberto Barbera</a>.
36 // AliITS. Inner Traking System base class.
37 // This class contains the base procedures for the Inner Tracking System
41 <img src="picts/ITS/AliITS_Class_Diagram.gif">
44 <font size=+2 color=red>
45 <p>This show the class diagram of the different elements that are part of
53 // Written by Rene Brun, Federico Carminati, and Roberto Barbera
56 // Modified and documented by Bjorn S. Nilsen
60 // Modified and documented by A. Bologna
63 // AliITS is the general base class for the ITS. Also see AliDetector for
64 // futher information.
66 ///////////////////////////////////////////////////////////////////////////////
68 #include <Riostream.h>
72 #include <TClonesArray.h>
75 #include <TObjectTable.h>
81 #include <TVirtualMC.h>
83 #include "AliConfig.h"
84 #include "AliHeader.h"
86 #include "AliITSClusterFinderSDD.h"
87 #include "AliITSClusterFinderSPD.h"
88 #include "AliITSClusterFinderSSD.h"
89 #include "AliITSDetType.h"
90 #include "AliITSLoader.h"
91 #include "AliITSRawClusterSPD.h"
92 #include "AliITSRawClusterSDD.h"
93 #include "AliITSRawClusterSSD.h"
94 #include "AliITSRecPoint.h"
95 #include "AliITSdigit.h"
96 #include "AliITSgeom.h"
97 #include "AliITShit.h"
98 #include "AliITSmodule.h"
99 #include "AliITSpList.h"
100 #include "AliITSresponseSDD.h"
101 #include "AliITSresponseSPD.h"
102 #include "AliITSresponseSSD.h"
103 #include "AliITSsegmentationSDD.h"
104 #include "AliITSsegmentationSPD.h"
105 #include "AliITSsegmentationSSD.h"
106 #include "AliITSsimulationSDD.h"
107 #include "AliITSsimulationSPD.h"
108 #include "AliITSsimulationSSD.h"
110 #include "AliITSDigitizer.h"
111 #include "AliITSclustererV2.h"
112 #include "AliITStrackerV2.h"
113 #include "AliITSpidESD.h"
114 #include "AliV0vertexer.h"
115 #include "AliCascadeVertexer.h"
120 //______________________________________________________________________
121 AliITS::AliITS() : AliDetector() {
122 // Default initializer for ITS
123 // The default constructor of the AliITS class. In addition to
124 // creating the AliITS class it zeros the variables fIshunt (a member
125 // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
126 // fITSpoints, fIdSens, and fIdName. The AliDetector default constructor
135 fIshunt = 0; // not zeroed in AliDetector.
142 // SetDetectors(); // default to fOpt="All". This variable not written out.
148 fNDetTypes = kNTYPES;
163 SetMarkerColor(kRed);
165 //______________________________________________________________________
166 AliITS::AliITS(const char *name, const char *title):AliDetector(name,title){
167 // The standard Constructor for the ITS class. In addition to
168 // creating the AliITS class, it allocates memory for the TClonesArrays
169 // fHits, fSDigits, fDigits, fITSpoints, and the TObjArray of fCtype
170 // (clusters). It also zeros the variables
171 // fIshunt (a member of AliDetector class), fEuclidOut, and fIdN, and zeros
172 // the pointers fIdSens and fIdName. To help in displaying hits via the
173 // ROOT macro display.C AliITS also sets the marker color to red. The
174 // variables passes with this constructor, const char *name and *title,
175 // are used by the constructor of AliDetector class. See AliDetector
176 // class for a description of these parameters and its constructor
179 // const char *name Detector name. Should always be "ITS"
180 // const char *title Detector title.
186 fIshunt = 0; // not zeroed in AliDetector
187 fHits = new TClonesArray("AliITShit", 1560);//not done in AliDetector
188 gAlice->GetMCApp()->AddHitList(fHits); // Not done in AliDetector.
193 SetDetectors(); // default to fOpt="All". This variable not written out.
199 fNDetTypes = kNTYPES;
200 fDetTypes = new TObjArray(fNDetTypes);
202 fSDigits = new TClonesArray("AliITSpListItem",1000);
205 fNdtype = new Int_t[fNDetTypes];
206 fDtype = new TObjArray(fNDetTypes);
208 fCtype = new TObjArray(fNDetTypes);
209 fNctype = new Int_t[fNDetTypes];
211 fRecPoints = new TClonesArray("AliITSRecPoint",1000);
215 for(i=0;i<fNDetTypes;i++) {
216 fDetTypes->AddAt(new AliITSDetType(),i);
221 SetMarkerColor(kRed);
224 //______________________________________________________________________
226 // Default destructor for ITS.
227 // The default destructor of the AliITS class. In addition to deleting
228 // the AliITS class it deletes the memory pointed to by the fHits, fDigits,
229 // fSDigits, fCtype, fITSmodules, fITSgeom, fRecPoints, fIdSens, fIdName,
230 // fITSpoints, fDetType and it's contents.
254 fRecPoints->Delete();
258 delete[] fIdName; // Array of TStrings
261 this->ClearModules();
264 }// end if fITSmodules!=0
282 } // end if fDetTypes
285 if (fITSgeom) delete fITSgeom;
287 //______________________________________________________________________
288 AliITS::AliITS(const AliITS &source) : AliDetector(source){
289 // Copy constructor. This is a function which is not allowed to be
290 // done to the ITS. It exits with an error.
292 // AliITS &source An AliITS class.
298 if(this==&source) return;
299 Error("Copy constructor",
300 "You are not allowed to make a copy of the AliITS");
303 //______________________________________________________________________
304 AliITS& AliITS::operator=(AliITS &source){
305 // Assignment operator. This is a function which is not allowed to be
306 // done to the ITS. It exits with an error.
308 // AliITS &source An AliITS class.
314 if(this==&source) return *this;
315 Error("operator=","You are not allowed to make a copy of the AliITS");
317 return *this; //fake return
319 //______________________________________________________________________
320 Int_t AliITS::DistancetoPrimitive(Int_t,Int_t) const{
321 // Distance from mouse to ITS on the screen. Dummy routine
322 // A dummy routine used by the ROOT macro display.C to allow for the
323 // use of the mouse (pointing device) in the macro. In general this should
324 // never be called. If it is it returns the number 9999 for any value of
327 // Int_t Dummy screen coordinate.
328 // Int_t Dummy screen coordinate.
332 // Int_t Dummy = 9999 distance to ITS.
336 //______________________________________________________________________
338 // Initializer ITS after it has been built
339 // This routine initializes the AliITS class. It is intended to be
340 // called from the Init function in AliITSv?. Besides displaying a banner
341 // indicating that it has been called it initializes the array fIdSens
342 // and sets the default segmentation, response, digit and raw cluster
343 // classes therefore it should be called after a call to CreateGeometry.
354 for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
356 //______________________________________________________________________
357 void AliITS::SetDefaults(){
358 // sets the default segmentation, response, digit and raw cluster classes.
366 if(fDebug) printf("%s: SetDefaults\n",ClassName());
368 AliITSDetType *iDetType;
372 if (!iDetType->GetSegmentationModel()) {
373 AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(fITSgeom);
374 SetSegmentationModel(0,seg0);
376 if (!iDetType->GetResponseModel()) {
377 SetResponseModel(0,new AliITSresponseSPD());
379 // set digit and raw cluster classes to be used
381 const char *kData0=(iDetType->GetResponseModel())->DataType();
382 if (strstr(kData0,"real")) {
383 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSPD");
384 } else iDetType->ClassNames("AliITSdigitSPD","AliITSRawClusterSPD");
388 if (!iDetType->GetResponseModel()) {
389 SetResponseModel(1,new AliITSresponseSDD("simulated"));
391 AliITSresponse *resp1=iDetType->GetResponseModel();
392 if (!iDetType->GetSegmentationModel()) {
393 AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(fITSgeom,resp1);
394 SetSegmentationModel(1,seg1);
396 const char *kData1=(iDetType->GetResponseModel())->DataType();
397 const char *kopt=iDetType->GetResponseModel()->ZeroSuppOption();
398 if((!strstr(kopt,"2D"))&&(!strstr(kopt,"1D")) || strstr(kData1,"real") ){
399 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
400 } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
404 if (!iDetType->GetSegmentationModel()) {
405 AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(fITSgeom);
406 SetSegmentationModel(2,seg2);
408 if (!iDetType->GetResponseModel()) {
409 SetResponseModel(2,new AliITSresponseSSD("simulated"));
411 const char *kData2=(iDetType->GetResponseModel())->DataType();
412 if (strstr(kData2,"real")) {
413 iDetType->ClassNames("AliITSdigit","AliITSRawClusterSSD");
414 } else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
417 Warning("SetDefaults",
418 "Only the three basic detector types are initialized!");
421 //______________________________________________________________________
422 void AliITS::SetDefaultSimulation(){
423 // sets the default simulation.
431 AliITSDetType *iDetType;
432 AliITSsimulation *sim;
434 sim = iDetType->GetSimulationModel();
436 AliITSsegmentation *seg0=
437 (AliITSsegmentation*)iDetType->GetSegmentationModel();
438 AliITSresponse *res0 = (AliITSresponse*)iDetType->GetResponseModel();
439 AliITSsimulationSPD *sim0=new AliITSsimulationSPD(seg0,res0);
440 SetSimulationModel(0,sim0);
441 }else{ // simulation exists, make sure it is set up properly.
442 ((AliITSsimulationSPD*)sim)->Init(
443 (AliITSsegmentationSPD*) iDetType->GetSegmentationModel(),
444 (AliITSresponseSPD*) iDetType->GetResponseModel());
445 // if(sim->GetResponseModel()==0) sim->SetResponseModel(
446 // (AliITSresponse*)iDetType->GetResponseModel());
447 // if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
448 // (AliITSsegmentation*)iDetType->GetSegmentationModel());
451 sim = iDetType->GetSimulationModel();
453 AliITSsegmentation *seg1=
454 (AliITSsegmentation*)iDetType->GetSegmentationModel();
455 AliITSresponse *res1 = (AliITSresponse*)iDetType->GetResponseModel();
456 AliITSsimulationSDD *sim1=new AliITSsimulationSDD(seg1,res1);
457 SetSimulationModel(1,sim1);
458 }else{ // simulation exists, make sure it is set up properly.
459 ((AliITSsimulationSDD*)sim)->Init(
460 (AliITSsegmentationSDD*) iDetType->GetSegmentationModel(),
461 (AliITSresponseSDD*) iDetType->GetResponseModel());
462 // if(sim->GetResponseModel()==0) sim->SetResponseModel(
463 // (AliITSresponse*)iDetType->GetResponseModel());
464 // if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
465 // (AliITSsegmentation*)iDetType->GetSegmentationModel());
468 sim = iDetType->GetSimulationModel();
470 AliITSsegmentation *seg2=
471 (AliITSsegmentation*)iDetType->GetSegmentationModel();
472 AliITSresponse *res2 = (AliITSresponse*)iDetType->GetResponseModel();
473 AliITSsimulationSSD *sim2=new AliITSsimulationSSD(seg2,res2);
474 SetSimulationModel(2,sim2);
475 }else{ // simulation exists, make sure it is set up properly.
476 ((AliITSsimulationSSD*)sim)->Init(
477 (AliITSsegmentationSSD*) iDetType->GetSegmentationModel(),
478 (AliITSresponseSSD*) iDetType->GetResponseModel());
479 // if(sim->GetResponseModel()==0) sim->SetResponseModel(
480 // (AliITSresponse*)iDetType->GetResponseModel());
481 // if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
482 // (AliITSsegmentation*)iDetType->GetSegmentationModel());
485 //______________________________________________________________________
486 void AliITS::SetDefaultClusterFinders(){
487 // Sets the default cluster finders. Used in finding RecPoints.
496 AliITSDetType *iDetType;
500 if (!iDetType->GetReconstructionModel()) {
501 AliITSsegmentation *seg0 =
502 (AliITSsegmentation*)iDetType->GetSegmentationModel();
503 TClonesArray *dig0=DigitsAddress(0);
504 TClonesArray *recp0=ClustersAddress(0);
505 AliITSClusterFinderSPD *rec0 = new AliITSClusterFinderSPD(seg0,dig0,
507 SetReconstructionModel(0,rec0);
512 if (!iDetType->GetReconstructionModel()) {
513 AliITSsegmentation *seg1 =
514 (AliITSsegmentation*)iDetType->GetSegmentationModel();
515 AliITSresponse *res1 = (AliITSresponse*)iDetType->GetResponseModel();
516 TClonesArray *dig1=DigitsAddress(1);
517 TClonesArray *recp1=ClustersAddress(1);
518 AliITSClusterFinderSDD *rec1 =
519 new AliITSClusterFinderSDD(seg1,res1,dig1,recp1);
520 SetReconstructionModel(1,rec1);
525 if (!iDetType->GetReconstructionModel()) {
526 AliITSsegmentation *seg2=
527 (AliITSsegmentation*)iDetType->GetSegmentationModel();
528 TClonesArray *dig2=DigitsAddress(2);
529 AliITSClusterFinderSSD *rec2= new AliITSClusterFinderSSD(seg2,dig2);
530 SetReconstructionModel(2,rec2);
533 //______________________________________________________________________
534 void AliITS::MakeBranch(Option_t* option){
535 // Creates Tree branches for the ITS.
537 // Option_t *option String of Tree types S,D, and/or R.
538 // const char *file String of the file name where these branches
539 // are to be stored. If blank then these branches
540 // are written to the same tree as the Hits were
546 Bool_t cH = (strstr(option,"H")!=0);
547 Bool_t cS = (strstr(option,"S")!=0);
548 Bool_t cD = (strstr(option,"D")!=0);
549 Bool_t cR = (strstr(option,"R")!=0);
550 Bool_t cRF = (strstr(option,"RF")!=0);
553 if(cH && (fHits == 0x0)) fHits = new TClonesArray("AliITShit", 1560);
555 AliDetector::MakeBranch(option);
557 if(cS) MakeBranchS(0);
558 if(cD) MakeBranchD(0);
559 if(cR) MakeBranchR(0);
560 if(cRF) MakeBranchRF(0);
562 //______________________________________________________________________
563 void AliITS::SetTreeAddress(){
564 // Set branch address for the Trees.
571 TTree *treeS = fLoader->TreeS();
572 TTree *treeD = fLoader->TreeD();
573 TTree *treeR = fLoader->TreeR();
574 if (fLoader->TreeH() && (fHits == 0x0)) fHits = new TClonesArray("AliITShit", 1560);
576 AliDetector::SetTreeAddress();
578 SetTreeAddressS(treeS);
579 SetTreeAddressD(treeD);
580 SetTreeAddressR(treeR);
582 //______________________________________________________________________
583 AliITSDetType* AliITS::DetType(Int_t id){
584 // Return pointer to id detector type.
586 // Int_t id detector id number.
590 // returned, a pointer to a AliITSDetType.
592 return ((AliITSDetType*) fDetTypes->At(id));
594 //______________________________________________________________________
595 void AliITS::SetResponseModel(Int_t id, AliITSresponse *response){
596 // Set the response model for the id detector type.
598 // Int_t id detector id number.
599 // AliITSresponse* a pointer containing an instance of AliITSresponse
600 // to be stored/owned b y AliITSDetType.
606 ((AliITSDetType*) fDetTypes->At(id))->ResponseModel(response);
608 //______________________________________________________________________
609 void AliITS::SetSegmentationModel(Int_t id, AliITSsegmentation *seg){
610 // Set the segmentation model for the id detector type.
612 // Int_t id detector id number.
613 // AliITSsegmentation* a pointer containing an instance of
614 // AliITSsegmentation to be stored/owned b y
621 ((AliITSDetType*) fDetTypes->At(id))->SegmentationModel(seg);
623 //______________________________________________________________________
624 void AliITS::SetSimulationModel(Int_t id, AliITSsimulation *sim){
625 // Set the simulation model for the id detector type.
627 // Int_t id detector id number.
628 // AliITSresponse* a pointer containing an instance of AliITSresponse
629 // to be stored/owned b y AliITSDetType.
635 ((AliITSDetType*) fDetTypes->At(id))->SimulationModel(sim);
638 //______________________________________________________________________
639 void AliITS::SetReconstructionModel(Int_t id, AliITSClusterFinder *reconst){
640 // Set the cluster finder model for the id detector type.
642 // Int_t id detector id number.
643 // AliITSClusterFinder* a pointer containing an instance of
644 // AliITSClusterFinder to be stored/owned b y
651 ((AliITSDetType*) fDetTypes->At(id))->ReconstructionModel(reconst);
653 //______________________________________________________________________
654 void AliITS::SetClasses(Int_t id, const char *digit, const char *cluster){
655 // Set the digit and cluster classes name to be used for the id detector
658 // Int_t id detector id number.
659 // const char *digit Digit class name for detector id.
660 // const char *cluster Cluster class name for detector id.
666 ((AliITSDetType*) fDetTypes->At(id))->ClassNames(digit,cluster);
668 //______________________________________________________________________
669 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
671 // The function to add information to the AliITShit class. See the
672 // AliITShit class for a full description. This function allocates the
673 // necessary new space for the hit information and passes the variable
674 // track, and the pointers *vol and *hits to the AliITShit constructor
677 // Int_t track Track number which produced this hit.
678 // Int_t *vol Array of Integer Hit information. See AliITShit.h
679 // Float_t *hits Array of Floating Hit information. see AliITShit.h
685 TClonesArray &lhits = *fHits;
686 new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
688 //______________________________________________________________________
689 void AliITS::InitModules(Int_t size,Int_t &nmodules){
690 // Initialize the modules array.
692 // Int_t size Size of array of the number of modules to be
693 // created. If size <=0 then the number of modules
694 // is gotten from AliITSgeom class kept in fITSgeom.
696 // Int_t &nmodules The number of modules existing.
701 fITSmodules->Delete();
703 } // end fir fITSmoudles
705 Int_t nl,indexMAX,index;
707 if(size<=0){ // default to using data stored in AliITSgeom
709 Error("InitModules","fITSgeom not defined");
711 } // end if fITSgeom==0
712 nl = fITSgeom->GetNlayers();
713 indexMAX = fITSgeom->GetModuleIndex(nl,fITSgeom->GetNladders(nl),
714 fITSgeom->GetNdetectors(nl))+1;
716 fITSmodules = new TObjArray(indexMAX);
717 for(index=0;index<indexMAX;index++){
718 fITSmodules->AddAt( new AliITSmodule(index),index);
721 fITSmodules = new TObjArray(size);
722 for(index=0;index<size;index++) {
723 fITSmodules->AddAt( new AliITSmodule(index),index);
729 //______________________________________________________________________
730 void AliITS::FillModules(Int_t evnt,Int_t bgrev,Int_t nmodules,
731 Option_t *option,Text_t *filename){
732 // fill the modules with the sorted by module hits; add hits from
733 // background if option=Add.
735 // Int_t evnt Event to be processed.
736 // Int_t bgrev Background Hit tree number.
737 // Int_t nmodules Not used.
738 // Option_t *option String indicating if merging hits or not. To
739 // merge hits set equal to "Add". Otherwise no
740 // background hits are considered.
741 // Test_t *filename File name containing the background hits..
746 static TTree *trH1; //Tree with background hits
747 static Bool_t first=kTRUE;
749 const char *addBgr = strstr(option,"Add");
751 evnt = nmodules; // Dummy use of variables to remove warnings
754 file=new TFile(filename);
759 // Get Hits Tree header from file
760 if(trH1) delete trH1;
764 sprintf(treeName,"TreeH%d",bgrev);
765 trH1 = (TTree*)gDirectory->Get(treeName);
767 Error("FillModules","cannot find Hits Tree for event:%d",bgrev);
769 // Set branch addresses
772 FillModules(fLoader->TreeH(),0); // fill from this file's tree.
775 FillModules(trH1,10000000); // Default mask 10M.
776 TTree *fAli=fLoader->GetRunLoader()->TreeK();
778 if (fAli) fileAli =fAli->GetCurrentFile();
782 //______________________________________________________________________
783 void AliITS::FillModules(TTree *treeH, Int_t mask) {
784 // fill the modules with the sorted by module hits;
785 // can be called many times to do a merging
787 // TTree *treeH The tree containing the hits to be copied into
789 // Int_t mask The track number mask to indecate which file
790 // this hits came from.
798 Error("FillModules","Tree is NULL");
800 Int_t lay,lad,det,index;
804 sprintf(branchname,"%s",GetName());
805 TBranch *branch = treeH->GetBranch(branchname);
807 Error("FillModules","%s branch in TreeH not found",branchname);
810 branch->SetAddress(&fHits);
811 Int_t nTracks =(Int_t) treeH->GetEntries();
813 for(iPrimTrack=0; iPrimTrack<nTracks; iPrimTrack++){
815 Int_t nBytes = treeH->GetEvent(iPrimTrack);
816 if (nBytes <= 0) continue;
817 Int_t nHits = fHits->GetEntriesFast();
818 for(h=0; h<nHits; h++){
819 itsHit = (AliITShit *)fHits->UncheckedAt(h);
820 itsHit->GetDetectorID(lay,lad,det);
822 index = fITSgeom->GetModuleIndex(lay,lad,det);
824 index=det-1; // This should not be used.
825 } // end if [You must have fITSgeom for this to work!]
826 mod = GetModule(index);
827 itsHit->SetTrack(itsHit->GetTrack()+mask); // Set track mask.
828 mod->AddHit(itsHit,iPrimTrack,h);
829 } // end loop over hits
830 } // end loop over tracks
832 //______________________________________________________________________
833 void AliITS::ClearModules(){
834 // Clear the modules TObjArray.
840 if(fITSmodules) fITSmodules->Delete();
842 //______________________________________________________________________
843 void AliITS::MakeBranchS(const char *fl){
844 // Creates Tree Branch for the ITS summable digits.
846 // cont char *fl File name where SDigits branch is to be written
847 // to. If blank it write the SDigits to the same
848 // file in which the Hits were found.
853 Int_t buffersize = 4000;
856 // only one branch for SDigits.
857 sprintf(branchname,"%s",GetName());
860 if(fLoader->TreeS()){
861 if (fSDigits == 0x0) fSDigits = new TClonesArray("AliITSpListItem",1000);
862 MakeBranchInTree(fLoader->TreeS(),branchname,&fSDigits,buffersize,fl);
865 //______________________________________________________________________
866 void AliITS::SetTreeAddressS(TTree *treeS){
867 // Set branch address for the ITS summable digits Trees.
869 // TTree *treeS Tree containing the SDigits.
877 if (fSDigits == 0x0) fSDigits = new TClonesArray("AliITSpListItem",1000);
879 sprintf(branchname,"%s",GetName());
880 branch = treeS->GetBranch(branchname);
881 if (branch) branch->SetAddress(&fSDigits);
883 //______________________________________________________________________
884 void AliITS::MakeBranchInTreeD(TTree *treeD,const char *file){
885 // Creates Tree branches for the ITS.
887 // TTree *treeD Pointer to the Digits Tree.
888 // cont char *file File name where Digits branch is to be written
889 // to. If blank it write the SDigits to the same
890 // file in which the Hits were found.
895 Int_t buffersize = 4000;
898 sprintf(branchname,"%s",GetName());
899 // one branch for digits per type of detector
900 const char *det[3] = {"SPD","SDD","SSD"};
904 for (i=0; i<kNTYPES ;i++) {
905 DetType(i)->GetClassNames(digclass,clclass);
907 if (fDtype == 0x0) fDtype = new TObjArray(fNDetTypes);
908 if(!(fDtype->At(i))) fDtype->AddAt(new TClonesArray(digclass,1000),i);
911 for (i=0; i<kNTYPES ;i++) {
912 if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
913 else sprintf(branchname,"%sDigits%d",GetName(),i+1);
914 if (fDtype && treeD) {
915 MakeBranchInTree(treeD, branchname, &((*fDtype)[i]),buffersize,file);
919 //______________________________________________________________________
920 void AliITS::SetTreeAddressD(TTree *treeD){
921 // Set branch address for the Trees.
923 // TTree *treeD Tree containing the Digits.
929 const char *det[3] = {"SPD","SDD","SSD"};
936 for (i=0; i<kNTYPES; i++) {
937 DetType(i)->GetClassNames(digclass,clclass);
939 if (fDtype == 0x0) fDtype = new TObjArray(fNDetTypes);
940 if(!(fDtype->At(i))) fDtype->AddAt(new TClonesArray(digclass,1000),i);
942 if (kNTYPES==3) sprintf(branchname,"%sDigits%s",GetName(),det[i]);
943 else sprintf(branchname,"%sDigits%d",GetName(),i+1);
945 branch = treeD->GetBranch(branchname);
946 if (branch) branch->SetAddress(&((*fDtype)[i]));
950 //______________________________________________________________________
951 void AliITS::Hits2SDigits(){
952 // Standard Hits to summable Digits function.
958 // return; // Using Hits in place of the larger sDigits.
959 fLoader->LoadHits("read");
960 fLoader->LoadSDigits("recreate");
961 AliRunLoader* rl = fLoader->GetRunLoader();
963 for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
964 // Do the Hits to Digits operation. Use Standard input values.
965 // Event number from file, no background hit merging , use size from
966 // AliITSgeom class, option="All", input from this file only.
967 rl->GetEvent(iEvent);
968 if (!fLoader->TreeS()) fLoader->MakeTree("S");
971 HitsToSDigits(iEvent,0,-1," ",fOpt," ");
974 fLoader->UnloadHits();
975 fLoader->UnloadSDigits();
977 //______________________________________________________________________
978 void AliITS::Hits2PreDigits(){
979 // Standard Hits to summable Digits function.
985 AliHeader *header=fLoader->GetRunLoader()->GetHeader(); // Get event number from this file.
986 // Do the Hits to Digits operation. Use Standard input values.
987 // Event number from file, no background hit merging , use size from
988 // AliITSgeom class, option="All", input from this file only.
989 HitsToPreDigits(header->GetEvent(),0,-1," ",fOpt," ");
991 //______________________________________________________________________
992 AliDigitizer* AliITS::CreateDigitizer(AliRunDigitizer* manager) const
994 return new AliITSDigitizer(manager);
996 //______________________________________________________________________
997 void AliITS::SDigitsToDigits(Option_t *opt){
998 // Standard Summable digits to Digits function.
1003 char name[20] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
1005 if(!GetITSgeom()) return; // need transformations to do digitization.
1006 AliITSgeom *geom = GetITSgeom();
1008 const char *all = strstr(opt,"All");
1009 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1011 if( !det[0] && !det[1] && !det[2] ) all = "All";
1013 static Bool_t setDef=kTRUE;
1014 if (setDef) SetDefaultSimulation();
1017 AliITSsimulation *sim = 0;
1018 AliITSDetType *iDetType = 0;
1019 TTree *trees = fLoader->TreeS();
1020 if( !(trees && this->GetSDigits()) ){
1021 Error("SDigits2Digits","Error: No trees or SDigits. Returning.");
1024 sprintf( name, "%s", this->GetName() );
1025 TBranch *brchSDigits = trees->GetBranch( name );
1028 for(module=0;module<geom->GetIndexMax();module++){
1029 id = geom->GetModuleType(module);
1030 if (!all && !det[id]) continue;
1031 iDetType = DetType(id);
1032 sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1034 Error("SDigit2Digits",
1035 "The simulation class was not instanciated!");
1038 sim->InitSimulationModule(module,gAlice->GetEvNumber());
1040 // add summable digits to module
1041 this->GetSDigits()->Clear();
1042 brchSDigits->GetEvent(module);
1043 sim->AddSDigitsToModule(GetSDigits(),0);
1045 // Digitise current module sum(SDigits)->Digits
1046 sim->FinishSDigitiseModule();
1048 // fills all branches - wasted disk space
1049 fLoader->TreeD()->Fill();
1050 this->ResetDigits();
1053 fLoader->TreeD()->GetEntries();
1055 fLoader->TreeD()->AutoSave();
1057 fLoader->TreeD()->Reset();
1060 //______________________________________________________________________
1061 void AliITS::Hits2Digits(){
1062 // Standard Hits to Digits function.
1068 fLoader->LoadHits("read");
1069 fLoader->LoadDigits("recreate");
1070 AliRunLoader* rl = fLoader->GetRunLoader();
1072 for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
1073 // Do the Hits to Digits operation. Use Standard input values.
1074 // Event number from file, no background hit merging , use size from
1075 // AliITSgeom class, option="All", input from this file only.
1076 rl->GetEvent(iEvent);
1077 if (!fLoader->TreeD()) fLoader->MakeTree("D");
1080 HitsToDigits(iEvent,0,-1," ",fOpt," ");
1083 fLoader->UnloadHits();
1084 fLoader->UnloadDigits();
1086 //______________________________________________________________________
1087 void AliITS::HitsToSDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1088 Option_t *option, Option_t *opt,Text_t *filename){
1089 // keep galice.root for signal and name differently the file for
1090 // background when add! otherwise the track info for signal will be lost !
1091 // the condition below will disappear when the geom class will be
1092 // initialized for all versions - for the moment it is only for v5 !
1093 // 7 is the SDD beam test version. Dummy routine. Hits are ITS's Summable
1096 // Int_t evnt Event to be processed.
1097 // Int_t bgrev Background Hit tree number.
1098 // Int_t nmodules Not used.
1099 // Option_t *option String indicating if merging hits or not. To
1100 // merge hits set equal to "Add". Otherwise no
1101 // background hits are considered.
1102 // Test_t *filename File name containing the background hits..
1107 // return; // using Hits instead of the larger sdigits.
1109 HitsToPreDigits(evNumber,bgrev,size,option,opt,filename);
1111 //______________________________________________________________________
1112 void AliITS::HitsToPreDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1113 Option_t *option, Option_t *opt,Text_t *filename){
1114 // Keep galice.root for signal and name differently the file for
1115 // background when add! otherwise the track info for signal will be lost !
1116 // the condition below will disappear when the geom class will be
1117 // initialized for all versions - for the moment it is only for v5 !
1118 // 7 is the SDD beam test version.
1120 // Int_t evnt Event to be processed.
1121 // Int_t bgrev Background Hit tree number.
1122 // Int_t nmodules Not used.
1123 // Option_t *option String indicating if merging hits or not. To
1124 // merge hits set equal to "Add". Otherwise no
1125 // background hits are considered.
1126 // Test_t *filename File name containing the background hits..
1132 if(!GetITSgeom()) return; // need transformations to do digitization.
1133 AliITSgeom *geom = GetITSgeom();
1135 const char *all = strstr(opt,"All");
1136 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1138 static Bool_t setDef=kTRUE;
1139 if (setDef) SetDefaultSimulation();
1143 InitModules(size,nmodules);
1144 FillModules(evNumber,bgrev,nmodules,option,filename);
1146 AliITSsimulation *sim = 0;
1147 AliITSDetType *iDetType = 0;
1148 AliITSmodule *mod = 0;
1150 for(module=0;module<geom->GetIndexMax();module++){
1151 id = geom->GetModuleType(module);
1152 if (!all && !det[id]) continue;
1153 iDetType = DetType(id);
1154 sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1156 Error("HitsToSDigits",
1157 "The simulation class was not instanciated!");
1160 mod = (AliITSmodule *)fITSmodules->At(module);
1161 sim->SDigitiseModule(mod,module,evNumber);
1162 // fills all branches - wasted disk space
1163 fLoader->TreeS()->Fill();
1169 fLoader->TreeS()->GetEntries();
1170 fLoader->TreeS()->AutoSave();
1171 fLoader->WriteSDigits("OVERWRITE");
1173 fLoader->TreeS()->Reset();
1175 //______________________________________________________________________
1176 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size,
1177 Option_t *option, Option_t *opt,Text_t *filename){
1178 // Keep galice.root for signal and name differently the file for
1179 // background when add! otherwise the track info for signal will be lost !
1180 // the condition below will disappear when the geom class will be
1181 // initialized for all versions - for the moment it is only for v5 !
1182 // 7 is the SDD beam test version.
1184 // Int_t evnt Event to be processed.
1185 // Int_t bgrev Background Hit tree number.
1186 // Int_t nmodules Not used.
1187 // Option_t *option String indicating if merging hits or not. To
1188 // merge hits set equal to "Add". Otherwise no
1189 // background hits are considered.
1190 // Test_t *filename File name containing the background hits..
1196 if(!GetITSgeom()) return; // need transformations to do digitization.
1197 AliITSgeom *geom = GetITSgeom();
1199 const char *all = strstr(opt,"All");
1200 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1202 static Bool_t setDef=kTRUE;
1203 if (setDef) SetDefaultSimulation();
1207 InitModules(size,nmodules);
1208 FillModules(evNumber,bgrev,nmodules,option,filename);
1210 AliITSsimulation *sim = 0;
1211 AliITSDetType *iDetType = 0;
1212 AliITSmodule *mod = 0;
1214 for(module=0;module<geom->GetIndexMax();module++){
1215 id = geom->GetModuleType(module);
1216 if (!all && !det[id]) continue;
1217 iDetType = DetType(id);
1218 sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1220 Error("HitsToDigits",
1221 "The simulation class was not instanciated!");
1224 mod = (AliITSmodule *)fITSmodules->At(module);
1225 sim->DigitiseModule(mod,module,evNumber);
1226 // fills all branches - wasted disk space
1227 fLoader->TreeD()->Fill();
1233 fLoader->TreeD()->GetEntries();
1234 fLoader->TreeD()->AutoSave();
1236 fLoader->TreeD()->Reset();
1238 //______________________________________________________________________
1239 void AliITS::ResetSDigits(){
1240 // Reset the Summable Digits array.
1246 if (fSDigits) fSDigits->Clear();
1249 //______________________________________________________________________
1250 void AliITS::ResetDigits(){
1251 // Reset number of digits and the digits array for the ITS detector.
1257 if (!fDtype) return;
1260 for (i=0;i<kNTYPES;i++ ) {
1261 if (fDtype->At(i)) ((TClonesArray*)fDtype->At(i))->Clear();
1262 if (fNdtype) fNdtype[i]=0;
1265 //______________________________________________________________________
1266 void AliITS::ResetDigits(Int_t i){
1267 // Reset number of digits and the digits array for this branch.
1273 if (fDtype->At(i)) ((TClonesArray*)fDtype->At(i))->Clear();
1274 if (fNdtype) fNdtype[i]=0;
1276 //______________________________________________________________________
1277 void AliITS::AddSumDigit(AliITSpListItem &sdig){
1278 // Adds the a module full of summable digits to the summable digits tree.
1280 // AliITSpListItem &sdig SDigit to be added to SDigits tree.
1286 TClonesArray &lsdig = *fSDigits;
1287 new(lsdig[fNSDigits++]) AliITSpListItem(sdig);
1289 //______________________________________________________________________
1290 void AliITS::AddRealDigit(Int_t id, Int_t *digits){
1291 // Add a real digit - as coming from data.
1293 // Int_t id Detector type number.
1294 // Int_t *digits Integer array containing the digits info. See
1301 TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1302 new(ldigits[fNdtype[id]++]) AliITSdigit(digits);
1304 //______________________________________________________________________
1305 void AliITS::AddSimDigit(Int_t id, AliITSdigit *d){
1306 // Add a simulated digit.
1308 // Int_t id Detector type number.
1309 // AliITSdigit *d Digit to be added to the Digits Tree. See
1316 TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1320 new(ldigits[fNdtype[id]++]) AliITSdigitSPD(*((AliITSdigitSPD*)d));
1323 new(ldigits[fNdtype[id]++]) AliITSdigitSDD(*((AliITSdigitSDD*)d));
1326 new(ldigits[fNdtype[id]++]) AliITSdigitSSD(*((AliITSdigitSSD*)d));
1330 //______________________________________________________________________
1331 void AliITS::AddSimDigit(Int_t id,Float_t phys,Int_t *digits,Int_t *tracks,
1332 Int_t *hits,Float_t *charges){
1333 // Add a simulated digit to the list.
1335 // Int_t id Detector type number.
1336 // Float_t phys Physics indicator. See AliITSdigits.h
1337 // Int_t *digits Integer array containing the digits info. See
1339 // Int_t *tracks Integer array [AliITSdigitS?D::GetNTracks()]
1340 // containing the track numbers that contributed to
1342 // Int_t *hits Integer array [AliITSdigitS?D::GetNTracks()]
1343 // containing the hit numbers, from AliITSmodule, that
1344 // contributed to this digit.
1345 // Float_t *charge Floating point array of the signals contributed
1346 // to this digit by each track.
1352 TClonesArray &ldigits = *((TClonesArray*)fDtype->At(id));
1355 new(ldigits[fNdtype[id]++]) AliITSdigitSPD(digits,tracks,hits);
1358 new(ldigits[fNdtype[id]++]) AliITSdigitSDD(phys,digits,tracks,
1362 new(ldigits[fNdtype[id]++]) AliITSdigitSSD(digits,tracks,hits);
1366 //______________________________________________________________________
1367 void AliITS::MakeTreeC(Option_t *option){
1368 // Create a separate tree to store the clusters.
1370 // Option_t *option string which must contain "C" otherwise
1371 // no Cluster Tree is created.
1377 AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;
1379 if (pITSLoader == 0x0) {
1380 Error("MakeTreeC","fLoader == 0x0 option=%s",option);
1383 if (pITSLoader->TreeC() == 0x0) pITSLoader->MakeTree("C");
1387 void AliITS::MakeBranchC()
1389 //Makes barnches in treeC
1390 AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;
1391 if (pITSLoader == 0x0)
1393 Error("MakeTreeC","fLoader == 0x0");
1396 TTree * lTC = pITSLoader->TreeC();
1399 Error("MakeTreeC","Can not get TreeC from Loader");
1403 Int_t buffersize = 4000;
1404 char branchname[30];
1405 const char *det[3] = {"SPD","SDD","SSD"};
1409 // one branch for Clusters per type of detector
1411 for (i=0; i<kNTYPES ;i++)
1413 AliITSDetType *iDetType=DetType(i);
1414 iDetType->GetClassNames(digclass,clclass);
1416 if (fCtype == 0x0) fCtype = new TObjArray(fNDetTypes);
1417 if(!ClustersAddress(i))
1419 fCtype->AddAt(new TClonesArray(clclass,1000),i);
1421 if (kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
1422 else sprintf(branchname,"%sClusters%d",GetName(),i+1);
1425 if (lTC->GetBranch(branchname))
1427 Warning("MakeBranchC","Branch %s alread exists in TreeC",branchname);
1431 Info("MakeBranchC","Creating branch %s in TreeC",branchname);
1432 lTC->Branch(branchname,&((*fCtype)[i]), buffersize);
1434 } // end if fCtype && lTC
1438 //______________________________________________________________________
1439 void AliITS::GetTreeC(Int_t event){
1440 // Get the clusters tree for this event and set the branch address.
1442 // Int_t event Event number for the cluster tree.
1447 char branchname[30];
1448 const char *det[3] = {"SPD","SDD","SSD"};
1450 AliITSLoader *pITSLoader = (AliITSLoader*)fLoader;
1451 TTree * lTC = pITSLoader->TreeC();
1455 pITSLoader->CleanRawClusters();
1465 for (i=0; i<kNTYPES; i++) {
1466 AliITSDetType *iDetType=DetType(i);
1467 iDetType->GetClassNames(digclass,clclass);
1469 if (fCtype == 0x0) fCtype = new TObjArray(fNDetTypes);
1470 if(!fCtype->At(i)) fCtype->AddAt(new TClonesArray(clclass,1000),i);
1471 if(kNTYPES==3) sprintf(branchname,"%sClusters%s",GetName(),det[i]);
1472 else sprintf(branchname,"%sClusters%d",GetName(),i+1);
1474 branch = lTC->GetBranch(branchname);
1475 if (branch) branch->SetAddress(&((*fCtype)[i]));
1479 Error("GetTreeC","cannot find Clusters Tree for event:%d",event);
1482 //______________________________________________________________________
1483 void AliITS::AddCluster(Int_t id, AliITSRawCluster *c){
1484 // Add a cluster to the list.
1486 // Int_t id Detector type number.
1487 // AliITSRawCluster *c Cluster class to be added to the tree of
1494 TClonesArray &lc = *((TClonesArray*)fCtype->At(id));
1498 new(lc[fNctype[id]++]) AliITSRawClusterSPD(*((AliITSRawClusterSPD*)c));
1501 new(lc[fNctype[id]++]) AliITSRawClusterSDD(*((AliITSRawClusterSDD*)c));
1504 new(lc[fNctype[id]++]) AliITSRawClusterSSD(*((AliITSRawClusterSSD*)c));
1508 //______________________________________________________________________
1509 void AliITS::ResetClusters(){
1510 // Reset number of clusters and the clusters array for ITS.
1517 for (i=0;i<kNTYPES;i++ ) ResetClusters(i);
1519 //______________________________________________________________________
1520 void AliITS::ResetClusters(Int_t i){
1521 // Reset number of clusters and the clusters array for this branch.
1523 // Int_t i Detector type number.
1529 if (fCtype->At(i)) ((TClonesArray*)fCtype->At(i))->Clear();
1530 if (fNctype) fNctype[i]=0;
1532 //______________________________________________________________________
1533 void AliITS::MakeBranchR(const char *file, Option_t *opt){
1534 // Creates Tree branches for the ITS Reconstructed points.
1536 // cont char *file File name where RecPoints branch is to be written
1537 // to. If blank it write the SDigits to the same
1538 // file in which the Hits were found.
1543 Int_t buffsz = 4000;
1544 char branchname[30];
1546 // only one branch for rec points for all detector types
1547 Bool_t oFast= (strstr(opt,"Fast")!=0);
1549 sprintf(branchname,"%sRecPointsF",GetName());
1551 sprintf(branchname,"%sRecPoints",GetName());
1555 if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1556 if (fLoader->TreeR()) {
1557 if (fRecPoints == 0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1558 MakeBranchInTree(fLoader->TreeR(),branchname,&fRecPoints,buffsz,file);
1561 //______________________________________________________________________
1562 void AliITS::SetTreeAddressR(TTree *treeR){
1563 // Set branch address for the Reconstructed points Trees.
1565 // TTree *treeR Tree containing the RecPoints.
1570 char branchname[30];
1573 if (fRecPoints == 0x0) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
1575 sprintf(branchname,"%sRecPoints",GetName());
1576 branch = treeR->GetBranch(branchname);
1578 branch->SetAddress(&fRecPoints);
1581 sprintf(branchname,"%sRecPointsF",GetName());
1582 branch = treeR->GetBranch(branchname);
1584 branch->SetAddress(&fRecPoints);
1588 //______________________________________________________________________
1589 void AliITS::AddRecPoint(const AliITSRecPoint &r){
1590 // Add a reconstructed space point to the list
1592 // const AliITSRecPoint &r RecPoint class to be added to the tree
1593 // of reconstructed points TreeR.
1599 TClonesArray &lrecp = *fRecPoints;
1600 new(lrecp[fNRecPoints++]) AliITSRecPoint(r);
1602 //______________________________________________________________________
1603 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
1604 Option_t *opt0,Option_t *opt1,Text_t *flnm){
1605 // keep galice.root for signal and name differently the file for
1606 // background when add! otherwise the track info for signal will be lost !
1607 // the condition below will disappear when the geom class will be
1608 // initialized for all versions - for the moment it is only for v5 !
1610 // Int_t evnt Event to be processed.
1611 // Int_t bgrev Background Hit tree number.
1612 // Int_t size Size used by InitModules. See InitModules.
1613 // Option_t *opt0 Option passed to FillModules. See FillModules.
1614 // Option_t *opt1 String indicating if merging hits or not. To
1615 // merge hits set equal to "Add". Otherwise no
1616 // background hits are considered.
1617 // Test_t *flnm File name containing the background hits..
1623 if(!GetITSgeom()) return;
1624 AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
1625 AliITSgeom *geom = GetITSgeom();
1627 const char *all = strstr(opt1,"All");
1628 const char *det[3] ={strstr(opt1,"SPD"),strstr(opt1,"SDD"),
1629 strstr(opt1,"SSD")};
1631 InitModules(size,nmodules);
1632 FillModules(evNumber,bgrev,nmodules,opt0,flnm);
1634 AliITSsimulation *sim = 0;
1635 AliITSDetType *iDetType = 0;
1636 AliITSmodule *mod = 0;
1639 //m.b. : this change is nothing but a nice way to make sure
1642 cout<<"HitsToFastRecPoints: N mod = "<<geom->GetIndexMax()<<endl;
1644 for(module=0;module<geom->GetIndexMax();module++)
1646 id = geom->GetModuleType(module);
1647 if (!all && !det[id]) continue;
1648 iDetType = DetType(id);
1649 sim = (AliITSsimulation*)iDetType->GetSimulationModel();
1652 Error("HitsToFastPoints","The simulation class was not instanciated!");
1655 mod = (AliITSmodule *)fITSmodules->At(module);
1656 sim->CreateFastRecPoints(mod,module,gRandom);
1657 cout<<module<<"\r";fflush(0);
1658 //gAlice->TreeR()->Fill();
1659 TTree *lTR = pITSloader->TreeR();
1660 TBranch *br=lTR->GetBranch("ITSRecPointsF");
1667 fLoader->WriteRecPoints("OVERWRITE");
1669 //______________________________________________________________________
1670 void AliITS::Digits2Reco(){
1671 // Find clusters and reconstruct space points.
1677 AliHeader *header=fLoader->GetRunLoader()->GetHeader();
1678 // to Digits to RecPoints for event in file, all digits in file, and
1679 // all ITS detectors.
1680 DigitsToRecPoints(header->GetEvent(),0,fOpt);
1682 //______________________________________________________________________
1683 void AliITS::DigitsToRecPoints(Int_t evNumber,Int_t lastentry,Option_t *opt){
1684 // cluster finding and reconstruction of space points
1685 // the condition below will disappear when the geom class will be
1686 // initialized for all versions - for the moment it is only for v5 !
1687 // 7 is the SDD beam test version
1689 // Int_t evNumber Event number to be processed.
1690 // Int_t lastentry Offset for module when not all of the modules
1692 // Option_t *opt String indicating which ITS sub-detectors should
1693 // be processed. If ="All" then all of the ITS
1694 // sub detectors are processed.
1700 if(!GetITSgeom()) return;
1701 AliITSgeom *geom = GetITSgeom();
1703 const char *all = strstr(opt,"All");
1704 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
1706 static Bool_t setRec=kTRUE;
1707 if (setRec) SetDefaultClusterFinders();
1710 AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
1711 TTree *treeC=pITSloader->TreeC();
1712 AliITSClusterFinder *rec = 0;
1713 AliITSDetType *iDetType = 0;
1714 Int_t id,module,first=0;
1715 for(module=0;module<geom->GetIndexMax();module++){
1716 id = geom->GetModuleType(module);
1717 if (!all && !det[id]) continue;
1718 if(det[id]) first = geom->GetStartDet(id);
1719 iDetType = DetType(id);
1720 rec = (AliITSClusterFinder*)iDetType->GetReconstructionModel();
1721 TClonesArray *itsDigits = this->DigitsAddress(id);
1723 Error("DigitsToRecPoints",
1724 "The reconstruction class was not instanciated! event=%d",
1728 this->ResetDigits();
1729 TTree *lTD = pITSloader->TreeD();
1731 lTD->GetEvent(lastentry+module);
1733 lTD->GetEvent(lastentry+(module-first));
1735 Int_t ndigits = itsDigits->GetEntriesFast();
1736 if (ndigits) rec->FindRawClusters(module);
1737 pITSloader->TreeR()->Fill();
1744 pITSloader->WriteRecPoints("OVERWRITE");
1745 pITSloader->WriteRawClusters("OVERWRITE");
1747 //______________________________________________________________________
1748 void AliITS::ResetRecPoints(){
1749 // Reset number of rec points and the rec points array.
1755 if (fRecPoints) fRecPoints->Clear();
1758 //______________________________________________________________________
1759 AliLoader* AliITS::MakeLoader(const char* topfoldername)
1761 //builds ITSgetter (AliLoader type)
1762 //if detector wants to use castomized getter, it must overload this method
1764 Info("MakeLoader","Creating AliITSLoader. Top folder is %s.",topfoldername);
1765 fLoader = new AliITSLoader(GetName(),topfoldername);
1770 //_____________________________________________________________________________
1771 void AliITS::Reconstruct() const
1773 // reconstruct clusters
1775 AliLoader* loader = GetLoader();
1776 loader->LoadRecPoints("recreate");
1777 loader->LoadDigits("read");
1779 AliITSclustererV2 clusterer(GetITSgeom());
1780 AliRunLoader* runLoader = loader->GetRunLoader();
1781 Int_t nEvents = runLoader->GetNumberOfEvents();
1783 for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
1784 runLoader->GetEvent(iEvent);
1786 TTree* treeClusters = loader->TreeR();
1787 if (!treeClusters) {
1788 loader->MakeTree("R");
1789 treeClusters = loader->TreeR();
1791 TTree* treeDigits = loader->TreeD();
1793 Error("Reconstruct", "Can't get digits tree !");
1797 clusterer.Digits2Clusters(treeDigits, treeClusters);
1799 loader->WriteRecPoints("OVERWRITE");
1802 loader->UnloadRecPoints();
1803 loader->UnloadDigits();
1806 //_____________________________________________________________________________
1807 AliTracker* AliITS::CreateTracker() const
1809 // create an ITS tracker
1811 return new AliITStrackerV2(GetITSgeom());
1814 //_____________________________________________________________________________
1815 void AliITS::FillESD(AliESD* esd) const
1817 // make PID, find V0s and cascades
1819 Double_t parITS[] = {34., 0.15, 10.};
1820 AliITSpidESD itsPID(parITS);
1821 itsPID.MakePID(esd);
1824 Double_t cuts[]={33, // max. allowed chi2
1825 0.16,// min. allowed negative daughter's impact parameter
1826 0.05,// min. allowed positive daughter's impact parameter
1827 0.080,// max. allowed DCA between the daughter tracks
1828 0.998,// max. allowed cosine of V0's pointing angle
1829 0.9, // min. radius of the fiducial volume
1830 2.9 // max. radius of the fiducial volume
1832 AliV0vertexer vtxer(cuts);
1833 Double_t vtx[3], cvtx[6];
1834 esd->GetVertex(vtx, cvtx);
1835 vtxer.SetVertex(vtx);
1836 vtxer.Tracks2V0vertices(esd);
1839 Double_t cts[]={33., // max. allowed chi2
1840 0.05, // min. allowed V0 impact parameter
1841 0.008, // window around the Lambda mass
1842 0.035, // min. allowed bachelor's impact parameter
1843 0.10, // max. allowed DCA between a V0 and a track
1844 0.9985, //max. allowed cosine of the cascade pointing angle
1845 0.9, // min. radius of the fiducial volume
1846 2.9 // max. radius of the fiducial volume
1848 AliCascadeVertexer cvtxer=AliCascadeVertexer(cts);
1849 cvtxer.V0sTracks2CascadeVertices(esd);