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