]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFDigitizer.cxx
Fix for DA's
[u/mrichter/AliRoot.git] / TOF / AliTOFDigitizer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2000, 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 //_________________________________________________________________________//
17 //                                                                         //
18 // This is a TTask that makes TOF-Digits out of TOF-SDigits.               //
19 // The simulation of the detector is performed at sdigits level:           //
20 // during digitization the unique task is the sum of all sdigits in the    //
21 // same pad.                                                               //
22 // Digits are written to TreeD in branch "TOF".                            //
23 //                                                                         //
24 // -- Author :  F. Pierella (Bologna University) pierella@bo.infn.it       //
25 //                                                                         //
26 //_________________________________________________________________________//
27
28 //#include "Riostream.h"
29
30 //#include "TFile.h"
31 #include "TMath.h"
32 #include "TH1F.h"
33 #include "TTree.h"
34 #include "TRandom.h"
35 #include "TObjArray.h"
36
37 #include "AliLoader.h"
38 #include "AliLog.h"
39 #include "AliRunDigitizer.h"
40 #include "AliRunLoader.h"
41 #include "AliRun.h"
42
43 #include "AliTOFcalib.h"
44 //#include "AliTOFChannelOnlineArray.h"
45 //#include "AliTOFChannelOnlineStatusArray.h"
46 //#include "AliTOFChannelOffline.h"
47 #include "AliTOFDigitizer.h"
48 #include "AliTOFdigit.h"
49 #include "AliTOFHitMap.h"
50 #include "AliTOFGeometry.h"
51 #include "AliTOFSDigit.h"
52 #include "AliTOF.h"
53
54 extern TDirectory *gDirectory;
55 //extern TFile *gFile;
56 extern TRandom *gRandom;
57
58 extern AliRun *gAlice;
59
60
61 ClassImp(AliTOFDigitizer)
62
63 //___________________________________________
64   AliTOFDigitizer::AliTOFDigitizer()  :
65     AliDigitizer(),
66     fDigits(new TClonesArray("AliTOFdigit",4000)),
67     fSDigitsArray(new TClonesArray("AliTOFSDigit",1000)),
68   fhitMap(0x0),
69   fCalib(new AliTOFcalib())
70 {
71   // Default ctor - don't use it
72   InitDecalibration();
73 }
74
75 //___________________________________________
76 AliTOFDigitizer::AliTOFDigitizer(AliRunDigitizer* manager): 
77   AliDigitizer(manager), 
78   fDigits(new TClonesArray("AliTOFdigit",4000)),
79   fSDigitsArray(new TClonesArray("AliTOFSDigit",1000)),
80   fhitMap(0x0),
81   fCalib(new AliTOFcalib())
82 {
83   //ctor with RunDigitizer
84   InitDecalibration();
85 }
86
87 //------------------------------------------------------------------------
88 AliTOFDigitizer::AliTOFDigitizer(const AliTOFDigitizer &source):
89   AliDigitizer(source),
90   fDigits(source.fDigits),
91   fSDigitsArray(source.fSDigitsArray),
92   fhitMap(source.fhitMap),
93   fCalib(source.fCalib)
94 {
95   // copy constructor
96 }
97
98 //------------------------------------------------------------------------
99   AliTOFDigitizer& AliTOFDigitizer::operator=(const AliTOFDigitizer &source)
100 {
101   // ass. op.
102   
103   if (this == &source)
104     return *this;
105
106   AliDigitizer::operator=(source);
107   fDigits=source.fDigits;
108   fSDigitsArray=source.fSDigitsArray;
109   fhitMap=source.fhitMap;
110   fCalib=source.fCalib;
111   return *this;
112
113 }
114
115 //------------------------------------------------------------------------
116 AliTOFDigitizer::~AliTOFDigitizer()
117 {
118   // Destructor
119   delete fCalib;
120   if (fDigits){
121     fDigits->Delete();
122     delete fDigits;
123     fDigits=0x0;
124   }
125   if (fSDigitsArray){
126     fSDigitsArray->Delete();
127     delete fSDigitsArray;
128     fSDigitsArray=0x0;
129   }
130 }
131
132 //---------------------------------------------------------------------
133
134 void AliTOFDigitizer::Exec(Option_t* /*option*/)
135 {
136   //
137   // Perform digitization and merging.
138   // The algorithm is the following:
139   // - a hitmap is created to check if a pad is already activated;
140   // - an sdigits container is created to collect all sdigits from
141   //   different files;
142   // - sdigits are summed using the hitmap;
143   // - the sdigits container is used to create the array of AliTOFdigit.
144   //
145
146   AliDebug(1, "");
147
148
149   // get the ptr to TOF detector
150   AliTOF * tof = (AliTOF *) gAlice->GetDetector("TOF") ;
151
152   //Make branches
153
154   const Int_t kSize = 20;
155   char branchname[kSize];
156   snprintf(branchname,kSize,"%s", tof->GetName ());
157  
158   AliRunLoader* outrl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
159   if (outrl == 0x0)
160    {
161      AliError("Can not find Run Loader in output folder.");
162      return;
163    }
164    
165   AliLoader* outgime = outrl->GetLoader("TOFLoader");
166   if (outgime == 0x0)
167    {
168      AliError("Can not get TOF Loader from Output Run Loader.");
169      return;
170    }
171   
172   TTree* treeD = outgime->TreeD();
173   if (treeD == 0x0)
174    {
175      outgime->MakeTree("D");
176      treeD = outgime->TreeD();
177    }
178   //Make branch for digits (to be created in Init())
179   tof->MakeBranchInTree(treeD,branchname,&fDigits,4000);
180
181   // container for all summed sdigits (to be created in Init())
182   //fSDigitsArray=new TClonesArray("AliTOFSDigit",1000);
183   
184   // create hit map (to be created in Init())
185   fhitMap = new AliTOFHitMap(fSDigitsArray);
186   
187   // Loop over files to digitize
188
189   for (Int_t inputFile=0; inputFile<fManager->GetNinputs();
190        inputFile++) {
191     ReadSDigit(inputFile);
192    }
193
194   // create digits
195   CreateDigits();
196
197   // free used memory for Hit Map in current event
198   delete fhitMap;
199   fSDigitsArray->Clear();
200
201   treeD->Fill();
202
203   AliDebug(2,"----------------------------------------");
204   AliDebug(1,Form("%d digits have been created", fDigits->GetEntriesFast()));
205   AliDebug(2,"----------------------------------------");
206
207   outgime->WriteDigits("OVERWRITE");
208   outgime->UnloadDigits();
209   fDigits->Clear();
210
211 }
212
213 //---------------------------------------------------------------------
214
215 void AliTOFDigitizer::CreateDigits()
216 {
217   // loop on sdigits container to fill the AliTOFdigit TClonesArray
218   // start digitizing all the collected sdigits 
219
220   Int_t ndump=0; // dump the first ndump created digits for each event
221
222   // get the total number of collected sdigits
223   Int_t ndig = fSDigitsArray->GetEntriesFast();
224
225   Int_t  vol[5]={-1,-1,-1,-1,-1};  // location for a digit
226   Int_t  digit[4] = {0,0,0,0};     // TOF digit variables
227   Int_t tracknum[AliTOFSDigit::kMAXDIGITS]; // contributing tracks for the current slot
228   for (Int_t aa=0; aa<AliTOFSDigit::kMAXDIGITS; aa++) tracknum[aa] = -1;
229
230   for (Int_t k = 0; k < ndig; k++) {
231     
232     for (Int_t i=0; i<5; i++) vol[i] = -1;
233     
234     // Get the information for this digit
235     AliTOFSDigit *tofsdigit = (AliTOFSDigit *) fSDigitsArray->UncheckedAt(k);
236     
237     Int_t nslot=tofsdigit->GetNDigits(); // get the number of slots
238     // for current sdigit
239     
240     // TOF sdigit volumes (always the same for all slots)
241     Int_t sector    = tofsdigit->GetSector(); // range [0-17]
242     Int_t plate     = tofsdigit->GetPlate();  // range [0- 4]
243     Int_t strip     = tofsdigit->GetStrip();  // range [0-14/18/19]
244     Int_t padz      = tofsdigit->GetPadz();   // range [0- 1]
245     Int_t padx      = tofsdigit->GetPadx();   // range [0-47]
246     
247     vol[0] = sector;
248     vol[1] = plate;
249     vol[2] = strip;
250     vol[3] = padx;
251     vol[4] = padz;
252     
253     //--------------------- QA section ----------------------
254     // in the while, I perform QA
255     Bool_t isSDigitBad = (sector<0 || sector>17 || plate<0 || plate >4 || padz<0 || padz>1 || padx<0 || padx>47);
256     
257     if (isSDigitBad)
258       AliFatal(Form("strange sdigit found   %2d  %1d  %2d  %1d %2d", sector, plate, strip, padz, padx));
259     //-------------------------------------------------------
260     
261     //------------------- Dump section ----------------------
262     if (k<ndump) {
263       AliInfo(Form("%2d-th digit: Sector %2d | Plate %1d | Strip %2d | PadZ %1d | PadX %2d ", k, sector, plate, strip, padz, padx));
264       AliInfo("----------------------------------------------------");
265     }
266     // ------------------------------------------------------
267     
268     // start loop on number of slots for current sdigit
269     for (Int_t islot = 0; islot < nslot; islot++) {
270       for (Int_t aa=0; aa<4; aa++) digit[aa] = 0; // TOF digit variables
271       for (Int_t aa=0; aa<AliTOFSDigit::kMAXDIGITS; aa++) tracknum[aa] = -1;
272       
273       Int_t tdc=tofsdigit->GetTdc(islot); digit[0]=tdc;
274       Int_t adc=tofsdigit->GetAdc(islot); digit[1]=adc;
275
276       //if (tdc>=8192) continue;//AdC
277
278       tracknum[0]=tofsdigit->GetTrack(islot,0);
279       tracknum[1]=tofsdigit->GetTrack(islot,1);
280       tracknum[2]=tofsdigit->GetTrack(islot,2);
281       
282       // new with placement must be used
283       // adding a TOF digit for each slot
284       TClonesArray &aDigits = *fDigits;
285       Int_t last=fDigits->GetEntriesFast();
286       new (aDigits[last]) AliTOFdigit(tracknum, vol, digit);
287
288     }
289     
290   } // end loop on sdigits - end digitizing all collected sdigits
291
292   //Insert Decalibration 
293   AliDebug(2,"in digitizer, create digits");
294   DecalibrateTOFSignal();
295 }
296
297 //---------------------------------------------------------------------
298
299 void AliTOFDigitizer::ReadSDigit(Int_t inputFile )
300 {
301   // Read sdigits for current event and inputFile; 
302   // store them into the sdigits container
303   // and update the hit map
304   // SDigits from different files are assumed to
305   // be created with the same simulation parameters.
306   
307   // creating the TClonesArray to store the digits
308   static TClonesArray sdigitsClonesArray("AliTOFSDigit",  1000); 
309   sdigitsClonesArray.Clear();
310
311   // get the treeS from manager
312   AliRunLoader* rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
313   if (rl == 0x0)
314    {
315      AliError(Form("Can not find Run Loader in input %d folder.",inputFile));
316      return;
317    }
318
319   AliLoader* gime = rl->GetLoader("TOFLoader");
320   if (gime == 0x0)
321    {
322      AliError(Form("Can not get TOF Loader from Input %d Run Loader.",inputFile));
323      return;
324    }
325
326   TTree* currentTreeS=gime->TreeS();
327   if (currentTreeS == 0x0)
328    {
329      Int_t retval = gime->LoadSDigits();
330      if (retval) 
331       {
332          AliError(Form("Error occured while loading S. Digits for Input %d",inputFile));
333          return;
334       }
335      currentTreeS=gime->TreeS();
336      if (currentTreeS == 0x0)
337       {
338          AliError(Form("Can not get S. Digits Tree for Input %d",inputFile));
339          return;
340       }
341    } 
342   // get the branch TOF inside the treeS
343   TClonesArray * sdigitsDummyContainer=&sdigitsClonesArray;
344   // check if the branch exist
345   TBranch* tofBranch=currentTreeS->GetBranch("TOF");
346
347   if(!tofBranch){
348     AliFatal(Form("TOF branch not found for input %d",inputFile));
349     return;
350   }
351   
352   tofBranch->SetAddress(&sdigitsDummyContainer);           
353   
354   Int_t nEntries = (Int_t)tofBranch->GetEntries();                                
355
356   // Loop through all entries in the tree
357   Int_t nbytes = 0;
358   
359   Int_t  vol[5]; // location for a sdigit
360   for (Int_t i=0; i<5; i++) vol[i] = -1;
361
362   for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
363     
364     // Import the tree
365     nbytes += tofBranch->GetEvent(iEntry);
366     
367     // Get the number of sdigits
368     Int_t ndig = sdigitsDummyContainer->GetEntriesFast();
369     
370     for (Int_t k=0; k<ndig; k++) {
371       AliTOFSDigit *tofSdigit= (AliTOFSDigit*) sdigitsDummyContainer->UncheckedAt(k);
372       
373       for (Int_t i=0; i<5; i++) vol[i] = -1;
374
375       // check the sdigit volume
376       vol[0] = tofSdigit->GetSector();
377       vol[1] = tofSdigit->GetPlate();
378       vol[2] = tofSdigit->GetStrip();
379       vol[3] = tofSdigit->GetPadx();
380       vol[4] = tofSdigit->GetPadz();
381       
382       if (fhitMap->TestHit(vol) != kEmpty) {
383         AliTOFSDigit *sdig = static_cast<AliTOFSDigit*>(fhitMap->GetHit(vol));
384         sdig->Update(tofSdigit);
385
386       } else {
387
388         CollectSDigit(tofSdigit); // collect the current sdigit
389         fhitMap->SetHit(vol);     // update the hitmap for location vol
390
391       } // if (hitMap->TestHit(vol) != kEmpty)
392       
393     } // for (Int_t k=0; k<ndig; k++)
394
395   } // end loop on entries
396
397 }
398
399
400 //_____________________________________________________________________________
401 void AliTOFDigitizer::CollectSDigit(const AliTOFSDigit * const sdigit)
402 {
403   //
404   // Add a TOF sdigit in container
405   // new with placement must be used
406   TClonesArray &aSDigitsArray = *fSDigitsArray;
407   Int_t last=fSDigitsArray->GetEntriesFast();
408   // make a copy of the current sdigit and
409   // put it into tmp array
410   new (aSDigitsArray[last]) AliTOFSDigit(*sdigit);
411 }
412
413 //_____________________________________________________________________________
414 void AliTOFDigitizer::InitDecalibration() const {
415   //
416   // Initialize TOF digits decalibration
417   //
418
419   fCalib->Init();
420   /*
421   fCalib->CreateCalArrays();
422   fCalib->ReadSimHistoFromCDB("TOF/Calib", -1); // use AliCDBManager's number
423   fCalib->ReadParOfflineFromCDB("TOF/Calib", -1); // use AliCDBManager's number
424   */
425 }
426 //---------------------------------------------------------------------
427 void AliTOFDigitizer::DecalibrateTOFSignal() {
428   //
429   // Decalibrate TOF signals according to OCDB parameters
430   //
431
432   Double_t time=0., tot=0., corr=0.;
433   Int_t deltaBC=0, l0l1=0, tdcBin=0;
434   Int_t index = -1;
435   Int_t detId[5] ={-1,-1,-1,-1,-1};
436   UInt_t timestamp=0;
437
438   Int_t ndigits = fDigits->GetEntriesFast();
439   // Loop on TOF Digits
440   for (Int_t i=0;i<ndigits;i++){
441     AliTOFdigit * dig = (AliTOFdigit*)fDigits->At(i);
442     detId[0] = dig->GetSector();
443     detId[1] = dig->GetPlate();
444     detId[2] = dig->GetStrip();
445     detId[3] = dig->GetPadz();
446     detId[4] = dig->GetPadx();
447     dig->SetTdcND(dig->GetTdc()); // save the non decalibrated time
448
449     index = AliTOFGeometry::GetIndex(detId); // The channel index    
450
451     // Read Calibration parameters from the CDB
452     // get digit info
453     time = dig->GetTdc() * AliTOFGeometry::TdcBinWidth(); /* ps */
454     tot = dig->GetToT() * AliTOFGeometry::ToTBinWidth() * 1.e-3; /* ns */
455     deltaBC = 0;//dig->GetDeltaBC();
456     l0l1 = 0;//dig->GetL0L1Latency();
457
458     // get correction
459     corr = fCalib->GetTimeCorrection(index, tot, deltaBC, l0l1, timestamp); /* ps */
460     AliDebug(2, Form("calibrate index %d: time=%f (ps) tot=%f (ns) deltaBC=%d l0l1=%d timestamp=%d corr=%f (ps)",
461                      index, time, tot, deltaBC, l0l1, timestamp, corr));
462
463     // apply time correction
464     time += corr;
465
466     // convert in TDC bins and set digit
467     //tdcBin = (Int_t)(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
468     tdcBin = TMath::Nint(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
469     dig->SetTdc(tdcBin);
470
471   }
472
473   AliDebug(1,"Simulating miscalibrated digits");
474
475   return;
476 }
477
478 //---------------------------------------------------------------------
479 /*
480 void AliTOFDigitizer::DecalibrateTOFSignal(){ // Old implementation
481
482   // Read Calibration parameters from the CDB
483
484   TObjArray * calOffline= fCalib->GetTOFCalArrayOffline();
485
486   AliDebug(2,Form("Size of array for Offline Calibration = %i",calOffline->GetEntries()));
487
488   // Initialize Quantities to Simulate ToT Spectra
489
490   TH1F * hToT= fCalib->GetTOFSimToT();
491   Int_t nbins = hToT->GetNbinsX();
492   Float_t delta = hToT->GetBinWidth(1);
493   Float_t maxch = hToT->GetBinLowEdge(nbins)+delta;
494   Float_t minch = hToT->GetBinLowEdge(1);
495   Float_t max=0,min=0; //maximum and minimum value of the distribution
496   Int_t maxbin=0,minbin=0; //maximum and minimum bin of the distribution
497
498   for (Int_t ii=nbins; ii>0; ii--){
499     if (hToT->GetBinContent(ii)!= 0) {
500       max = maxch - (nbins-ii-1)*delta;
501       maxbin = ii; 
502       break;}
503   }
504   for (Int_t j=1; j<nbins; j++){
505     if (hToT->GetBinContent(j)!= 0) {
506       min = minch + (j-1)*delta;
507       minbin = j; 
508       break;}
509   }
510
511   Float_t maxToT=max;
512   Float_t minToT=min;
513  
514   Float_t maxToTDistr=hToT->GetMaximum();
515
516   AliDebug (1, Form(" The minimum ToT = %f", minToT)); 
517   AliDebug (1, Form(" The maximum ToT = %f", maxToT)); 
518   AliDebug (1, Form(" The maximum peak in ToT = %f", maxToTDistr)); 
519   
520   // Loop on TOF Digits
521
522   Bool_t isToTSimulated=kFALSE;
523   Bool_t misCalibPars=kFALSE;
524   if(hToT->GetEntries()>0)isToTSimulated=kTRUE;  
525   Int_t ndigits = fDigits->GetEntriesFast();
526   for (Int_t i=0;i<ndigits;i++){
527     AliTOFdigit * dig = (AliTOFdigit*)fDigits->At(i);
528     Int_t detId[5];
529     detId[0] = dig->GetSector();
530     detId[1] = dig->GetPlate();
531     detId[2] = dig->GetStrip();
532     detId[3] = dig->GetPadz();
533     detId[4] = dig->GetPadx();
534     dig->SetTdcND(dig->GetTdc()); // save the non decalibrated time
535     if(isToTSimulated){  
536
537       //A realistic ToT Spectrum was found in input, 
538       //decalibrated TOF Digits likely to be simulated....
539  
540       Int_t index = AliTOFGeometry::GetIndex(detId); // The channel index    
541       AliTOFChannelOffline *calChannelOffline = (AliTOFChannelOffline *)calOffline->At(index); //retrieve the info for time slewing 
542       Double_t par[6];  // time slewing parameters
543   
544       //check whether we actually ask for miscalibration
545
546       for (Int_t j = 0; j<6; j++){
547         par[j]=(Double_t)calChannelOffline->GetSlewPar(j);
548         if(par[j]!=0)misCalibPars=kTRUE;
549       }
550       AliDebug(2,Form(" Calib Pars = %f (0-th parameter for time slewing + time delay), %f, %f, %f, %f, %f ",par[0],par[1],par[2],par[3],par[4],par[5]));
551
552       // Now generate Realistic ToT distribution from TestBeam Data. 
553       // Tot is in ns, assuming a Matching Window of 10 ns.
554
555       Float_t simToT = 0;
556       Float_t trix = 0;
557       Float_t triy = 0;
558       Double_t timeCorr;
559       Double_t tToT;
560       while (simToT <= triy){
561         trix = gRandom->Rndm(i);
562         triy = gRandom->Rndm(i);
563         trix = (maxToT-minToT)*trix + minToT; 
564         triy = maxToTDistr*triy;
565         Int_t binx=hToT->FindBin(trix);
566         simToT=hToT->GetBinContent(binx);
567       }
568       // the generated ToT (ns)
569       tToT= (Double_t) trix; // to apply slewing we start from ns..
570       // transform TOF signal in ns
571       AliDebug(2,Form(" The Initial Time (counts): %i: ",dig->GetTdc()));
572       AliDebug(2,Form(" Time before miscalibration (ps) %e: ",dig->GetTdc()*(Double_t)AliTOFGeometry::TdcBinWidth()));
573       // add slewing effect
574       timeCorr=par[0] + tToT*(par[1] +tToT*(par[2] +tToT*(par[3] +tToT*(par[4] +tToT*par[5])))); 
575       AliDebug(2,Form(" The Time slewing + delay (ns): %f: ",timeCorr));
576       // add global time shift
577       //convert to ps
578       timeCorr*=1E3;
579       Double_t timeMis = (Double_t)(dig->GetTdc())*(Double_t)AliTOFGeometry::TdcBinWidth();
580       timeMis = timeMis+timeCorr;
581       AliDebug(2,Form(" The Miscalibrated time (ps): %e: ",timeMis));
582
583       // now update the digit info
584  
585       Int_t tdcCorr= (Int_t)(timeMis/AliTOFGeometry::TdcBinWidth());
586       AliDebug(2,Form(" Final Time (counts): %i: ",tdcCorr));
587       // Setting Decalibrated Time signal (TDC counts)    
588       dig->SetTdc(tdcCorr);   
589       // Setting realistic ToT signal (TDC counts) 
590       tToT*=1E3; //back to ps  
591       Int_t tot=(Int_t)(tToT/AliTOFGeometry::ToTBinWidth());//(factor 1E3 as input ToT is in ns)
592       dig->SetToT(tot); 
593       AliDebug(2,Form(" Final Time and ToT (counts): %d: , %d:",dig->GetTdc(),dig->GetToT()));
594       if(tdcCorr<0){
595         AliWarning (Form(" The bad Slewed Time(TDC counts)= %d ", tdcCorr)); 
596         AliWarning(Form(" The bad ToT (TDC counts)= %d ", tot)); 
597       }
598     }
599     else{
600     // For Data with no Miscalibration, set ToT signal == Adc
601       dig->SetToT((Int_t)(dig->GetAdc()/AliTOFGeometry::ToTBinWidth())); //remove the factor 10^3 just to have a reasonable ToT range for raw data simulation even in the case of non-realistic ToT distribution (n.b. fAdc is practically an arbitrary quantity, and ToT has no impact on the TOF reco for non-miscalibrated digits)
602     }
603   }
604
605   if(!isToTSimulated)
606     AliDebug(1,"Standard Production, no miscalibrated digits");
607   else
608     if(!misCalibPars)
609       AliDebug(1,"Standard Production, no miscalibrated digits");
610     else
611       AliDebug(1,"Simulating miscalibrated digits");
612
613   return;
614 }
615 */