]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONSDigitizerV2.cxx
a/ AliTRDCalibTask.cxx .h: one histo more to quantify the event selection if any...
[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
21 #include "AliMUONSDigitizerV2.h"
22
23 #include "AliMUON.h"
24 #include "AliMUONChamber.h"
25 #include "AliMUONVDigit.h"
26 #include "AliMUONHit.h"
27 #include "AliMUONVDigitStore.h"
28 #include "AliMUONVHitStore.h"
29 #include "AliMUONCalibrationData.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 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 : TTask("AliMUONSDigitizerV2","From Hits to SDigits for MUON")
71 {
72   ///
73   /// ctor.
74   ///
75
76   // Load mapping
77   if ( ! AliMpCDB::LoadMpSegmentation() ) {
78     AliFatal("Could not access mapping from OCDB !");
79   }
80 }
81
82 //_____________________________________________________________________________
83 AliMUONSDigitizerV2::~AliMUONSDigitizerV2()
84 {
85   ///
86   /// dtor.
87   ///
88 }
89
90 //_____________________________________________________________________________
91 void
92 AliMUONSDigitizerV2::Exec(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   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   AliMUONCalibrationData *calibrationData = 0x0;
113
114   if(muon->GetTriggerEffCells()){
115     Int_t runnumber = AliCDBManager::Instance()->GetRun();
116     calibrationData = new AliMUONCalibrationData(runnumber);
117     for (Int_t chamber = 10; chamber < 14; chamber++) {
118       ((AliMUONResponseTrigger *) (muon->Chamber(chamber).ResponseModel()))->InitTriggerEfficiency(calibrationData->TriggerEfficiency()); // Init trigger efficiency
119     }
120   }
121   
122   Int_t nofEvents(runLoader->GetNumberOfEvents());
123   
124   TString classname = muon->DigitStoreClassName();
125   
126   AliMUONVDigitStore* sDigitStore = AliMUONVDigitStore::Create(classname.Data());
127   
128   if (!sDigitStore)
129   {
130     AliFatal(Form("Could not create digitstore of class %s",classname.Data()));
131   }
132   
133   AliDebug(1,Form("Will use digitStore of type %s",sDigitStore->ClassName()));
134           
135   // average arrival time to chambers, for pileup studies
136
137   for ( Int_t iEvent = 0; iEvent < nofEvents; ++iEvent ) 
138   {    
139     // Loop over events.
140     TObjArray tdlist;
141     tdlist.SetOwner(kTRUE);
142     
143     AliDebug(1,Form("iEvent=%d",iEvent));
144     runLoader->GetEvent(iEvent);
145   
146     // for pile up studies
147     float t0=fgkMaxIntTime;  int aa=0;
148     AliHeader* header = runLoader->GetHeader();   
149     AliGenCocktailEventHeader* cocktailHeader =
150       dynamic_cast<AliGenCocktailEventHeader*>(header->GenEventHeader());
151     if (cocktailHeader) {
152       AliGenCocktailEventHeader* genEventHeader = (AliGenCocktailEventHeader*) (header->GenEventHeader());
153       TList* headers = genEventHeader->GetHeaders();
154       TIter nextH(headers);
155       AliGenEventHeader *entry; 
156       while((entry = (AliGenEventHeader*)nextH())) {
157         float t = entry->InteractionTime();     
158         if (TMath::Abs(t)<TMath::Abs(t0)) t0 = t;      
159         aa++;
160       }
161     } else {
162       AliGenEventHeader* evtHeader = 
163         (AliGenEventHeader*)(header->GenEventHeader());
164       float t = evtHeader->InteractionTime();           
165       if (TMath::Abs(t)<TMath::Abs(t0)) t0 = t;           
166       aa++;
167     }
168
169     loader->MakeSDigitsContainer();
170
171     TTree* treeS = loader->TreeS();
172
173     if ( !treeS )
174     {
175       AliFatal("");
176     }
177
178     sDigitStore->Connect(*treeS);
179     
180     TTree* treeH = loader->TreeH();
181
182     AliMUONVHitStore* hitStore = AliMUONVHitStore::Create(*treeH);
183     hitStore->Connect(*treeH);
184     
185     Long64_t nofTracks = treeH->GetEntries();
186     
187     for ( Long64_t iTrack = 0; iTrack < nofTracks; ++iTrack )
188     {
189       // Loop over the tracks of this event.
190       treeH->GetEvent(iTrack);
191
192       AliMUONHit* hit;
193       TIter next(hitStore->CreateIterator());
194       Int_t ihit(0);
195       
196       while ( ( hit = static_cast<AliMUONHit*>(next()) ) )       
197       {
198         Int_t chamberId = hit->Chamber()-1;
199         Float_t age = hit->Age()-t0;
200
201         AliMUONChamber& chamber = muon->Chamber(chamberId);
202         AliMUONResponse* response = chamber.ResponseModel();
203         
204         // This is the heart of this method : the dis-integration
205         TList digits;        
206         if (aa>1){  // if there are pileup events
207           Float_t chamberTime = AliMUONConstants::AverageChamberT(chamberId);
208           Float_t timeDif=age-chamberTime;        
209           if (timeDif>fgkMaxPosTimeDif || timeDif<fgkMaxNegTimeDif) {
210             continue;
211           }
212           if(TMath::Abs(timeDif)>fgkMinTimeDif){
213             response->DisIntegrate(*hit,digits,timeDif);
214           }
215           else{
216             response->DisIntegrate(*hit,digits,0.);
217           }
218         }
219         else{
220           response->DisIntegrate(*hit,digits,0.);
221         }
222         
223         TIter nextd(&digits);
224         AliMUONVDigit* d;
225         while ( ( d = (AliMUONVDigit*)nextd() ) )
226         {
227           // Update some sdigit information that could not be known
228           // by the DisIntegrate method
229           d->SetHit(ihit);
230           d->SetTime(age);
231           d->AddTrack(hit->GetTrack(),d->Charge());
232           tdlist.Add(d);
233         }
234         ++ihit;
235       }
236       hitStore->Clear();
237     } // end of loop on tracks within an event
238     
239     TIter next(&tdlist);
240     AliMUONVDigit* d;
241     
242     while ( ( d = static_cast<AliMUONVDigit*>(next()) ) )
243     {
244       if ( d->Charge() > 0 ) // that check would be better in the disintegrate
245         // method, but to compare with old sdigitizer, it has to be there.
246       {
247         AliMUONVDigit* added = sDigitStore->Add(*d,AliMUONVDigitStore::kMerge);
248         if (!added)
249         {
250           AliError("Could not add digit to digitStore");
251         }
252       }
253     }
254
255     treeS->Fill();
256     
257     loader->WriteSDigits("OVERWRITE");
258     
259     sDigitStore->Clear();
260     
261     loader->UnloadSDigits();
262     
263     delete hitStore;
264     
265   } // loop on events
266   
267   loader->UnloadHits();  
268   
269   delete sDigitStore;
270
271   if(calibrationData) delete calibrationData;
272 }