]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ACORDE/AliACORDEDigitizer.cxx
Added const members to contain the IDs of FMD files used by the preprocessor and...
[u/mrichter/AliRoot.git] / ACORDE / AliACORDEDigitizer.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 ///_________________________________________________________________________
19 ///
20 /// This class constructs Digits out of Hits
21 ///
22 ///
23
24 // --- Standard library ---
25
26 // --- ROOT system ---
27 #include <TMath.h>
28 #include <TTree.h>
29 #include <TRandom.h>
30
31 // --- AliRoot header files ---
32 #include "AliRun.h"
33 #include "AliACORDE.h"
34 #include "AliACORDEhit.h"
35 #include "AliRunLoader.h"
36 #include "AliLoader.h"
37 #include "AliRunDigitizer.h"
38 #include "AliCDBManager.h"
39 #include "AliCDBStorage.h"
40 #include "AliCDBEntry.h"
41 #include "AliACORDECalibData.h"
42 #include "AliACORDEConstants.h"
43
44 #include "AliACORDEdigit.h"
45 #include "AliACORDEDigitizer.h"
46
47 ClassImp(AliACORDEDigitizer)
48
49 AliACORDEDigitizer::AliACORDEDigitizer()
50   :AliDigitizer(),
51    fCalibData(GetCalibData()),
52    fNdigits(0),
53    fDigits(0)
54   
55 {
56   // default constructor
57 }
58
59 AliACORDEDigitizer::AliACORDEDigitizer(AliRunDigitizer* manager)
60   :AliDigitizer(manager),
61    fCalibData(GetCalibData()),
62    fNdigits(0),
63    fDigits(0)
64
65 {
66   // constructor
67
68 }
69
70 AliACORDEDigitizer::~AliACORDEDigitizer()
71 {
72   // destructor
73   
74   if (fDigits) {
75     fDigits->Delete();
76     delete fDigits;
77     fDigits=0;
78   }
79 }
80
81
82 Bool_t AliACORDEDigitizer::Init()
83 {
84   // Initialises the digitizer
85   
86   // Initialises the Digit array
87   fDigits = new TClonesArray ("AliACORDEdigit", 1000);
88   
89   return kTRUE;
90 }
91
92 void AliACORDEDigitizer::Exec(Option_t* /*option*/)
93 {
94
95   // Creates digits from hits
96
97   // 1.- create and initialize temporal variables
98   // 2.- get loaders to access hots and digits
99   // 3.- loop over all input. 
100   //     3.1 for each input loop over all track
101   //     3.2 for each track loop over all hits
102   //     3.3 for each hit, check 
103   //         if energy above threshold
104   //         if survives efficiency
105   //         then store the plastic, the time and track in temp arrays
106   //         if a hit already survived for this plastic take the hit
107   //         with the lowest time
108   // 4.- loop over temporal array
109   //     if both plastic have a surviving hit and the time
110   //     difference is below the limit, add a new digit
111   // 
112
113
114   // 1.- temporal variables
115   Float_t emin = AliACORDEConstants::Instance()->HitEnergyThreshold();
116   Float_t td = AliACORDEConstants::Instance()->MaxHitTimeDifference();
117   Float_t PlasticTimes[2][60]; 
118   Int_t PlasticTracks[2][60];
119   for (Int_t i=0;i<60;i++) {
120     PlasticTimes[0][i]=-1.0;
121     PlasticTimes[1][i]=-1.0;
122     PlasticTracks[0][i]=-1;
123     PlasticTracks[1][i]=-1;
124   }
125   
126   
127   // 2.- get loaders
128   AliRunLoader* outRunLoader =
129     AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
130   if (!outRunLoader) {
131     Error("Exec", "Can not get output Run Loader");
132     return;}
133   
134   AliLoader* outLoader = outRunLoader->GetLoader("ACORDELoader");
135   if (!outLoader) {
136     Error("Exec", "Can not get output ACORDE Loader");
137     return;}
138   
139   outLoader->LoadDigits("update");
140   if (!outLoader->TreeD()) outLoader->MakeTree("D");
141   outLoader->MakeDigitsContainer();
142   TTree* treeD  = outLoader->TreeD();
143   Int_t bufsize = 16000;
144   treeD->Branch("ACORDEdigit", &fDigits, bufsize);
145   
146   // 3. loop over inputs
147   for (Int_t iInput = 0; iInput < fManager->GetNinputs(); iInput++) {
148     AliRunLoader* runLoader =
149       AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
150     AliLoader* loader = runLoader->GetLoader("ACORDELoader");
151     if (!loader) {
152       Error("Exec", "Can not get ACORDE Loader for input %d", iInput);
153       continue;}
154     
155     if (!runLoader->GetAliRun()) runLoader->LoadgAlice();
156     
157     AliACORDE* acorde = (AliACORDE*) runLoader->GetAliRun()->GetDetector("ACORDE");
158     if (!acorde) {
159       Error("Exec", "No ACORDE detector for input %d", iInput);
160       continue;}
161     
162     loader->LoadHits();
163     TTree* treeH = loader->TreeH();
164     if (!treeH) {
165       Error("Exec", "Cannot get TreeH for input %d", iInput);
166       continue; }
167     
168     TClonesArray* hits = acorde->Hits();
169     
170     // 3.1 loop over all tracks
171     Int_t nTracks = (Int_t) treeH->GetEntries();
172     for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
173       acorde->ResetHits();
174       treeH->GetEvent(iTrack);
175       Int_t nHits = hits->GetEntriesFast();
176       // 3.2 loop over hits
177       for (Int_t iHit = 0; iHit < nHits; iHit++) {
178         AliACORDEhit* hit = (AliACORDEhit *)hits->UncheckedAt(iHit);
179         // 3.3 select hit
180         // get hit info
181         Float_t eloss_mev = hit->Eloss()*1000.0;
182         Int_t module = hit->GetModule();
183         Int_t plastic = hit->GetPlastic();
184         Float_t time_ns = hit->GetTime()*1e9;
185         Float_t eff = TMath::Sqrt(fCalibData->GetEfficiency(module));
186
187         // if enough energy and efficiency
188         if( eloss_mev > emin && gRandom->Uniform() < eff ) {
189           // if first hit or earlier track
190           if ((PlasticTimes[plastic-1][module-1] == -1.0) ||
191               (PlasticTimes[plastic-1][module-1] > time_ns) ) {
192             PlasticTimes[plastic-1][module-1]= time_ns;
193             PlasticTracks[plastic-1][module-1]= hit->GetTrack();
194           } 
195         }
196       } // end of hit   loop
197     } // end of track loop
198     
199     loader->UnloadHits();
200     
201   }  // end of input loop
202
203   // 4.- loop over temporal arrays to add hits
204
205   Int_t tracks[3]={-1,-1,-1};
206   for (Int_t i=0; i<60; i++) {
207     // if both modules have a hit
208     if (PlasticTimes[0][i] == -1) continue;
209     if (PlasticTimes[1][i] == -1) continue;
210     // if time diff small enough
211     Float_t diff = TMath::Abs(PlasticTimes[0][i]-PlasticTimes[1][i]);
212     if (diff < td) {
213       tracks[0] = PlasticTracks[0][i];
214       if (PlasticTracks[0][i] != PlasticTracks[1][i]) 
215         tracks[1] = PlasticTracks[1][i];
216       Int_t module = i+1;
217       Float_t module_time = (PlasticTimes[0][i] > PlasticTimes[1][i] ? 
218                            PlasticTimes[0][i] : PlasticTimes[1][i]);
219       AddDigit(tracks, module, module_time);
220     }
221   }
222   
223   treeD->Fill();
224   outLoader->WriteDigits("OVERWRITE");
225   outLoader->UnloadDigits();
226   ResetDigit();
227 }
228
229 //____________________________________________________________________________
230
231 void AliACORDEDigitizer::AddDigit(Int_t* track, Int_t module, Float_t time)
232 {
233   
234   // Adds Digit
235   
236   TClonesArray &ldigits = *fDigits;
237   new(ldigits[fNdigits++]) AliACORDEdigit(track,module,time);
238 }
239
240
241 void AliACORDEDigitizer::ResetDigit()
242 {
243 //
244 // Clears Digits
245 //
246   fNdigits = 0;
247   if (fDigits) fDigits->Delete();
248 }
249
250
251 AliACORDECalibData* AliACORDEDigitizer::GetCalibData() const
252
253 {
254   AliCDBManager *man = AliCDBManager::Instance();
255
256   AliCDBEntry *entry=0;
257
258   entry = man->Get("ACORDE/Calib/Data");
259
260   if(!entry){
261     AliWarning("Load of calibration data from default storage failed!");
262     AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
263     Int_t runNumber = man->GetRun();
264     entry = man->GetStorage("local://$ALICE_ROOT")
265       ->Get("ACORDE/Calib/Data",runNumber);
266
267   }
268
269   // Retrieval of data in directory ACORDE/Calib/Data:
270
271
272   AliACORDECalibData *calibdata = 0;
273
274   if (entry) calibdata = (AliACORDECalibData*) entry->GetObject();
275   if (!calibdata)  AliError("No calibration data from calibration database !");
276
277
278   return calibdata;
279
280 }
281