Line data Source code
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"
43 : #include "AliGenCocktailEventHeader.h"
44 :
45 : //-----------------------------------------------------------------------------
46 : /// The sdigitizer performs the transformation from hits (energy deposits by
47 : /// the transport code) to sdigits (equivalent of charges on pad).
48 : ///
49 : /// It does so by converting the energy deposit into a charge and then spreading
50 : /// this charge over several pads, according to the response function (a
51 : /// Mathieson distribution, basically).
52 : ///
53 : /// See also AliMUONResponseV0, which is doing the real job (in DisIntegrate
54 : /// method), while this sdigitizer is just "steering" the process.
55 : ///
56 : /// Please note that we do *not* merge sdigits after creation, which means
57 : /// that after sdigitization, a given pad can have several sdigits. This
58 : /// merging is taken care of later on by the digitizer(V3).
59 : //-----------------------------------------------------------------------------
60 :
61 16 : 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 : //_____________________________________________________________________________
69 : AliMUONSDigitizerV2::AliMUONSDigitizerV2()
70 1 : : TNamed("AliMUONSDigitizerV2","From Hits to SDigits for MUON")
71 5 : {
72 : ///
73 : /// ctor.
74 : ///
75 :
76 : // Load mapping
77 2 : if ( ! AliMpCDB::LoadMpSegmentation() ) {
78 0 : AliFatal("Could not access mapping from OCDB !");
79 : }
80 2 : }
81 :
82 : //_____________________________________________________________________________
83 : AliMUONSDigitizerV2::~AliMUONSDigitizerV2()
84 2 : {
85 : ///
86 : /// dtor.
87 : ///
88 3 : }
89 :
90 : //_____________________________________________________________________________
91 : void
92 : AliMUONSDigitizerV2::Digitize(Option_t*)
93 : {
94 : ///
95 : /// Go from hits to sdigits.
96 : ///
97 : /// In the code below, apart from the loop itself (which look complicated
98 : /// but is really only a loop on each hit in the input file) the main
99 : /// work is done in AliMUONResponse::DisIntegrate method, which converts
100 : /// a single hit in (possibly) several sdigits.
101 : ///
102 :
103 4 : AliDebug(1,"");
104 :
105 1 : AliRunLoader* runLoader = AliRunLoader::Instance();
106 1 : AliLoader* loader = runLoader->GetDetectorLoader("MUON");
107 :
108 1 : loader->LoadHits("READ");
109 :
110 1 : AliMUON* muon = static_cast<AliMUON*>(gAlice->GetModule("MUON"));
111 :
112 1 : Int_t nofEvents(runLoader->GetNumberOfEvents());
113 :
114 1 : TString classname = muon->DigitStoreClassName();
115 :
116 2 : AliMUONVDigitStore* sDigitStore = AliMUONVDigitStore::Create(classname.Data());
117 :
118 1 : if (!sDigitStore)
119 : {
120 0 : AliFatal(Form("Could not create digitstore of class %s",classname.Data()));
121 : }
122 :
123 5 : AliDebug(1,Form("Will use digitStore of type %s",sDigitStore->ClassName()));
124 :
125 : // average arrival time to chambers, for pileup studies
126 :
127 10 : for ( Int_t iEvent = 0; iEvent < nofEvents; ++iEvent )
128 : {
129 : // Loop over events.
130 4 : TObjArray tdlist;
131 4 : tdlist.SetOwner(kTRUE);
132 :
133 20 : AliDebug(1,Form("iEvent=%d",iEvent));
134 4 : runLoader->GetEvent(iEvent);
135 :
136 : // for pile up studies
137 4 : float t0=fgkMaxIntTime; int aa=0;
138 4 : AliHeader* header = runLoader->GetHeader();
139 : AliGenCocktailEventHeader* cocktailHeader =
140 16 : dynamic_cast<AliGenCocktailEventHeader*>(header->GenEventHeader());
141 4 : if (cocktailHeader) {
142 8 : AliGenCocktailEventHeader* genEventHeader = (AliGenCocktailEventHeader*) (header->GenEventHeader());
143 4 : TList* headers = genEventHeader->GetHeaders();
144 4 : TIter nextH(headers);
145 : AliGenEventHeader *entry;
146 60 : while((entry = (AliGenEventHeader*)nextH())) {
147 16 : float t = entry->InteractionTime();
148 20 : if (TMath::Abs(t)<TMath::Abs(t0)) t0 = t;
149 16 : aa++;
150 : }
151 4 : } else {
152 : AliGenEventHeader* evtHeader =
153 0 : (AliGenEventHeader*)(header->GenEventHeader());
154 0 : if (evtHeader)
155 : {
156 0 : float t = evtHeader->InteractionTime();
157 0 : if (TMath::Abs(t)<TMath::Abs(t0)) t0 = t;
158 : aa++;
159 0 : }
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 4 : loader->MakeSDigitsContainer();
169 :
170 4 : TTree* treeS = loader->TreeS();
171 :
172 4 : if ( !treeS )
173 : {
174 0 : AliFatal("");
175 : }
176 :
177 4 : sDigitStore->Connect(*treeS);
178 :
179 4 : TTree* treeH = loader->TreeH();
180 :
181 4 : AliMUONVHitStore* hitStore = AliMUONVHitStore::Create(*treeH);
182 4 : hitStore->Connect(*treeH);
183 :
184 4 : Long64_t nofTracks = treeH->GetEntries();
185 :
186 232 : for ( Long64_t iTrack = 0; iTrack < nofTracks; ++iTrack )
187 : {
188 : // Loop over the tracks of this event.
189 112 : treeH->GetEvent(iTrack);
190 :
191 : AliMUONHit* hit;
192 224 : TIter next(hitStore->CreateIterator());
193 : Int_t ihit(0);
194 :
195 566 : while ( ( hit = static_cast<AliMUONHit*>(next()) ) )
196 : {
197 230 : Int_t chamberId = hit->Chamber()-1;
198 115 : Float_t age = hit->Age()-t0;
199 :
200 115 : AliMUONChamber& chamber = muon->Chamber(chamberId);
201 230 : AliMUONResponse* response = chamber.ResponseModel();
202 :
203 : // This is the heart of this method : the dis-integration
204 115 : TList digits;
205 115 : if (aa>1){ // if there are pileup events
206 115 : Float_t chamberTime = AliMUONConstants::AverageChamberT(chamberId);
207 115 : Float_t timeDif=age-chamberTime;
208 230 : if (timeDif>fgkMaxPosTimeDif || timeDif<fgkMaxNegTimeDif) {
209 0 : continue;
210 : }
211 230 : if(TMath::Abs(timeDif)>fgkMinTimeDif){
212 115 : response->DisIntegrate(*hit,digits,timeDif);
213 : }
214 : else{
215 115 : response->DisIntegrate(*hit,digits,0.);
216 : }
217 115 : }
218 : else{
219 0 : response->DisIntegrate(*hit,digits,0.);
220 : }
221 :
222 115 : TIter nextd(&digits);
223 : AliMUONVDigit* d;
224 2293 : while ( ( d = (AliMUONVDigit*)nextd() ) )
225 : {
226 : // Update some sdigit information that could not be known
227 : // by the DisIntegrate method
228 974 : d->SetHit(ihit);
229 974 : d->SetTime(age);
230 2922 : d->AddTrack(hit->GetTrack(),d->Charge());
231 974 : tdlist.Add(d);
232 : }
233 115 : ++ihit;
234 230 : }
235 112 : hitStore->Clear();
236 112 : } // end of loop on tracks within an event
237 :
238 4 : TIter next(&tdlist);
239 : AliMUONVDigit* d;
240 :
241 2934 : while ( ( d = static_cast<AliMUONVDigit*>(next()) ) )
242 : {
243 974 : d->ChargeInFC(kTRUE);
244 :
245 974 : AliMUONVDigit* added = sDigitStore->Add(*d,AliMUONVDigitStore::kMerge);
246 974 : if (!added)
247 : {
248 0 : AliError("Could not add digit to digitStore");
249 : }
250 : }
251 :
252 4 : treeS->Fill();
253 :
254 4 : loader->WriteSDigits("OVERWRITE");
255 :
256 4 : sDigitStore->Clear();
257 :
258 4 : loader->UnloadSDigits();
259 :
260 8 : delete hitStore;
261 :
262 4 : } // loop on events
263 :
264 1 : loader->UnloadHits();
265 :
266 2 : delete sDigitStore;
267 :
268 1 : }
|