]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONPedestalEventGenerator.cxx
Adding comments only
[u/mrichter/AliRoot.git] / MUON / AliMUONPedestalEventGenerator.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 "AliMUONPedestalEventGenerator.h"
19
20 #include "AliDAQ.h"
21 #include "AliHeader.h"
22 #include "AliLog.h"
23 #include "AliMUONCalibrationData.h"
24 #include "AliMUONData.h"
25 #include "AliMUONDigit.h"
26 #include "AliMUONRawWriter.h"
27 #include "AliMUONVCalibParam.h"
28 #include "AliMpIntPair.h"
29 #include "AliMpManuList.h"
30 #include "AliRunLoader.h"
31
32 #include "TClonesArray.h"
33 #include "TRandom.h"
34 #include "TStopwatch.h"
35 #include "TSystem.h"
36 #include "TMath.h"
37
38 ///
39 /// \class AliMUONPedestalEventGenerator
40 ///
41 /// Generate simulated pedestal events for MUON TRK, to be able to e.g. test
42 /// online calibration routines.
43 ///
44 /// The pedestals themselves are taken from the CDB. What we get from the CDB
45 /// is, per channel, the mean and the sigma of the pedestal. We then use
46 /// those informations to randomly get the pedestals for each channel, for
47 /// each event (picking in a gaus(mean,sigma)).
48 ///
49 /// Output can be just digits, or digits + raw (ddl), or digits + raw (ddl)
50 /// + raw (date files, one per LDC), depending of ctor and MakeDDL() method.
51 ///
52 /// \author L. Aphecetche
53 ///
54
55 /// \cond CLASSIMP
56 ClassImp(AliMUONPedestalEventGenerator)
57 /// \endcond
58
59 Int_t AliMUONPedestalEventGenerator::fgCounter(0);
60
61 //_____________________________________________________________________________
62 AliMUONPedestalEventGenerator::AliMUONPedestalEventGenerator(Int_t runNumber,
63                                                              Int_t nevents,
64                                                              const char* filename)
65 : TTask("AliMUONPedestalEventGenerator","Generate fake pedestal events"), 
66 fManuList(AliMpManuList::ManuList()),
67 fCalibrationData(new AliMUONCalibrationData(runNumber)),
68 fDateFileName(filename),
69 fGAliceFileName("galice.root"),
70 fMakeDDL(kTRUE)
71 {
72   /// Will generate pedestals according to (mean,sigma)s found in CDB
73   /// for run runNumber.
74   /// Will generate nevents events
75   /// If filename is != "", it will be the basename of the output LDC files
76   ///
77   if (!gSystem->IsAbsoluteFileName(fGAliceFileName)) 
78   {
79     char* absFileName = gSystem->ConcatFileName(gSystem->WorkingDirectory(),
80                                                 fGAliceFileName);
81     fGAliceFileName = absFileName;
82     delete[] absFileName;
83   }
84   
85   AliRunLoader* runLoader = LoadRun("recreate");
86   
87   runLoader->SetNumberOfEventsPerFile(nevents);
88   
89   if (!runLoader)
90   {
91     AliError("Could not create RunLoader");
92     return;
93   }
94   
95   AliConfig::Instance()
96     ->CreateDetectorFolders(runLoader->GetEventFolder(), 
97                             "MUON", "MUON");
98   
99   AliLoader* loader = new AliLoader("MUON",runLoader->GetEventFolder());
100   
101         runLoader->AddLoader(loader);
102   
103   // Initialize event headers.
104   runLoader->MakeTree("E");
105   for ( Int_t iEvent = 0; iEvent < nevents; ++iEvent )
106   {
107     runLoader->SetEventNumber(iEvent);
108     runLoader->GetHeader()->Reset(runNumber,iEvent);
109     runLoader->TreeE()->Fill();
110   }
111   runLoader->WriteHeader("OVERWRITE");
112   runLoader->CdGAFile();
113   runLoader->Write(0, TObject::kOverwrite);  
114
115   delete runLoader;
116 }
117
118
119 //_____________________________________________________________________________
120 AliMUONPedestalEventGenerator::~AliMUONPedestalEventGenerator()
121 {
122   /// dtor
123   delete fManuList;
124   delete fCalibrationData;
125   AliInfo(Form("make a digit counter %d",fgCounter));
126 }
127
128 //_____________________________________________________________________________
129 Bool_t 
130 AliMUONPedestalEventGenerator::ConvertRawFilesToDate()
131 {
132   /// convert raw data DDL files to DATE files with the program "dateStream".
133   /// we make one file per LDC
134   
135   const Int_t kIDet = AliDAQ::DetectorID("MUONTRK");
136   
137   const Int_t kNLDCs = TMath::CeilNint(AliDAQ::NumberOfLdcs(kIDet));
138   
139   char* path = gSystem->Which(gSystem->Getenv("PATH"), "dateStream");
140   if (!path) 
141   {
142     AliError("the program dateStream was not found");
143     return kFALSE;
144   }
145   
146   delete[] path;
147   
148   AliRunLoader* runLoader = LoadRun("read");
149   if (!runLoader) return kFALSE;
150   
151   AliInfo(Form("converting raw data DDL files to DATE files %s", fDateFileName.Data()));
152   FILE** pipe = new FILE*[kNLDCs];
153   
154   for ( Int_t iFile = 0; iFile < kNLDCs; ++iFile)
155   {
156     char command[256];
157     // Note the option -s. It is used in order to avoid
158     // the generation of SOR/EOR events.
159     sprintf(command, "dateStream -s -D -o %s.LDC%d -# %d -C", 
160             fDateFileName.Data(), iFile, runLoader->GetNumberOfEvents());
161     pipe[iFile] = gSystem->OpenPipe(command, "w");
162   }
163   
164   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); ++iEvent) 
165   {
166     Float_t ldc = 0;
167     Int_t prevLDC = -1;
168
169     for (Int_t iDDL = 0; iDDL < AliDAQ::NumberOfDdls(kIDet); ++iDDL) 
170     {        
171       Int_t ddlID = AliDAQ::DdlID(kIDet,iDDL);
172       Int_t ldcID = Int_t(ldc + 0.0001);
173       ldc += AliDAQ::NumberOfLdcs(kIDet) / AliDAQ::NumberOfDdls(kIDet);
174       
175       char rawFileName[256];
176       sprintf(rawFileName, "raw%d/%s", 
177               iEvent, AliDAQ::DdlFileName(kIDet,iDDL));
178       
179       // check existence and size of raw data file
180       FILE* file = fopen(rawFileName, "rb");
181       if (!file) continue;
182       fseek(file, 0, SEEK_END);
183       unsigned long size = ftell(file);
184       fclose(file);
185       if (!size) continue;
186       
187       if (ldcID != prevLDC) {
188         fprintf(pipe[ldcID], " LDC Id %d\n", ldcID);
189         prevLDC = ldcID;
190       }
191       fprintf(pipe[ldcID], "  Equipment Id %d Payload %s\n", ddlID, rawFileName);
192     }
193   }
194   
195   Int_t result(0);
196   
197   for ( Int_t iFile = 0; iFile < kNLDCs; ++iFile)
198   {
199     result += gSystem->ClosePipe(pipe[iFile]);
200   }
201   
202   delete [] pipe;
203   delete runLoader;
204   return (result == 0);
205 }
206
207 //_____________________________________________________________________________
208 void
209 AliMUONPedestalEventGenerator::Exec(Option_t*)
210 {  
211   /// Main steering method
212   
213   TStopwatch timer;
214   timer.Start(kTRUE);
215   
216   TStopwatch writeTimer;
217   writeTimer.Start(kTRUE); writeTimer.Stop();
218   TStopwatch fillTimer;
219   fillTimer.Start(kTRUE); fillTimer.Stop();
220   TStopwatch resetTimer;
221   resetTimer.Start(kTRUE); resetTimer.Stop();
222   TStopwatch getEventTimer;
223   getEventTimer.Start(kTRUE); getEventTimer.Stop();
224   TStopwatch branchTimer;
225   branchTimer.Start(kTRUE); branchTimer.Stop();
226   TStopwatch generateDigitsTimer;
227   generateDigitsTimer.Start(kTRUE); generateDigitsTimer.Stop();
228   TStopwatch digits2RawTimer;
229   digits2RawTimer.Start(kTRUE); digits2RawTimer.Stop();
230   TStopwatch convertRawFilesToDateTimer;
231   convertRawFilesToDateTimer.Start(kTRUE); convertRawFilesToDateTimer.Stop();
232   TStopwatch getTimer;
233
234   getTimer.Start(kTRUE); 
235   AliMUONData* data = GetDataAccess("update");
236   AliRunLoader* runLoader = data->GetLoader()->GetRunLoader();  
237   getTimer.Stop();
238   
239   for ( Int_t i = 0; i < runLoader->GetNumberOfEvents(); ++i )
240   {
241     getEventTimer.Start(kFALSE);
242     runLoader->GetEvent(i);
243     getEventTimer.Stop();
244     
245     branchTimer.Start(kFALSE);
246     if ( data->TreeD() == 0x0 )
247     {
248       AliDebug(1,"Calling MakeDigitsContainer");
249       data->GetLoader()->MakeDigitsContainer();
250     }
251     data->MakeBranch("D,GLT");
252     data->SetTreeAddress("D,GLT");
253     branchTimer.Stop();
254     
255     generateDigitsTimer.Start(kFALSE);
256     GenerateDigits(data);
257     generateDigitsTimer.Stop();
258     
259     fillTimer.Start(kTRUE);
260     // Fill the output treeD
261     data->Fill("D,GLT");
262     fillTimer.Stop();
263     
264     // Write to the output tree(D).
265     // Please note that as GlobalTrigger, LocalTrigger and Digits are in the same
266     // tree (=TreeD) in different branches, this WriteDigits in fact writes all of 
267     // the 3 branches.
268     
269     writeTimer.Start(kFALSE);
270     data->GetLoader()->WriteDigits("OVERWRITE");
271     writeTimer.Stop();
272     
273     // Finally, we clean up after ourselves.
274     resetTimer.Start(kTRUE);
275     data->ResetDigits();
276     data->ResetTrigger();
277     data->GetLoader()->UnloadDigits();
278     resetTimer.Stop();
279   }
280   
281   TStopwatch endTimer;
282   endTimer.Start(kTRUE);
283   runLoader->WriteRunLoader("OVERWRITE");
284   delete data;
285   delete runLoader;
286   endTimer.Stop();
287   
288   if ( fMakeDDL ) 
289   {
290     // Now convert the digits.root file(s) to DDL files
291     digits2RawTimer.Start(kFALSE);
292     Digits2Raw();
293     digits2RawTimer.Stop();
294   }
295   
296   // Finally, if instructed to do so, convert DDL files to DATE file(s)
297   if ( fMakeDDL && fDateFileName.Length() > 0 ) 
298   {
299     convertRawFilesToDateTimer.Start(kFALSE);
300     AliDebug(1,"Converting to DATE file(s)");
301     
302     Bool_t dateOutput = ConvertRawFilesToDate();
303     convertRawFilesToDateTimer.Stop();
304     if (!dateOutput) 
305     {
306       AliError("DATE output failed. Aborting.");
307       return;
308     }    
309   }
310   
311   AliInfo(Form("Execution time Exec : R:%.2fs C:%.2fs",
312                timer.RealTime(),timer.CpuTime()));
313
314   AliInfo(Form("  Execution time branch : R:%.2fs C:%.2fs",
315                branchTimer.RealTime(),branchTimer.CpuTime()));
316   AliInfo(Form("  Execution time getEvent : R:%.2fs C:%.2fs",
317                getEventTimer.RealTime(),getEventTimer.CpuTime()));
318   AliInfo(Form("  Execution time fill digits : R:%.2fs C:%.2fs",
319                fillTimer.RealTime(),fillTimer.CpuTime()));
320   AliInfo(Form("  Execution time write digits : R:%.2fs C:%.2fs",
321                writeTimer.RealTime(),writeTimer.CpuTime()));
322   AliInfo(Form("  Execution time reset digits : R:%.2fs C:%.2fs",
323                resetTimer.RealTime(),resetTimer.CpuTime()));
324   AliInfo(Form("  Execution time for GenerateDigits : R:%.2fs C:%.2fs",
325                generateDigitsTimer.RealTime(),generateDigitsTimer.CpuTime()));
326   AliInfo(Form("  Execution time for Digits2Raw : R:%.2fs C:%.2fs",
327                digits2RawTimer.RealTime(),digits2RawTimer.CpuTime()));
328   AliInfo(Form("  Execution time for ConvertRawFilesToDate : R:%.2fs C:%.2fs",
329                convertRawFilesToDateTimer.RealTime(),convertRawFilesToDateTimer.CpuTime()));
330   AliInfo(Form("  Execution time for get : R:%.2fs C:%.2fs",
331                getTimer.RealTime(),getTimer.CpuTime()));
332   AliInfo(Form("  Execution time for end : R:%.2fs C:%.2fs",
333                endTimer.RealTime(),endTimer.CpuTime()));
334 }
335
336 //_____________________________________________________________________________
337 void
338 AliMUONPedestalEventGenerator::Digits2Raw()
339 {
340   /// Converts digits (from MUON.Digits.root file) to Raw DDL ascii files.
341   
342   AliMUONData* data = GetDataAccess("read");
343   AliRunLoader* runLoader = data->GetLoader()->GetRunLoader();
344   AliDebug(1,Form("runLoader=%p",runLoader));
345   
346   AliMUONRawWriter rawWriter(data);
347   
348   // Generate RAW data from the digits
349   // Be carefull to create&change to the correct directory first...
350
351   TString baseDir = gSystem->WorkingDirectory();
352
353   for ( Int_t i = 0; i < runLoader->GetNumberOfEvents(); ++i )
354   {
355     runLoader->GetEvent(i);
356     
357     AliDebug(1,Form("processing event %d", i));
358     
359     char dirName[256];
360     sprintf(dirName, "raw%d", i);
361     gSystem->MakeDirectory(dirName);
362     if (!gSystem->ChangeDirectory(dirName)) 
363     {
364       AliError(Form("couldn't change to directory %s", dirName));
365       return;
366     }
367
368     rawWriter.Digits2Raw();
369     
370     gSystem->ChangeDirectory(baseDir);
371   }
372   
373   delete data;
374   delete runLoader;
375   
376   AliDebug(1,Form("DDL files written successfully"));    
377 }
378
379 //_____________________________________________________________________________
380 void
381 AliMUONPedestalEventGenerator::GenerateDigits(AliMUONData* data)
382 {  
383   /// Generate digits (where ADC is set to pedestal value) for all MUON TRK
384   /// and for 1 event.
385   
386   TIter next(fManuList);
387   AliMpIntPair* p;
388   Int_t ngenerated(0);
389   
390   while ( ( p = static_cast<AliMpIntPair*>(next())) )
391   {
392     Int_t detElemId = p->GetFirst();
393
394     Int_t chamber = detElemId/100-1;
395         
396     TClonesArray* pdigits = data->Digits(chamber);
397     if (!pdigits)
398     {
399       AliError(Form("No digits for chamber %d",1));
400       continue;
401     }
402     TClonesArray& digits = *pdigits;
403     
404     Int_t manuId = p->GetSecond();
405     
406     AliMUONVCalibParam* pedestals = fCalibrationData->Pedestals(detElemId,manuId);
407     if (!pedestals)
408     {
409       AliError(Form("Could not find pedestals for (DE,MANU)=(%d,%d)",detElemId,
410                     manuId));
411       return;
412     }
413     for ( Int_t manuChannel = 0; manuChannel < pedestals->Size(); ++manuChannel )
414     {
415       Float_t mean = pedestals->ValueAsFloat(manuChannel,0);
416       if (mean == AliMUONVCalibParam::InvalidFloatValue())
417       {
418         // This is a poor's man way of knowing if that channel really exists.
419         // Better and safer way (but much slower too) would be to check pad existence
420         // using AliMpVSegmentation::PadByLocation(AliMpIntPair(manuId,manuChannel))
421         continue;
422       }
423       else
424       {
425         Float_t sigma = pedestals->ValueAsFloat(manuChannel,1);
426         AliMUONDigit d;
427         Float_t ped = gRandom->Gaus(mean,sigma);
428         Int_t pedADC = TMath::FloorNint(ped);
429         d.SetElectronics(manuId, manuChannel);
430         d.SetADC(pedADC);
431         d.SetSignal(ped);
432         d.SetDetElemId(detElemId); 
433         // we do not set the remaining parts of the digit, as in principle
434         // this is all we need : manuId, manuChannel and ADC, as far as
435         // real data is concerned.
436         new (digits[digits.GetLast()+1]) AliMUONDigit(d);
437         ++fgCounter;
438         ++ngenerated;
439       }
440     }
441   }
442   AliDebug(1,Form("ngenerated=%d",ngenerated));
443 }
444
445 //_____________________________________________________________________________
446 AliMUONData*
447 AliMUONPedestalEventGenerator::GetDataAccess(const char* mode)
448 {
449   /// Get the pointer to AliMUONData object
450   AliRunLoader* runLoader = LoadRun(mode);
451   if (!runLoader)
452   {
453     AliError("Could not get RunLoader");
454     return 0x0;
455   }
456   AliLoader* loader = static_cast<AliLoader*>(runLoader->GetLoader("MUONLoader"));
457   if (!loader)
458   {
459     AliError("Could not get MuonLoader");
460     return 0x0;
461   }
462   
463   return new AliMUONData(loader,"MUON","MUONData");
464 }
465
466 //_____________________________________________________________________________
467 AliRunLoader*
468 AliMUONPedestalEventGenerator::LoadRun(const char* mode)
469 {
470   /// Get access to AliRunLoader object
471   while (AliRunLoader::GetRunLoader()) 
472   {
473     AliDebug(1,Form("Deleting AliRunLoader %p",AliRunLoader::GetRunLoader()));
474     delete AliRunLoader::GetRunLoader();
475   }
476   
477   AliRunLoader* runLoader = 
478     AliRunLoader::Open(fGAliceFileName,AliConfig::GetDefaultEventFolderName(), 
479                        mode);
480   
481   AliDebug(1,Form("AliRunLoader(%s)=%p",mode,runLoader));
482            
483   if (!runLoader) 
484   {
485     AliError("No run loader found in file galice.root");
486   }
487   return runLoader;
488 }