1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///_________________________________________________________________________
20 /// This class constructs Digits out of Hits
24 // --- Standard library ---
26 // --- ROOT system ---
31 // --- AliRoot header files ---
33 #include "AliACORDE.h"
34 #include "AliACORDEhit.h"
35 #include "AliRunLoader.h"
36 #include "AliLoader.h"
37 #include "AliDigitizationInput.h"
38 #include "AliCDBManager.h"
39 #include "AliCDBStorage.h"
40 #include "AliCDBEntry.h"
41 #include "AliACORDECalibData.h"
42 #include "AliACORDEConstants.h"
44 #include "AliACORDEdigit.h"
45 #include "AliACORDEDigitizer.h"
47 ClassImp(AliACORDEDigitizer)
49 AliACORDEDigitizer::AliACORDEDigitizer()
51 fCalibData(GetCalibData()),
56 // default constructor
59 AliACORDEDigitizer::AliACORDEDigitizer(AliDigitizationInput* digInput)
60 :AliDigitizer(digInput),
61 fCalibData(GetCalibData()),
70 AliACORDEDigitizer::~AliACORDEDigitizer()
82 Bool_t AliACORDEDigitizer::Init()
84 // Initialises the digitizer
86 // Initialises the Digit array
87 fDigits = new TClonesArray ("AliACORDEdigit", 1000);
92 void AliACORDEDigitizer::Digitize(Option_t* /*option*/)
95 // Creates digits from hits
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
114 // 1.- temporal variables
115 Float_t emin = AliACORDEConstants::Instance()->HitEnergyThreshold();
116 Float_t td = AliACORDEConstants::Instance()->MaxHitTimeDifference();
117 Int_t modules[60]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
118 Int_t moduls[60]={0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
120 Float_t PlasticTimes[2][60];
121 Int_t PlasticTracks[2][60];
122 for (Int_t i=0;i<60;i++) {
123 PlasticTimes[0][i]=-1.0;
124 PlasticTimes[1][i]=-1.0;
125 PlasticTracks[0][i]=-1;
126 PlasticTracks[1][i]=-1;
130 AliRunLoader* outRunLoader =
131 AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
133 Error("Exec", "Can not get output Run Loader");
136 AliLoader* outLoader = outRunLoader->GetLoader("ACORDELoader");
138 Error("Exec", "Can not get output ACORDE Loader");
141 outLoader->LoadDigits("update");
142 if (!outLoader->TreeD()) outLoader->MakeTree("D");
143 outLoader->MakeDigitsContainer();
144 TTree* treeD = outLoader->TreeD();
145 Int_t bufsize = 16000;
146 treeD->Branch("ACORDEdigit", &fDigits, bufsize);
148 // 3. loop over inputs
149 for (Int_t iInput = 0; iInput < fDigInput->GetNinputs(); iInput++) {
150 AliRunLoader* runLoader =
151 AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
152 AliLoader* loader = runLoader->GetLoader("ACORDELoader");
154 Error("Exec", "Can not get ACORDE Loader for input %d", iInput);
157 if (!runLoader->GetAliRun()) runLoader->LoadgAlice();
159 AliACORDE* acorde = (AliACORDE*) runLoader->GetAliRun()->GetDetector("ACORDE");
161 Error("Exec", "No ACORDE detector for input %d", iInput);
165 TTree* treeH = loader->TreeH();
167 Error("Exec", "Cannot get TreeH for input %d", iInput);
170 TClonesArray* hits = acorde->Hits();
172 // 3.1 loop over all tracks
173 Int_t nTracks = (Int_t) treeH->GetEntries();
174 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
176 treeH->GetEvent(iTrack);
177 Int_t nHits = hits->GetEntriesFast();
178 // 3.2 loop over hits
179 for (Int_t iHit = 0; iHit < nHits; iHit++) {
180 AliACORDEhit* hit = (AliACORDEhit *)hits->UncheckedAt(iHit);
183 Float_t eloss_mev = hit->Eloss()*1000.0;
184 Int_t module = hit->GetModule();
186 Int_t plastic = hit->GetPlastic();
187 Float_t time_ns = hit->GetTime()*1e9;
188 Float_t eff = TMath::Sqrt(fCalibData->GetEfficiency(module));
189 // if enough energy and efficiency
190 if( eloss_mev > emin && gRandom->Uniform() < eff ) {
191 // if first hit or earlier track
192 if ((PlasticTimes[plastic-1][module-1] == -1.0) ||
193 (PlasticTimes[plastic-1][module-1] > time_ns) ) {
194 PlasticTimes[plastic-1][module-1]= time_ns;
195 PlasticTracks[plastic-1][module-1]= hit->GetTrack();
199 } // end of track loop
200 for(Int_t i=0;i<60;i++){moduls[i]=modules[i];}
202 loader->UnloadHits();
204 } // end of input loop
206 // 4.- loop over temporal arrays to add hits
207 Int_t tracks[3]={-1,-1,-1};
208 for (Int_t i=0; i<60; i++) {
209 // if both modules have a hit
210 // if time diff small enough
211 Float_t diff = TMath::Abs(PlasticTimes[0][i]-PlasticTimes[1][i]);
213 tracks[0] = PlasticTracks[0][i];
214 if (PlasticTracks[0][i] != PlasticTracks[1][i])
215 tracks[1] = PlasticTracks[1][i];
218 // Float_t module_time = (PlasticTimes[0][i] > PlasticTimes[1][i] ?
219 // PlasticTimes[0][i] : PlasticTimes[1][i]);
220 // AddDigit(tracks, mods, module_time);
221 AddDigit(tracks, mods, 0);
226 outLoader->WriteDigits("OVERWRITE");
227 outLoader->UnloadDigits();
231 //____________________________________________________________________________
233 void AliACORDEDigitizer::AddDigit(Int_t* track, Int_t module, Float_t time)
238 TClonesArray &ldigits = *fDigits;
239 new(ldigits[fNdigits++]) AliACORDEdigit(track,module,time);
242 void AliACORDEDigitizer::AddDigit(Int_t* modul,Float_t time)
245 TClonesArray &ldigits = *fDigits;
246 new(ldigits[fNdigits++]) AliACORDEdigit(modul,time);
250 void AliACORDEDigitizer::ResetDigit()
256 if (fDigits) fDigits->Delete();
260 AliACORDECalibData* AliACORDEDigitizer::GetCalibData() const
263 AliCDBManager *man = AliCDBManager::Instance();
265 AliCDBEntry *entry=0;
267 entry = man->Get("ACORDE/Calib/Data");
270 AliWarning("Load of calibration data from default storage failed!");
271 AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
272 Int_t runNumber = man->GetRun();
273 entry = man->GetStorage("local://$ALICE_ROOT/OCDB")
274 ->Get("ACORDE/Calib/Data",runNumber);
278 // Retrieval of data in directory ACORDE/Calib/Data:
281 AliACORDECalibData *calibdata = 0;
283 if (entry) calibdata = (AliACORDECalibData*) entry->GetObject();
284 if (!calibdata) AliError("No calibration data from calibration database !");