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 ///////////////////////////////////////////////////////////////////////////////
69 #include <TClonesArray.h>
71 #include <TParticle.h>
74 #include <TVirtualMC.h>
75 #include "AliDetector.h"
77 #include "AliITSDetTypeSim.h"
78 #include "AliITSDDLRawData.h"
79 #include "AliITSLoader.h"
80 #include "AliITShit.h"
81 #include "AliITSmodule.h"
82 #include "AliITSpListItem.h"
83 #include "AliITSsimulation.h"
84 #include "AliITSsimulationFastPoints.h"
86 #include "AliITSRecPoint.h"
87 #include "AliITSsegmentationSPD.h"
88 #include "AliITSsegmentationSDD.h"
89 #include "AliITSsimulationSDD.h"
90 #include "AliITSCalibrationSDD.h"
91 #include "AliITSCalibrationSSD.h"
92 #include "AliITSsegmentationSSD.h"
93 #include "AliITSRawStreamSPD.h"
94 #include "AliITSRawStreamSSD.h"
95 #include "AliITSRawStreamSDD.h"
96 #include "AliRawReader.h"
99 #include "AliITSInitGeometry.h"
100 #include "AliITSFOSignalsSPD.h"
104 //______________________________________________________________________
105 AliITS::AliITS() : AliDetector(),
118 // Default initializer for ITS
119 // The default constructor of the AliITS class. In addition to
120 // creating the AliITS class it zeros the variables fIshunt (a member
121 // of AliDetector class), fEuclidOut, and fIdN, and zeros the pointers
122 // fIdSens, and fIdName. The AliDetector default constructor
125 // SetDetectors(); // default to fOpt="All". This variable not written out.
126 //PH SetMarkerColor(kRed);
128 //______________________________________________________________________
129 AliITS::AliITS(const Char_t *title):
130 AliDetector("ITS",title),
143 // The standard Constructor for the ITS class.
144 // It also zeros the variables
145 // fIshunt (a member of AliDetector class), fEuclidOut, and zeros
146 // the pointers fIdSens and fIdName. To help in displaying hits via the
147 // ROOT macro display.C AliITS also sets the marker color to red. The
148 // variables passes with this constructor, const char *name and *title,
149 // are used by the constructor of AliDetector class. See AliDetector
150 // class for a description of these parameters and its constructor
153 // Char_t *title Simulation title for the ITS
159 fHits = new TClonesArray("AliITShit",1560); // from AliDetector
160 if(gAlice->GetMCApp()) gAlice->GetMCApp()->AddHitList(fHits);
161 //fNhits=0; //done in AliDetector(name,title)
162 SetDetectors(); // default to fOpt="All". This variable not written out.
163 fDetTypeSim = new AliITSDetTypeSim();
164 //PH SetMarkerColor(kRed);
165 if(!fLoader) MakeLoader(AliConfig::GetDefaultEventFolderName());
166 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
168 //______________________________________________________________________
169 AliITS::AliITS(const char *name, const char *title):
170 AliDetector(name,title),
183 // The standard Constructor for the ITS class.
184 // It also zeros the variables
185 // fIshunt (a member of AliDetector class), fEuclidOut, and zeros
186 // the pointers fIdSens and fIdName. To help in displaying hits via the
187 // ROOT macro display.C AliITS also sets the marker color to red. The
188 // variables passes with this constructor, const char *name and *title,
189 // are used by the constructor of AliDetector class. See AliDetector
190 // class for a description of these parameters and its constructor
193 fHits = new TClonesArray("AliITShit",1560);
194 if(gAlice->GetMCApp()) gAlice->GetMCApp()->AddHitList(fHits);
195 //fNhits=0; //done in AliDetector(name,title)
197 SetDetectors(); // default to fOpt="All". This variable not written out.
199 fDetTypeSim = new AliITSDetTypeSim();
200 //PH SetMarkerColor(kRed);
201 if(!fLoader) MakeLoader(AliConfig::GetDefaultEventFolderName());
202 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
205 //______________________________________________________________________
207 // Default destructor for ITS.
208 // The default destructor of the AliITS class. In addition to deleting
209 // the AliITS class it deletes the memory pointed to by
210 // fIdSens, fIdName, fDetTypeSim and it's contents.
224 this->ClearModules();
227 }// end if fITSmodules!=0
229 delete[] fIdName; // Array of TStrings
231 Int_t size = AliITSgeomTGeo::GetNModules();
242 for(Int_t j=0; j<size; j++){
255 //______________________________________________________________________
257 // Initializer ITS after it has been built
258 // This routine initializes the AliITS class. It is intended to be
259 // called from the Init function in AliITSv?. Besides displaying a banner
260 // indicating that it has been called it initializes the array fIdSens
261 // and sets the default segmentation, response, digit and raw cluster
262 // classes therefore it should be called after a call to CreateGeometry.
271 if(gMC) for(i=0;i<fIdN;i++) fIdSens[i] = gMC->VolId(fIdName[i]);
274 //______________________________________________________________________
275 void AliITS::SetDefaults(){
276 // sets the default segmentation, response, digit and raw cluster classes.
283 AliInfoClass("AliITS::Setting Defaults");
285 Error("SetDefaults()","fDetTypeSim is 0!");
289 fDetTypeSim->SetDefaults();
290 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
293 //______________________________________________________________________
294 void AliITS::SetDefaultSimulation(){
295 // sets the default simulation.
303 Error("SetDefaultSimulation()","fDetTypeSim is 0!");
307 fDetTypeSim->SetDefaultSimulation();
308 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
313 //______________________________________________________________________
314 void AliITS::MakeBranch(Option_t* option){
315 // Creates Tree branches for the ITS.
317 // Option_t *option String of Tree types S,D, and/or R.
318 // const char *file String of the file name where these branches
319 // are to be stored. If blank then these branches
320 // are written to the same tree as the Hits were
327 Error("MakeBranch","fDetTypeSim is 0!");
331 Bool_t cH = (strstr(option,"H")!=0);
332 Bool_t cS = (strstr(option,"S")!=0);
333 Bool_t cD = (strstr(option,"D")!=0);
335 if(cH && (fHits == 0x0)) fHits = new TClonesArray("AliITShit", 1560);
336 AliDetector::MakeBranch(option);
338 if(cS) MakeBranchS(0);
339 if(cD) MakeBranchD(0);
343 //___________________________________________________________________
344 void AliITS::MakeBranchS(const char* fl){
346 // Creates Tree Branch for the ITS summable digits.
348 // cont char *fl File name where SDigits branch is to be written
349 // to. If blank it write the SDigits to the same
350 // file in which the Hits were found.
354 Error("MakeBranchS","fDetTypeSim is 0!");
356 Int_t buffersize = 4000;
359 // only one branch for SDigits.
360 snprintf(branchname,30,"%s",GetName());
362 if(fLoader->TreeS()){
363 TClonesArray* sdig = (TClonesArray*)fDetTypeSim->GetSDigits();
364 MakeBranchInTree(fLoader->TreeS(),branchname,&sdig,buffersize,fl);
367 //______________________________________________________________________
368 void AliITS::MakeBranchD(const char* file){
370 //Make branch for digits
372 Warning("MakeBranchD","fDetTypeSim is 0!");
375 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
376 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
377 MakeBranchInTreeD(fLoader->TreeD(),file);
380 //___________________________________________________________________
381 void AliITS:: MakeBranchInTreeD(TTree* treeD, const char* file){
382 // Creates Tree branches for the ITS.
385 Error("MakeBranchS","fDetTypeSim is 0!");
387 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
388 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
390 const Char_t *det[3] = {"SPD","SDD","SSD"};
391 const Char_t* digclass;
392 Int_t buffersize = 4000;
393 Char_t branchname[31];
395 if(!fDetTypeSim->GetDigits()){
396 fDetTypeSim->SetDigits(new TObjArray(fgkNTYPES));
398 for(Int_t i=0;i<fgkNTYPES;i++){
399 digclass = fDetTypeSim->GetDigitClassName(i);
400 TString classn = digclass;
401 if(!((fDetTypeSim->GetDigits())->At(i))){
402 (fDetTypeSim->GetDigits())->AddAt(new TClonesArray(classn.Data(),1000),i);
405 if(fgkNTYPES==3) snprintf(branchname,30,"%sDigits%s",GetName(),det[i]);
406 else sprintf(branchname,"%sDigits%d",GetName(),i+1);
407 TObjArray* dig = DigitsAddress(i);
408 if(GetDigits() && treeD) AliDetector::MakeBranchInTree(treeD,branchname, &dig,buffersize,file);
412 //______________________________________________________________________
413 void AliITS::SetTreeAddress(){
414 // Set branch address for the Trees.
423 Error("SetTreeAddress","fDetTypeSim is 0!");
427 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
428 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
430 TTree *treeS = fLoader->TreeS();
431 TTree *treeD = fLoader->TreeD();
432 if (fLoader->TreeH() && (fHits == 0x0)) {
433 fHits = new TClonesArray("AliITShit", 1560);
435 AliDetector::SetTreeAddress();
437 fDetTypeSim->SetTreeAddressS(treeS, (Char_t*)GetName());
438 fDetTypeSim->SetTreeAddressD(treeD, (Char_t*)GetName());
440 //______________________________________________________________________
441 void AliITS::AddHit(Int_t track, Int_t *vol, Float_t *hits){
443 // The function to add information to the AliITShit class. See the
444 // AliITShit class for a full description. This function allocates the
445 // necessary new space for the hit information and passes the variable
446 // track, and the pointers *vol and *hits to the AliITShit constructor
449 // Int_t track Track number which produced this hit.
450 // Int_t *vol Array of Integer Hit information. See AliITShit.h
451 // Float_t *hits Array of Floating Hit information. see AliITShit.h
456 TClonesArray &lhits = *fHits;
457 new(lhits[fNhits++]) AliITShit(fIshunt,track,vol,hits);
460 //______________________________________________________________________
461 void AliITS::FillModules(Int_t /* evnt */,Int_t bgrev,Int_t /* nmodules */,
462 Option_t *option, const char *filename){
463 // fill the modules with the sorted by module hits; add hits from
464 // background if option=Add.
466 static TTree *trH1; //Tree with background hits
467 static Bool_t first=kTRUE;
469 const char *addBgr = strstr(option,"Add");
473 file=new TFile(filename);
478 // Get Hits Tree header from file
479 if(trH1) delete trH1;
483 snprintf(treeName,20,"TreeH%d",bgrev);
484 trH1 = (TTree*)gDirectory->Get(treeName);
486 Error("FillModules","cannot find Hits Tree for event:%d",bgrev);
488 // Set branch addresses
491 FillModules(fLoader->TreeH(),0); // fill from this file's tree.
494 FillModules(trH1,10000000); // Default mask 10M.
495 TTree *fAli=fLoader->GetRunLoader()->TreeK();
498 fileAli =fAli->GetCurrentFile();
505 //______________________________________________________________________
506 void AliITS::FillModules(TTree *treeH, Int_t mask) {
507 // fill the modules with the sorted by module hits;
508 // can be called many times to do a merging
510 // TTree *treeH The tree containing the hits to be copied into
512 // Int_t mask The track number mask to indecate which file
513 // this hits came from.
521 AliError("Tree H is NULL");
524 Int_t lay,lad,det,index;
528 snprintf(branchname,20,"%s",GetName());
529 TBranch *branch = treeH->GetBranch(branchname);
531 Error("FillModules","%s branch in TreeH not found",branchname);
534 branch->SetAddress(&fHits);
535 Int_t nTracks =(Int_t) treeH->GetEntries();
537 for(iPrimTrack=0; iPrimTrack<nTracks; iPrimTrack++){
539 Int_t nBytes = treeH->GetEvent(iPrimTrack);
540 if (nBytes <= 0) continue;
541 Int_t nHits = fHits->GetEntriesFast();
542 for(h=0; h<nHits; h++){
543 itsHit = (AliITShit *)fHits->UncheckedAt(h);
544 itsHit->GetDetectorID(lay,lad,det);
546 index = GetITSgeom()->GetModuleIndex(lay,lad,det);
548 index=det-1; // This should not be used.
549 } // end if [You must have fITSgeom for this to work!]
550 mod = GetModule(index);
551 itsHit->SetTrack(itsHit->GetTrack()+mask); // Set track mask.
552 mod->AddHit(itsHit,iPrimTrack,h);
553 } // end loop over hits
554 } // end loop over tracks
557 //______________________________________________________________________
558 Bool_t AliITS::InitModules(Int_t size,Int_t &nmodules){
559 // Initialize the modules array.
561 // Int_t size Size of array of the number of modules to be
562 // created. If size <=0 then the number of modules
563 // is gotten from AliITSgeom class kept in fITSgeom.
565 // Int_t &nmodules The number of modules existing.
570 fITSmodules->Delete();
572 } // end fir fITSmoudles
575 Error("InitModules","fDetTypeSim is null!");
578 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
580 Int_t nl,indexMAX,index;
582 if(size<=0){ // default to using data stored in AliITSgeom
583 if(fDetTypeSim->GetITSgeom()==0) {
584 Error("InitModules","fITSgeom not defined");
586 } // end if fITSgeom==0
587 nl = fDetTypeSim->GetITSgeom()->GetNlayers();
588 indexMAX = fDetTypeSim->GetITSgeom()->GetIndexMax();
590 fITSmodules = new TObjArray(indexMAX);
591 for(index=0;index<indexMAX;index++){
592 fITSmodules->AddAt( new AliITSmodule(index),index);
595 fITSmodules = new TObjArray(size);
596 for(index=0;index<size;index++) {
597 fITSmodules->AddAt( new AliITSmodule(index),index);
604 //______________________________________________________________________
605 void AliITS::Hits2SDigits(){
606 // Standard Hits to summable Digits function.
614 Error("Hits2SDigits","fDetTypeSim is null!");
619 fLoader->LoadHits("read");
620 fLoader->LoadSDigits("recreate");
621 AliRunLoader* rl = fLoader->GetRunLoader();
622 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
623 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
625 for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
626 // Do the Hits to Digits operation. Use Standard input values.
627 // Event number from file, no background hit merging , use size from
628 // AliITSgeom class, option="All", input from this file only.
629 rl->GetEvent(iEvent);
630 if (!fLoader->TreeS()) fLoader->MakeTree("S");
633 HitsToPreDigits(iEvent,0,-1," ",fOpt," ");
636 fLoader->UnloadHits();
637 fLoader->UnloadSDigits();
640 //______________________________________________________________________
641 void AliITS::Hits2Digits(){
643 //Conversion from hits to digits
645 Error("Hits2SDigits","fDetTypeSim is 0!");
649 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
650 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
653 fLoader->LoadHits("read");
654 fLoader->LoadDigits("recreate");
655 AliRunLoader* rl = fLoader->GetRunLoader();
656 for (Int_t iEvent = 0; iEvent < rl->GetNumberOfEvents(); iEvent++) {
657 rl->GetEvent(iEvent);
658 if (!fLoader->TreeD()) fLoader->MakeTree("D");
661 HitsToDigits(iEvent,0,-1," ",fOpt," ");
664 fLoader->UnloadHits();
665 fLoader->UnloadDigits();
669 //______________________________________________________________________
670 void AliITS::HitsToDigits(Int_t evNumber,Int_t bgrev,Int_t size,
671 Option_t *option,Option_t *opt,
672 const char *filename){
673 // Keep galice.root for signal and name differently the file for
674 // background when add! otherwise the track info for signal will be lost !
675 // the condition below will disappear when the geom class will be
676 // initialized for all versions - for the moment it is only for v5 !
677 // 7 is the SDD beam test version.
679 // Int_t evnt Event to be processed.
680 // Int_t bgrev Background Hit tree number.
681 // Int_t nmodules Not used.
682 // Option_t *option String indicating if merging hits or not. To
683 // merge hits set equal to "Add". Otherwise no
684 // background hits are considered.
685 // Test_t *filename File name containing the background hits..
692 Error("HitsToDigits","fDetTypeSim is null!");
695 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
696 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
697 if(!GetITSgeom()) return; // need transformations to do digitization.
698 AliITSgeom *geom = GetITSgeom();
700 const char *all = strstr(opt,"All");
701 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
703 static Bool_t setDef=kTRUE;
704 if (setDef) SetDefaultSimulation();
708 InitModules(size,nmodules);
709 FillModules(evNumber,bgrev,nmodules,option,filename);
711 // Reset Fast-OR signals for this event
712 fDetTypeSim->ResetFOSignals();
714 AliITSsimulation *sim = 0;
715 AliITSmodule *mod = 0;
717 for(Int_t module=0;module<geom->GetIndexMax();module++){
718 id = geom->GetModuleType(module);
719 if (!all && !det[id]) continue;
720 sim = (AliITSsimulation*)fDetTypeSim->GetSimulationModel(id);
722 Error("HitsToDigits","The simulation class was not "
723 "instanciated for module %d type %s!",module,
724 geom->GetModuleTypeName(module));
727 mod = (AliITSmodule *)fITSmodules->At(module);
728 sim->DigitiseModule(mod,module,evNumber);
729 // fills all branches - wasted disk space
730 fLoader->TreeD()->Fill();
736 // Add random noise to FO signals
737 if (all || det[0]) { // SPD present
738 fDetTypeSim->ProcessNoiseForFastOr();
741 // Add Fast-OR signals to event (only one object per event)
742 if (all || det[0]) { // SPD present
743 fDetTypeSim->WriteFOSignals();
747 fLoader->TreeD()->GetEntries();
748 fLoader->TreeD()->AutoSave();
750 fLoader->TreeD()->Reset();
752 //_____________________________________________________________________
753 void AliITS::Hits2PreDigits(){
754 // Turn hits into SDigits
757 Error("Hits2SDigits","fDetTypeSim is 0!");
761 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
762 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
765 HitsToPreDigits(fLoader->GetRunLoader()->GetEventNumber(),
769 //______________________________________________________________________
770 void AliITS::HitsToPreDigits(Int_t evNumber,Int_t bgrev,Int_t size,
771 Option_t *option,Option_t *opt,
772 const char *filename){
773 // Keep galice.root for signal and name differently the file for
774 // background when add! otherwise the track info for signal will be lost !
775 // the condition below will disappear when the geom class will be
776 // initialized for all versions - for the moment it is only for v5 !
777 // 7 is the SDD beam test version.
779 // Int_t evnt Event to be processed.
780 // Int_t bgrev Background Hit tree number.
781 // Int_t nmodules Not used.
782 // Option_t *option String indicating if merging hits or not. To
783 // merge hits set equal to "Add". Otherwise no
784 // background hits are considered.
785 // Test_t *filename File name containing the background hits..
793 Error("HitsToPreDigits","fDetTypeSim is null!");
796 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
797 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
800 Error("HitsToPreDigits","fGeom is null!");
801 return; // need transformations to do digitization.
803 AliITSgeom *geom = GetITSgeom();
805 const char *all = strstr(opt,"All");
806 const char *det[3] = {strstr(opt,"SPD"),strstr(opt,"SDD"),
808 static Bool_t setDef=kTRUE;
809 if (setDef) SetDefaultSimulation();
813 InitModules(size,nmodules);
814 FillModules(evNumber,bgrev,nmodules,option,filename);
817 AliITSsimulation *sim = 0;
818 AliITSmodule *mod = 0;
820 for(module=0;module<geom->GetIndexMax();module++){
821 id = geom->GetModuleType(module);
822 if (!all && !det[id]) continue;
823 sim = (AliITSsimulation*)GetSimulationModel(id);
825 Error("HitsToPreDigits","The simulation class was not "
826 "instanciated for module %d type %s!",module,
827 geom->GetModuleTypeName(module));
830 mod = (AliITSmodule *)fITSmodules->At(module);
831 sim->SDigitiseModule(mod,module,evNumber);
832 // fills all branches - wasted disk space
833 fLoader->TreeS()->Fill();
834 fDetTypeSim->ResetSDigits();
840 fLoader->TreeS()->GetEntries();
841 fLoader->TreeS()->AutoSave();
842 fLoader->WriteSDigits("OVERWRITE");
844 fLoader->TreeS()->Reset();
847 //_____________________________________________________________________
848 void AliITS::HitsToFastRecPoints(Int_t evNumber,Int_t bgrev,Int_t size,
849 Option_t *opt0,Option_t *opt1,
851 // keep galice.root for signal and name differently the file for
852 // background when add! otherwise the track info for signal will be lost !
853 // the condition below will disappear when the geom class will be
854 // initialized for all versions - for the moment it is only for v5 !
856 // Int_t evnt Event to be processed.
857 // Int_t bgrev Background Hit tree number.
858 // Int_t size Size used by InitModules. See InitModules.
859 // Option_t *opt0 Option passed to FillModules. See FillModules.
860 // Option_t *opt1 String indicating if merging hits or not. To
861 // merge hits set equal to "Add". Otherwise no
862 // background hits are considered.
863 // Test_t *flnm File name containing the background hits..
872 Error("HitsToPreDigits","fGeom is null!");
873 return; // need transformations to do digitization.
875 AliITSgeom *geom = GetITSgeom();
877 AliITSLoader *pITSloader = (AliITSLoader*)fLoader;
879 const char *all = strstr(opt1,"All");
880 const char *det[3] ={strstr(opt1,"SPD"),strstr(opt1,"SDD"),
883 InitModules(size,nmodules);
884 FillModules(evNumber,bgrev,nmodules,opt0,flnm);
886 AliITSsimulation *sim = 0;
887 AliITSmodule *mod = 0;
890 TTree *lTR = pITSloader->TreeR();
892 pITSloader->MakeTree("R");
893 lTR = pITSloader->TreeR();
896 TClonesArray* ptarray = new TClonesArray("AliITSRecPoint",1000);
897 TBranch* branch = (TBranch*)lTR->Branch("ITSRecPointsF",&ptarray);
898 branch->SetAddress(&ptarray);
899 //m.b. : this change is nothing but a nice way to make sure
901 for(module=0;module<geom->GetIndexMax();module++){
902 id = geom->GetModuleType(module);
903 if (!all && !det[id]) continue;
904 sim = (AliITSsimulation*)GetSimulationModel(id);
906 Error("HitsToFastPoints","The simulation class was not "
907 "instantiated for module %d type %s!",module,
908 geom->GetModuleTypeName(module));
911 mod = (AliITSmodule *)fITSmodules->At(module);
912 sim->CreateFastRecPoints(mod,module,gRandom,ptarray);
918 fLoader->WriteRecPoints("OVERWRITE");
921 //_____________________________________________________________________
922 Int_t AliITS::Hits2Clusters(TTree *hTree, TTree *cTree) {
923 //------------------------------------------------------------
924 // This function creates ITS clusters
925 //------------------------------------------------------------
927 Error("HitsToPreDigits","fGeom is null!");
928 return 1; // need transformations to do digitization.
930 AliITSgeom *geom=GetITSgeom();
931 Int_t mmax=geom->GetIndexMax();
933 InitModules(-1,mmax);
934 FillModules(hTree,0);
936 TClonesArray *points = new TClonesArray("AliITSRecPoint",1000);
937 TBranch *branch=cTree->GetBranch("ITSRecPoints");
938 if (!branch) cTree->Branch("ITSRecPoints",&points);
939 else branch->SetAddress(&points);
941 AliITSsimulationFastPoints sim;
943 for (Int_t m=0; m<mmax; m++) {
944 AliITSmodule *mod=GetModule(m);
945 sim.CreateFastRecPoints(mod,m,gRandom,points);
946 ncl+=points->GetEntriesFast();
951 AliDebug(1,Form("Number of found fast clusters : %d",ncl));
959 //_____________________________________________________________________
960 void AliITS::CheckLabels(Int_t lab[3]) const {
961 //------------------------------------------------------------
962 // Tries to find mother's labels
963 //------------------------------------------------------------
965 if(lab[0]<0 && lab[1]<0 && lab[2]<0) return; // In case of no labels just exit
967 Int_t ntracks = gAlice->GetMCApp()->GetNtrack();
968 for (Int_t i=0;i<3;i++){
969 Int_t label = lab[i];
970 if (label>=0 && label<ntracks) {
971 TParticle *part=(TParticle*)gAlice->GetMCApp()->Particle(label);
972 if (part->P() < 0.005) {
973 Int_t m=part->GetFirstMother();
977 if (part->GetStatusCode()>0) {
987 //______________________________________________________________________
988 void AliITS::SDigitsToDigits(Option_t *opt){
989 // Standard Summable digits to Digits function.
995 AliError("fDetTypeSim is 0!");
999 const char *all = strstr(opt,"All");
1000 const char *det[3] ={strstr(opt,"SPD"),strstr(opt,"SDD"),
1003 // Reset Fast-OR signals for this event
1004 fDetTypeSim->ResetFOSignals();
1006 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
1008 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
1009 fDetTypeSim->SDigitsToDigits(opt,(Char_t*)GetName());
1011 // Add random noise to FO signals
1012 if (all || det[0]) { // SPD present
1013 fDetTypeSim->ProcessNoiseForFastOr();
1015 // Add Fast-OR signals to event (only one object per event)
1016 if (all || det[0]) { // SPD present
1017 fDetTypeSim->WriteFOSignals();
1021 //______________________________________________________________________
1022 void AliITS::ResetDigits(){
1023 // Reset number of digits and the digits array for the ITS detector.
1029 Error("ResetDigits","fDetTypeSim is 0!");
1033 fDetTypeSim->ResetDigits();
1037 //______________________________________________________________________
1038 void AliITS::ResetDigits(Int_t branch){
1039 // Reset number of digits and the digits array for this branch.
1046 Error("ResetDigits","fDetTypeSim is 0!");
1050 fDetTypeSim->ResetDigits(branch);
1053 //______________________________________________________________________
1054 void AliITS::AddSumDigit(AliITSpListItem &sdig){
1055 // Adds the a module full of summable digits to the summable digits tree.
1057 // AliITSpListItem &sdig SDigit to be added to SDigits tree.
1064 Error("AddSumDigit","fDetTypeSim is 0!");
1067 fDetTypeSim->AddSumDigit(sdig);
1070 //______________________________________________________________________
1071 void AliITS::AddSimDigit(Int_t branch, AliITSdigit *d){
1072 // Add a simulated digit.
1074 // Int_t id Detector type number.
1075 // AliITSdigit *d Digit to be added to the Digits Tree. See
1083 Error("AddSimDigit","fDetTypeSim is 0!");
1086 fDetTypeSim->AddSimDigit(branch,d);
1089 //______________________________________________________________________
1090 void AliITS::AddSimDigit(Int_t branch,Float_t phys,Int_t *digits,Int_t *tracks,
1091 Int_t *hits,Float_t *charges, Int_t sigexpanded){
1092 // Add a simulated digit to the list.
1094 // Int_t id Detector type number.
1095 // Float_t phys Physics indicator. See AliITSdigits.h
1096 // Int_t *digits Integer array containing the digits info. See
1098 // Int_t *tracks Integer array [AliITSdigitS?D::GetNTracks()]
1099 // containing the track numbers that contributed to
1101 // Int_t *hits Integer array [AliITSdigitS?D::GetNTracks()]
1102 // containing the hit numbers, from AliITSmodule, that
1103 // contributed to this digit.
1104 // Float_t *charge Floating point array of the signals contributed
1105 // to this digit by each track.
1112 Error("AddSimDigit","fDetTypeSim is 0!");
1115 fDetTypeSim->AddSimDigit(branch,phys,digits,tracks,hits,charges,sigexpanded);
1118 //______________________________________________________________________
1119 void AliITS::Digits2Raw(){
1120 // convert digits of the current event to raw data
1123 Error("Digits2Raw","fDetTypeSim is 0!");
1126 fDetTypeSim->SetLoader((AliITSLoader*)fLoader);
1128 if(fSimuParam) fDetTypeSim->SetSimuParam(fSimuParam);
1129 fDetTypeSim->GetLoader()->LoadDigits();
1130 TTree* digits = fDetTypeSim->GetLoader()->TreeD();
1132 Error("Digits2Raw", "no digits tree");
1135 fDetTypeSim->SetTreeAddressD(digits,(Char_t*)GetName());
1137 // Get the FO signals for this event
1138 AliITSFOSignalsSPD* foSignals = NULL;
1139 AliRunLoader* runLoader = AliRunLoader::Instance();
1140 AliITSLoader* itsLoader = (AliITSLoader*) runLoader->GetLoader("ITSLoader");
1142 AliError("ITS loader is NULL.");
1145 if(!itsLoader->TreeD()) AliError(" !!! No TreeD available !!!");
1146 foSignals = (AliITSFOSignalsSPD*)itsLoader->TreeD()->GetUserInfo()->FindObject("AliITSFOSignalsSPD");
1147 if(!foSignals) AliError("FO signals not retrieved");
1150 Bool_t deleteFOsignalsLater = kFALSE;
1152 AliError("FO signals not available. No FO bits will be written.");
1153 foSignals = new AliITSFOSignalsSPD(); // make a temporary dummy signals object
1154 deleteFOsignalsLater = kTRUE;
1158 AliITSDDLModuleMapSDD* ddlsdd=fDetTypeSim->GetDDLModuleMapSDD();
1159 Char_t rawSDD=fDetTypeSim->GetSimuParam()->GetSDDRawDataFormat();
1160 AliITSDDLRawData rawWriter;
1162 rawWriter.SetSDDRawFormat(rawSDD);
1166 // 2: txt files with digits
1167 //BE CAREFUL, verbose level 2 MUST be used only for debugging and
1168 //it is highly suggested to use this mode only for debugging digits files
1169 //reasonably small, because otherwise the size of the txt files can reach
1170 //quickly several MB wasting time and disk space.
1171 rawWriter.SetVerbose(0);
1173 //SILICON PIXEL DETECTOR
1174 AliDebug(1,"Formatting raw data for SPD");
1175 rawWriter.RawDataSPD(digits->GetBranch("ITSDigitsSPD"),foSignals);
1176 if(deleteFOsignalsLater) delete foSignals;
1178 //SILICON DRIFT DETECTOR
1179 AliDebug(1,Form("Formatting raw data for SDD - Format code =%d",rawSDD));
1180 rawWriter.RawDataSDD(digits->GetBranch("ITSDigitsSDD"),ddlsdd);
1182 //SILICON STRIP DETECTOR
1183 AliDebug(1,"Formatting raw data for SSD");
1184 rawWriter.RawDataSSD(digits->GetBranch("ITSDigitsSSD"));
1186 fLoader->UnloadDigits();
1188 //______________________________________________________________________
1189 AliLoader* AliITS::MakeLoader(const char* topfoldername){
1190 //builds ITSgetter (AliLoader type)
1191 //if detector wants to use castomized getter, it must overload this method
1193 AliDebug(1,Form("Creating AliITSLoader. Top folder is %s.",
1195 fLoader = new AliITSLoader(GetName(),topfoldername);
1198 //______________________________________________________________________
1199 Bool_t AliITS::Raw2SDigits(AliRawReader* rawReader)
1202 // Converts RAW data to SDigits
1207 Int_t size = AliITSgeomTGeo::GetNModules();
1209 fModA = new TClonesArray*[size];
1210 for (Int_t mod = 0; mod < size; mod++) fModA[mod] = new TClonesArray("AliITSpListItem", 10000);
1212 AliLoader* loader = (AliRunLoader::Instance())->GetLoader("ITSLoader");
1214 Error("Open","Can not get ITS loader from Run Loader");
1219 tree = loader->TreeS();
1221 loader->MakeTree("S");
1222 tree = loader->TreeS();
1225 // Array for SDigits
1228 fpSDigits = new TClonesArray("AliITSpListItem",10000);
1230 TClonesArray& aSDigits = *fpSDigits;
1231 Int_t bufsize = 32000;
1232 tree->Branch("ITS", &fpSDigits, bufsize);
1237 AliITSsegmentationSPD* segSPD = (AliITSsegmentationSPD*) fDetTypeSim->GetSegmentationModel(0);
1239 AliWarning("Set AliITS defaults");
1241 segSPD = (AliITSsegmentationSPD*) fDetTypeSim->GetSegmentationModel(0);
1243 npx = segSPD->Npx();
1244 Double_t thr, sigma;
1246 AliITSRawStreamSPD inputSPD(rawReader);
1248 Bool_t next = inputSPD.Next();
1251 Int_t module = inputSPD.GetModuleID();
1252 Int_t column = inputSPD.GetColumn();
1253 Int_t row = inputSPD.GetRow();
1254 Int_t index = npx * column + row;
1256 if (module >= size) continue;
1258 last = (fModA[module])->GetEntries();
1259 TClonesArray& dum = *fModA[module];
1260 fDetTypeSim->GetSimuParam()->SPDThresholds(module,thr,sigma);
1262 new (dum[last]) AliITSpListItem(-1, -1, module, index, thr);
1269 AliITSsegmentationSDD* segSDD = (AliITSsegmentationSDD*) fDetTypeSim->GetSegmentationModel(1);
1270 npx = segSDD->Npx();
1271 Int_t scalef=AliITSsimulationSDD::ScaleFourier(segSDD);
1272 Int_t firstSDD=AliITSgeomTGeo::GetModuleIndex(3,1,1);
1273 Int_t firstSSD=AliITSgeomTGeo::GetModuleIndex(5,1,1);
1275 AliITSRawStream* inputSDD=AliITSRawStreamSDD::CreateRawStreamSDD(rawReader);
1276 for(Int_t iMod=firstSDD; iMod<firstSSD; iMod++){
1277 AliITSCalibrationSDD* cal = (AliITSCalibrationSDD*)fDetTypeSim->GetCalibrationModel(iMod);
1278 Bool_t isZeroSupp=cal->GetZeroSupp();
1280 for(Int_t iSid=0; iSid<2; iSid++) inputSDD->SetZeroSuppLowThreshold(iMod-firstSDD,iSid,cal->GetZSLowThreshold(iSid));
1282 for(Int_t iSid=0; iSid<2; iSid++) inputSDD->SetZeroSuppLowThreshold(iMod-firstSDD,iSid,0);
1286 AliITSDDLModuleMapSDD* ddlmap=fDetTypeSim->GetDDLModuleMapSDD();
1287 inputSDD->SetDDLModuleMap(ddlmap);
1288 while(inputSDD->Next()){
1289 if(inputSDD->IsCompletedModule()==kFALSE &&
1290 inputSDD->IsCompletedDDL()==kFALSE){
1292 Int_t module = inputSDD->GetModuleID();
1293 Int_t anode = inputSDD->GetCoord1()+segSDD->NpzHalf()*inputSDD->GetChannel();
1294 Int_t time = inputSDD->GetCoord2();
1295 Int_t signal10 = inputSDD->GetSignal();
1296 Int_t index = AliITSpList::GetIndex(anode,time,scalef*npx);
1298 if (module >= size) continue;
1299 last = fModA[module]->GetEntries();
1300 TClonesArray& dum = *fModA[module];
1301 new (dum[last]) AliITSpListItem(-1, -1, module, index, Double_t(signal10));
1302 ((AliITSpListItem*) dum.At(last))->AddSignalAfterElect(module, index, Double_t(signal10));
1311 AliITSsegmentationSSD* segSSD = (AliITSsegmentationSSD*) fDetTypeSim->GetSegmentationModel(2);
1312 npx = segSSD->Npx();
1313 AliITSRawStreamSSD inputSSD(rawReader);
1315 Bool_t next = inputSSD.Next();
1318 Int_t module = inputSSD.GetModuleID();
1319 if(module<0)AliError(Form("Invalid SSD module %d \n",module));
1320 if(module<0)continue;
1321 Int_t side = inputSSD.GetSideFlag();
1322 Int_t strip = inputSSD.GetStrip();
1323 Int_t signal = inputSSD.GetSignal();
1324 Int_t index = npx * side + strip;
1326 if (module >= size) continue;
1328 last = fModA[module]->GetEntries();
1329 TClonesArray& dum = *fModA[module];
1330 new (dum[last]) AliITSpListItem(-1, -1, module, index, Double_t(signal));
1333 AliITSpListItem* sdig = 0;
1335 Int_t firstssd = GetITSgeom()->GetStartDet(kSSD);
1336 Double_t adcToEv = 1.;
1337 for (Int_t mod = 0; mod < size; mod++)
1340 AliITSCalibrationSSD* calssd = (AliITSCalibrationSSD*)fDetTypeSim->GetCalibrationModel(mod);
1341 adcToEv = 1./calssd->GetSSDDEvToADC(1.);
1343 Int_t nsdig = fModA[mod]->GetEntries();
1344 for (Int_t ie = 0; ie < nsdig; ie++) {
1345 sdig = (AliITSpListItem*) (fModA[mod]->At(ie));
1346 Double_t digsig = sdig->GetSignal();
1347 if(mod>=firstssd) digsig*=adcToEv; // for SSD: convert back charge from ADC to electron
1348 new (aSDigits[ie]) AliITSpListItem(-1, -1, mod, sdig->GetIndex(), digsig);
1349 Float_t sig = sdig->GetSignalAfterElect();
1350 if(mod>=firstssd) sig*=adcToEv;
1352 sdig = (AliITSpListItem*)aSDigits[ie];
1353 sdig->AddSignalAfterElect(mod, sdig->GetIndex(), Double_t(sig));
1359 fModA[mod]->Clear();
1361 loader->WriteSDigits("OVERWRITE");
1365 //______________________________________________________________________
1366 void AliITS::UpdateInternalGeometry(){
1367 //reads new geometry from TGeo
1368 // AliDebug(1,"Delete ITSgeom and create a new one reading TGeo");
1370 AliITSVersion_t version = (AliITSVersion_t)IsVersion();
1372 AliITSInitGeometry initgeom;
1373 AliITSgeom* geom = initgeom.CreateAliITSgeom(version,minor);
1376 //______________________________________________________________________
1377 AliTriggerDetector* AliITS::CreateTriggerDetector() const {
1378 // create an AliITSTrigger object (and set trigger conditions as input)
1379 return new AliITSTrigger(fDetTypeSim->GetTriggerConditions());