]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONSDigitizerV2.cxx
Adding 2 triggers to list of known unsupported triggers
[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 "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   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       float t = evtHeader->InteractionTime();           
155       if (TMath::Abs(t)<TMath::Abs(t0)) t0 = t;           
156       aa++;
157     }
158
159     loader->MakeSDigitsContainer();
160
161     TTree* treeS = loader->TreeS();
162
163     if ( !treeS )
164     {
165       AliFatal("");
166     }
167
168     sDigitStore->Connect(*treeS);
169     
170     TTree* treeH = loader->TreeH();
171
172     AliMUONVHitStore* hitStore = AliMUONVHitStore::Create(*treeH);
173     hitStore->Connect(*treeH);
174     
175     Long64_t nofTracks = treeH->GetEntries();
176     
177     for ( Long64_t iTrack = 0; iTrack < nofTracks; ++iTrack )
178     {
179       // Loop over the tracks of this event.
180       treeH->GetEvent(iTrack);
181
182       AliMUONHit* hit;
183       TIter next(hitStore->CreateIterator());
184       Int_t ihit(0);
185       
186       while ( ( hit = static_cast<AliMUONHit*>(next()) ) )       
187       {
188         Int_t chamberId = hit->Chamber()-1;
189         Float_t age = hit->Age()-t0;
190
191         AliMUONChamber& chamber = muon->Chamber(chamberId);
192         AliMUONResponse* response = chamber.ResponseModel();
193         
194         // This is the heart of this method : the dis-integration
195         TList digits;        
196         if (aa>1){  // if there are pileup events
197           Float_t chamberTime = AliMUONConstants::AverageChamberT(chamberId);
198           Float_t timeDif=age-chamberTime;        
199           if (timeDif>fgkMaxPosTimeDif || timeDif<fgkMaxNegTimeDif) {
200             continue;
201           }
202           if(TMath::Abs(timeDif)>fgkMinTimeDif){
203             response->DisIntegrate(*hit,digits,timeDif);
204           }
205           else{
206             response->DisIntegrate(*hit,digits,0.);
207           }
208         }
209         else{
210           response->DisIntegrate(*hit,digits,0.);
211         }
212         
213         TIter nextd(&digits);
214         AliMUONVDigit* d;
215         while ( ( d = (AliMUONVDigit*)nextd() ) )
216         {
217           // Update some sdigit information that could not be known
218           // by the DisIntegrate method
219           d->SetHit(ihit);
220           d->SetTime(age);
221           d->AddTrack(hit->GetTrack(),d->Charge());
222           tdlist.Add(d);
223         }
224         ++ihit;
225       }
226       hitStore->Clear();
227     } // end of loop on tracks within an event
228     
229     TIter next(&tdlist);
230     AliMUONVDigit* d;
231     
232     while ( ( d = static_cast<AliMUONVDigit*>(next()) ) )
233     {
234       d->ChargeInFC(kTRUE);
235       
236       AliMUONVDigit* added = sDigitStore->Add(*d,AliMUONVDigitStore::kMerge);
237       if (!added)
238       {
239         AliError("Could not add digit to digitStore");
240       }
241     }
242   
243     treeS->Fill();
244     
245     loader->WriteSDigits("OVERWRITE");
246     
247     sDigitStore->Clear();
248     
249     loader->UnloadSDigits();
250     
251     delete hitStore;
252     
253   } // loop on events
254   
255   loader->UnloadHits();  
256   
257   delete sDigitStore;
258
259 }