AliRoot Core  edcc906 (edcc906)
AliMUONSDigitizerV2.cxx
Go to the documentation of this file.
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15 
16 // $Id$
17 
18 #include <TClonesArray.h>
19 #include <TParticle.h>
20 #include <TTree.h>
21 
22 #include "AliMUONSDigitizerV2.h"
23 
24 #include "AliMUON.h"
25 #include "AliMUONChamber.h"
26 #include "AliMUONVDigit.h"
27 #include "AliMUONHit.h"
28 #include "AliMUONVDigitStore.h"
29 #include "AliMUONVHitStore.h"
30 #include "AliMUONResponseTrigger.h"
31 #include "AliMUONConstants.h"
32 
33 #include "AliMpCDB.h"
34 #include "AliMpDEManager.h"
35 
36 #include "AliLog.h"
37 #include "AliCDBManager.h"
38 #include "AliLoader.h"
39 #include "AliRun.h"
40 #include "AliRunLoader.h"
41 
42 #include "AliHeader.h"
44 
45 //-----------------------------------------------------------------------------
59 //-----------------------------------------------------------------------------
60 
61 ClassImp(AliMUONSDigitizerV2)
62 
63 Float_t AliMUONSDigitizerV2::fgkMaxIntTime = 10.0;
64 Float_t AliMUONSDigitizerV2::fgkMaxPosTimeDif = 1.22E-6;
65 Float_t AliMUONSDigitizerV2::fgkMaxNegTimeDif = -3.5E-6;
66 Float_t AliMUONSDigitizerV2::fgkMinTimeDif = 25E-9;
67 
68 //_____________________________________________________________________________
70 : TNamed("AliMUONSDigitizerV2","From Hits to SDigits for MUON")
71 {
75 
76  // Load mapping
77  if ( ! AliMpCDB::LoadMpSegmentation() ) {
78  AliFatal("Could not access mapping from OCDB !");
79  }
80 }
81 
82 //_____________________________________________________________________________
84 {
88 }
89 
90 //_____________________________________________________________________________
91 void
93 {
102 
103  AliDebug(1,"");
104 
105  AliRunLoader* runLoader = AliRunLoader::Instance();
106  AliLoader* loader = runLoader->GetDetectorLoader("MUON");
107 
108  loader->LoadHits("READ");
109 
110  AliMUON* muon = static_cast<AliMUON*>(gAlice->GetModule("MUON"));
111 
112  Int_t nofEvents(runLoader->GetNumberOfEvents());
113 
114  TString classname = muon->DigitStoreClassName();
115 
116  AliMUONVDigitStore* sDigitStore = AliMUONVDigitStore::Create(classname.Data());
117 
118  if (!sDigitStore)
119  {
120  AliFatal(Form("Could not create digitstore of class %s",classname.Data()));
121  }
122 
123  AliDebug(1,Form("Will use digitStore of type %s",sDigitStore->ClassName()));
124 
125  // average arrival time to chambers, for pileup studies
126 
127  for ( Int_t iEvent = 0; iEvent < nofEvents; ++iEvent )
128  {
129  // Loop over events.
130  TObjArray tdlist;
131  tdlist.SetOwner(kTRUE);
132 
133  AliDebug(1,Form("iEvent=%d",iEvent));
134  runLoader->GetEvent(iEvent);
135 
136  // for pile up studies
137  float t0=fgkMaxIntTime; int aa=0;
138  AliHeader* header = runLoader->GetHeader();
139  AliGenCocktailEventHeader* cocktailHeader =
140  dynamic_cast<AliGenCocktailEventHeader*>(header->GenEventHeader());
141  if (cocktailHeader) {
142  AliGenCocktailEventHeader* genEventHeader = (AliGenCocktailEventHeader*) (header->GenEventHeader());
143  TList* headers = genEventHeader->GetHeaders();
144  TIter nextH(headers);
145  AliGenEventHeader *entry;
146  while((entry = (AliGenEventHeader*)nextH())) {
147  float t = entry->InteractionTime();
148  if (TMath::Abs(t)<TMath::Abs(t0)) t0 = t;
149  aa++;
150  }
151  } else {
152  AliGenEventHeader* evtHeader =
153  (AliGenEventHeader*)(header->GenEventHeader());
154  if (evtHeader)
155  {
156  float t = evtHeader->InteractionTime();
157  if (TMath::Abs(t)<TMath::Abs(t0)) t0 = t;
158  aa++;
159  }
160  else
161  {
162  // some generators may not offer a header, if which
163  // case we cannot get the interaction time, so we assume zero
164  t0 = 0.;
165  }
166  }
167 
168  loader->MakeSDigitsContainer();
169 
170  TTree* treeS = loader->TreeS();
171 
172  if ( !treeS )
173  {
174  AliFatal("");
175  }
176 
177  sDigitStore->Connect(*treeS);
178 
179  TTree* treeH = loader->TreeH();
180 
181  AliMUONVHitStore* hitStore = AliMUONVHitStore::Create(*treeH);
182  hitStore->Connect(*treeH);
183 
184  Long64_t nofTracks = treeH->GetEntries();
185 
186  for ( Long64_t iTrack = 0; iTrack < nofTracks; ++iTrack )
187  {
188  // Loop over the tracks of this event.
189  treeH->GetEvent(iTrack);
190 
191  AliMUONHit* hit;
192  TIter next(hitStore->CreateIterator());
193  Int_t ihit(0);
194 
195  while ( ( hit = static_cast<AliMUONHit*>(next()) ) )
196  {
197  Int_t chamberId = hit->Chamber()-1;
198  Float_t age = hit->Age()-t0;
199 
200  AliMUONChamber& chamber = muon->Chamber(chamberId);
201  AliMUONResponse* response = chamber.ResponseModel();
202 
203  // This is the heart of this method : the dis-integration
204  TList digits;
205  if (aa>1){ // if there are pileup events
206  Float_t chamberTime = AliMUONConstants::AverageChamberT(chamberId);
207  Float_t timeDif=age-chamberTime;
208  if (timeDif>fgkMaxPosTimeDif || timeDif<fgkMaxNegTimeDif) {
209  continue;
210  }
211  if(TMath::Abs(timeDif)>fgkMinTimeDif){
212  response->DisIntegrate(*hit,digits,timeDif);
213  }
214  else{
215  response->DisIntegrate(*hit,digits,0.);
216  }
217  }
218  else{
219  response->DisIntegrate(*hit,digits,0.);
220  }
221 
222  TIter nextd(&digits);
223  AliMUONVDigit* d;
224  while ( ( d = (AliMUONVDigit*)nextd() ) )
225  {
226  // Update some sdigit information that could not be known
227  // by the DisIntegrate method
228  d->SetHit(ihit);
229  d->SetTime(age);
230  d->AddTrack(hit->GetTrack(),d->Charge());
231  tdlist.Add(d);
232  }
233  ++ihit;
234  }
235  hitStore->Clear();
236  } // end of loop on tracks within an event
237 
238  TIter next(&tdlist);
239  AliMUONVDigit* d;
240 
241  while ( ( d = static_cast<AliMUONVDigit*>(next()) ) )
242  {
243  d->ChargeInFC(kTRUE);
244 
245  AliMUONVDigit* added = sDigitStore->Add(*d,AliMUONVDigitStore::kMerge);
246  if (!added)
247  {
248  AliError("Could not add digit to digitStore");
249  }
250  }
251 
252  treeS->Fill();
253 
254  loader->WriteSDigits("OVERWRITE");
255 
256  sDigitStore->Clear();
257 
258  loader->UnloadSDigits();
259 
260  delete hitStore;
261 
262  } // loop on events
263 
264  loader->UnloadHits();
265 
266  delete sDigitStore;
267 
268 }
static AliRunLoader * Instance()
Definition: AliRunLoader.h:176
const TString DigitStoreClassName() const
Return digit store class name.
Definition: AliMUON.h:144
Interface for a digit container.
Virtual store to hold digit.
#define TObjArray
MUON tracking chamber class.
Int_t GetTrack() const
Definition: AliHit.h:19
virtual void ChargeInFC(Bool_t value=kTRUE)=0
Set the unit value (see note 1 in AliMUONVDigit.cxx)
Float_t Age() const
Return Particle Age.
Definition: AliMUONHit.h:53
AliLoader * GetDetectorLoader(const char *detname)
virtual AliMUONChamber & Chamber(Int_t id)
Return reference to Chamber id.
Definition: AliMUON.h:135
AliDetector class for MUON subsystem providing simulation data management.
Definition: AliMUON.h:37
TTree * TreeS() const
Definition: AliLoader.h:85
static Float_t fgkMaxIntTime
maximum time of interaction
virtual void AddTrack(Int_t, Float_t)
Add a track (and its charge) to the list of tracks we handle.
virtual Float_t Charge() const =0
The charge of this digit, calibrated or not depending on IsCalibrated()
MUON SDigitizer (from Hits to SDigits).
virtual Float_t InteractionTime() const
static Float_t fgkMaxPosTimeDif
maximum event time after the triggered event for a hit to be digitized
virtual void MakeSDigitsContainer() const
Definition: AliLoader.h:74
virtual Int_t WriteSDigits(Option_t *opt="") const
Definition: AliLoader.cxx:380
static Float_t fgkMaxNegTimeDif
maximum event time before the triggered event for a hit to be digitized
virtual AliMUONVDigitStore * Create() const =0
Create an (empty) object of the same concrete class as *this.
virtual AliMUONVStore * Create() const =0
Create an empty copy of this.
AliModule * GetModule(const char *name) const
Definition: AliRun.cxx:191
AliHeader * GetHeader() const
static Float_t fgkMinTimeDif
minimum time difference for the reduction factor to be applied
Int_t Chamber() const
Definition: AliMUONHit.cxx:160
virtual void SetTime(Float_t)
Set hit age.
AliRun * gAlice
Definition: AliRun.cxx:62
virtual AliMUONResponse *& ResponseModel()
Get pointer to response model.
Int_t GetNumberOfEvents()
static Float_t AverageChamberT(Int_t i)
Return average arrival time to chamber i.
virtual void Digitize(Option_t *opt="")
#define AliFatal(message)
Definition: AliLog.h:640
Int_t LoadHits(Option_t *opt="")
Definition: AliLoader.h:96
Int_t GetEvent(Int_t evno)
#define AliDebug(logLevel, message)
Definition: AliLog.h:300
virtual AliGenEventHeader * GenEventHeader() const
Definition: AliHeader.cxx:244
Chamber response base class.
AliMUON * muon()
MonteCarlo hit.
Definition: AliMUONHit.h:24
void UnloadHits() const
Definition: AliLoader.h:129
TTree * TreeH() const
Definition: AliLoader.h:83
virtual TIterator * CreateIterator() const =0
Return an iterator to loop over hits.
ABC of a MUON digit.
Definition: AliMUONVDigit.h:18
void UnloadSDigits() const
Definition: AliLoader.h:130
virtual void Clear(Option_t *opt="")=0
Clear ourselves (i.e. Reset)
#define AliError(message)
Definition: AliLog.h:591
virtual void DisIntegrate(const AliMUONHit &hit, TList &digits, Float_t timeDif)
for(Int_t itree=0;itree< arrInputTreesDistortionCalib->GetEntriesFast();++itree)
virtual Bool_t Connect(TTree &tree, Bool_t alone=kTRUE) const
Connect us to a TTree (only valid if CanConnect()==kTRUE)
static Bool_t LoadMpSegmentation(Bool_t warn=false)
Definition: AliMpCDB.cxx:113
virtual Bool_t Add(TObject *object)
Add an object, if it is of the right class.
virtual void SetHit(Int_t)
Set the hit number.