]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFDigitizer.cxx
adding separate primary vertex and V0 finder components (Timur)
[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   char branchname[20];
154   sprintf (branchname, "%s", tof->GetName ());
155  
156   AliRunLoader* outrl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
157   if (outrl == 0x0)
158    {
159      AliError("Can not find Run Loader in output folder.");
160      return;
161    }
162    
163   AliLoader* outgime = outrl->GetLoader("TOFLoader");
164   if (outgime == 0x0)
165    {
166      AliError("Can not get TOF Loader from Output Run Loader.");
167      return;
168    }
169   
170   TTree* treeD = outgime->TreeD();
171   if (treeD == 0x0)
172    {
173      outgime->MakeTree("D");
174      treeD = outgime->TreeD();
175    }
176   //Make branch for digits (to be created in Init())
177   tof->MakeBranchInTree(treeD,branchname,&fDigits,4000);
178
179   // container for all summed sdigits (to be created in Init())
180   //fSDigitsArray=new TClonesArray("AliTOFSDigit",1000);
181   
182   // create hit map (to be created in Init())
183   fhitMap = new AliTOFHitMap(fSDigitsArray);
184   
185   // Loop over files to digitize
186
187   for (Int_t inputFile=0; inputFile<fManager->GetNinputs();
188        inputFile++) {
189     ReadSDigit(inputFile);
190    }
191
192   // create digits
193   CreateDigits();
194
195   // free used memory for Hit Map in current event
196   delete fhitMap;
197   fSDigitsArray->Clear();
198
199   treeD->Fill();
200
201   AliDebug(2,"----------------------------------------");
202   AliDebug(1,Form("%d digits have been created", fDigits->GetEntriesFast()));
203   AliDebug(2,"----------------------------------------");
204
205   outgime->WriteDigits("OVERWRITE");
206   outgime->UnloadDigits();
207   fDigits->Clear();
208
209 }
210
211 //---------------------------------------------------------------------
212
213 void AliTOFDigitizer::CreateDigits()
214 {
215   // loop on sdigits container to fill the AliTOFdigit TClonesArray
216   // start digitizing all the collected sdigits 
217
218   Int_t ndump=0; // dump the first ndump created digits for each event
219
220   // get the total number of collected sdigits
221   Int_t ndig = fSDigitsArray->GetEntriesFast();
222
223   Int_t  vol[5]={-1,-1,-1,-1,-1};  // location for a digit
224   Int_t  digit[4] = {0,0,0,0};     // TOF digit variables
225   Int_t tracknum[AliTOFSDigit::kMAXDIGITS]; // contributing tracks for the current slot
226   for (Int_t aa=0; aa<AliTOFSDigit::kMAXDIGITS; aa++) tracknum[aa] = -1;
227
228   for (Int_t k = 0; k < ndig; k++) {
229     
230     for (Int_t i=0; i<5; i++) vol[i] = -1;
231     
232     // Get the information for this digit
233     AliTOFSDigit *tofsdigit = (AliTOFSDigit *) fSDigitsArray->UncheckedAt(k);
234     
235     Int_t nslot=tofsdigit->GetNDigits(); // get the number of slots
236     // for current sdigit
237     
238     // TOF sdigit volumes (always the same for all slots)
239     Int_t sector    = tofsdigit->GetSector(); // range [0-17]
240     Int_t plate     = tofsdigit->GetPlate();  // range [0- 4]
241     Int_t strip     = tofsdigit->GetStrip();  // range [0-14/18/19]
242     Int_t padz      = tofsdigit->GetPadz();   // range [0- 1]
243     Int_t padx      = tofsdigit->GetPadx();   // range [0-47]
244     
245     vol[0] = sector;
246     vol[1] = plate;
247     vol[2] = strip;
248     vol[3] = padx;
249     vol[4] = padz;
250     
251     //--------------------- QA section ----------------------
252     // in the while, I perform QA
253     Bool_t isSDigitBad = (sector<0 || sector>17 || plate<0 || plate >4 || padz<0 || padz>1 || padx<0 || padx>47);
254     
255     if (isSDigitBad)
256       AliFatal(Form("strange sdigit found   %2d  %1d  %2d  %1d %2d", sector, plate, strip, padz, padx));
257     //-------------------------------------------------------
258     
259     //------------------- Dump section ----------------------
260     if (k<ndump) {
261       AliInfo(Form("%2d-th digit: Sector %2d | Plate %1d | Strip %2d | PadZ %1d | PadX %2d ", k, sector, plate, strip, padz, padx));
262       AliInfo("----------------------------------------------------");
263     }
264     // ------------------------------------------------------
265     
266     // start loop on number of slots for current sdigit
267     for (Int_t islot = 0; islot < nslot; islot++) {
268       for (Int_t aa=0; aa<4; aa++) digit[aa] = 0; // TOF digit variables
269       for (Int_t aa=0; aa<AliTOFSDigit::kMAXDIGITS; aa++) tracknum[aa] = -1;
270       
271       Int_t tdc=tofsdigit->GetTdc(islot); digit[0]=tdc;
272       Int_t adc=tofsdigit->GetAdc(islot); digit[1]=adc;
273
274       //if (tdc>=8192) continue;//AdC
275
276       tracknum[0]=tofsdigit->GetTrack(islot,0);
277       tracknum[1]=tofsdigit->GetTrack(islot,1);
278       tracknum[2]=tofsdigit->GetTrack(islot,2);
279       
280       // new with placement must be used
281       // adding a TOF digit for each slot
282       TClonesArray &aDigits = *fDigits;
283       Int_t last=fDigits->GetEntriesFast();
284       new (aDigits[last]) AliTOFdigit(tracknum, vol, digit);
285
286     }
287     
288   } // end loop on sdigits - end digitizing all collected sdigits
289
290   //Insert Decalibration 
291   AliDebug(2,"in digitizer, create digits");
292   DecalibrateTOFSignal();
293 }
294
295 //---------------------------------------------------------------------
296
297 void AliTOFDigitizer::ReadSDigit(Int_t inputFile )
298 {
299   // Read sdigits for current event and inputFile; 
300   // store them into the sdigits container
301   // and update the hit map
302   // SDigits from different files are assumed to
303   // be created with the same simulation parameters.
304   
305   // creating the TClonesArray to store the digits
306   static TClonesArray sdigitsClonesArray("AliTOFSDigit",  1000); 
307   sdigitsClonesArray.Clear();
308
309   // get the treeS from manager
310   AliRunLoader* rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
311   if (rl == 0x0)
312    {
313      AliError(Form("Can not find Run Loader in input %d folder.",inputFile));
314      return;
315    }
316
317   AliLoader* gime = rl->GetLoader("TOFLoader");
318   if (gime == 0x0)
319    {
320      AliError(Form("Can not get TOF Loader from Input %d Run Loader.",inputFile));
321      return;
322    }
323
324   TTree* currentTreeS=gime->TreeS();
325   if (currentTreeS == 0x0)
326    {
327      Int_t retval = gime->LoadSDigits();
328      if (retval) 
329       {
330          AliError(Form("Error occured while loading S. Digits for Input %d",inputFile));
331          return;
332       }
333      currentTreeS=gime->TreeS();
334      if (currentTreeS == 0x0)
335       {
336          AliError(Form("Can not get S. Digits Tree for Input %d",inputFile));
337          return;
338       }
339    } 
340   // get the branch TOF inside the treeS
341   TClonesArray * sdigitsDummyContainer=&sdigitsClonesArray;
342   // check if the branch exist
343   TBranch* tofBranch=currentTreeS->GetBranch("TOF");
344
345   if(!tofBranch){
346     AliFatal(Form("TOF branch not found for input %d",inputFile));
347   }
348   
349   tofBranch->SetAddress(&sdigitsDummyContainer);           
350   
351   Int_t nEntries = (Int_t)tofBranch->GetEntries();                                
352
353   // Loop through all entries in the tree
354   Int_t nbytes = 0;
355   
356   Int_t  vol[5]; // location for a sdigit
357   for (Int_t i=0; i<5; i++) vol[i] = -1;
358
359   for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
360     
361     // Import the tree
362     nbytes += tofBranch->GetEvent(iEntry);
363     
364     // Get the number of sdigits
365     Int_t ndig = sdigitsDummyContainer->GetEntriesFast();
366     
367     for (Int_t k=0; k<ndig; k++) {
368       AliTOFSDigit *tofSdigit= (AliTOFSDigit*) sdigitsDummyContainer->UncheckedAt(k);
369       
370       for (Int_t i=0; i<5; i++) vol[i] = -1;
371
372       // check the sdigit volume
373       vol[0] = tofSdigit->GetSector();
374       vol[1] = tofSdigit->GetPlate();
375       vol[2] = tofSdigit->GetStrip();
376       vol[3] = tofSdigit->GetPadx();
377       vol[4] = tofSdigit->GetPadz();
378       
379       if (fhitMap->TestHit(vol) != kEmpty) {
380         AliTOFSDigit *sdig = static_cast<AliTOFSDigit*>(fhitMap->GetHit(vol));
381         sdig->Update(tofSdigit);
382
383       } else {
384
385         CollectSDigit(tofSdigit); // collect the current sdigit
386         fhitMap->SetHit(vol);     // update the hitmap for location vol
387
388       } // if (hitMap->TestHit(vol) != kEmpty)
389       
390     } // for (Int_t k=0; k<ndig; k++)
391
392   } // end loop on entries
393
394 }
395
396
397 //_____________________________________________________________________________
398 void AliTOFDigitizer::CollectSDigit(const AliTOFSDigit * const sdigit)
399 {
400   //
401   // Add a TOF sdigit in container
402   // new with placement must be used
403   TClonesArray &aSDigitsArray = *fSDigitsArray;
404   Int_t last=fSDigitsArray->GetEntriesFast();
405   // make a copy of the current sdigit and
406   // put it into tmp array
407   new (aSDigitsArray[last]) AliTOFSDigit(*sdigit);
408 }
409
410 //_____________________________________________________________________________
411 void AliTOFDigitizer::InitDecalibration() const {
412   //
413   // Initialize TOF digits decalibration
414   //
415
416   fCalib->Init();
417   /*
418   fCalib->CreateCalArrays();
419   fCalib->ReadSimHistoFromCDB("TOF/Calib", -1); // use AliCDBManager's number
420   fCalib->ReadParOfflineFromCDB("TOF/Calib", -1); // use AliCDBManager's number
421   */
422 }
423 //---------------------------------------------------------------------
424 void AliTOFDigitizer::DecalibrateTOFSignal() {
425   //
426   // Decalibrate TOF signals according to OCDB parameters
427   //
428
429   Double_t time=0., tot=0., corr=0.;
430   Int_t deltaBC=0, l0l1=0, tdcBin=0;
431   Int_t index = -1;
432   Int_t detId[5] ={-1,-1,-1,-1,-1};
433   UInt_t timestamp=0;
434
435   Int_t ndigits = fDigits->GetEntriesFast();
436   // Loop on TOF Digits
437   for (Int_t i=0;i<ndigits;i++){
438     AliTOFdigit * dig = (AliTOFdigit*)fDigits->At(i);
439     detId[0] = dig->GetSector();
440     detId[1] = dig->GetPlate();
441     detId[2] = dig->GetStrip();
442     detId[3] = dig->GetPadz();
443     detId[4] = dig->GetPadx();
444     dig->SetTdcND(dig->GetTdc()); // save the non decalibrated time
445
446     index = AliTOFGeometry::GetIndex(detId); // The channel index    
447
448     // Read Calibration parameters from the CDB
449     // get digit info
450     time = dig->GetTdc() * AliTOFGeometry::TdcBinWidth(); /* ps */
451     tot = dig->GetToT() * AliTOFGeometry::ToTBinWidth() * 1.e-3; /* ns */
452     deltaBC = 0;//dig->GetDeltaBC();
453     l0l1 = 0;//dig->GetL0L1Latency();
454
455     // get correction
456     corr = fCalib->GetTimeCorrection(index, tot, deltaBC, l0l1, timestamp); /* ps */
457     AliDebug(2, Form("calibrate index %d: time=%f (ps) tot=%f (ns) deltaBC=%d l0l1=%d timestamp=%d corr=%f (ps)",
458                      index, time, tot, deltaBC, l0l1, timestamp, corr));
459
460     // apply time correction
461     time += corr;
462
463     // convert in TDC bins and set digit
464     //tdcBin = (Int_t)(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
465     tdcBin = TMath::Nint(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
466     dig->SetTdc(tdcBin);
467
468   }
469
470   AliDebug(1,"Simulating miscalibrated digits");
471
472   return;
473 }
474
475 //---------------------------------------------------------------------
476 /*
477 void AliTOFDigitizer::DecalibrateTOFSignal(){ // Old implementation
478
479   // Read Calibration parameters from the CDB
480
481   TObjArray * calOffline= fCalib->GetTOFCalArrayOffline();
482
483   AliDebug(2,Form("Size of array for Offline Calibration = %i",calOffline->GetEntries()));
484
485   // Initialize Quantities to Simulate ToT Spectra
486
487   TH1F * hToT= fCalib->GetTOFSimToT();
488   Int_t nbins = hToT->GetNbinsX();
489   Float_t delta = hToT->GetBinWidth(1);
490   Float_t maxch = hToT->GetBinLowEdge(nbins)+delta;
491   Float_t minch = hToT->GetBinLowEdge(1);
492   Float_t max=0,min=0; //maximum and minimum value of the distribution
493   Int_t maxbin=0,minbin=0; //maximum and minimum bin of the distribution
494
495   for (Int_t ii=nbins; ii>0; ii--){
496     if (hToT->GetBinContent(ii)!= 0) {
497       max = maxch - (nbins-ii-1)*delta;
498       maxbin = ii; 
499       break;}
500   }
501   for (Int_t j=1; j<nbins; j++){
502     if (hToT->GetBinContent(j)!= 0) {
503       min = minch + (j-1)*delta;
504       minbin = j; 
505       break;}
506   }
507
508   Float_t maxToT=max;
509   Float_t minToT=min;
510  
511   Float_t maxToTDistr=hToT->GetMaximum();
512
513   AliDebug (1, Form(" The minimum ToT = %f", minToT)); 
514   AliDebug (1, Form(" The maximum ToT = %f", maxToT)); 
515   AliDebug (1, Form(" The maximum peak in ToT = %f", maxToTDistr)); 
516   
517   // Loop on TOF Digits
518
519   Bool_t isToTSimulated=kFALSE;
520   Bool_t misCalibPars=kFALSE;
521   if(hToT->GetEntries()>0)isToTSimulated=kTRUE;  
522   Int_t ndigits = fDigits->GetEntriesFast();
523   for (Int_t i=0;i<ndigits;i++){
524     AliTOFdigit * dig = (AliTOFdigit*)fDigits->At(i);
525     Int_t detId[5];
526     detId[0] = dig->GetSector();
527     detId[1] = dig->GetPlate();
528     detId[2] = dig->GetStrip();
529     detId[3] = dig->GetPadz();
530     detId[4] = dig->GetPadx();
531     dig->SetTdcND(dig->GetTdc()); // save the non decalibrated time
532     if(isToTSimulated){  
533
534       //A realistic ToT Spectrum was found in input, 
535       //decalibrated TOF Digits likely to be simulated....
536  
537       Int_t index = AliTOFGeometry::GetIndex(detId); // The channel index    
538       AliTOFChannelOffline *calChannelOffline = (AliTOFChannelOffline *)calOffline->At(index); //retrieve the info for time slewing 
539       Double_t par[6];  // time slewing parameters
540   
541       //check whether we actually ask for miscalibration
542
543       for (Int_t j = 0; j<6; j++){
544         par[j]=(Double_t)calChannelOffline->GetSlewPar(j);
545         if(par[j]!=0)misCalibPars=kTRUE;
546       }
547       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]));
548
549       // Now generate Realistic ToT distribution from TestBeam Data. 
550       // Tot is in ns, assuming a Matching Window of 10 ns.
551
552       Float_t simToT = 0;
553       Float_t trix = 0;
554       Float_t triy = 0;
555       Double_t timeCorr;
556       Double_t tToT;
557       while (simToT <= triy){
558         trix = gRandom->Rndm(i);
559         triy = gRandom->Rndm(i);
560         trix = (maxToT-minToT)*trix + minToT; 
561         triy = maxToTDistr*triy;
562         Int_t binx=hToT->FindBin(trix);
563         simToT=hToT->GetBinContent(binx);
564       }
565       // the generated ToT (ns)
566       tToT= (Double_t) trix; // to apply slewing we start from ns..
567       // transform TOF signal in ns
568       AliDebug(2,Form(" The Initial Time (counts): %i: ",dig->GetTdc()));
569       AliDebug(2,Form(" Time before miscalibration (ps) %e: ",dig->GetTdc()*(Double_t)AliTOFGeometry::TdcBinWidth()));
570       // add slewing effect
571       timeCorr=par[0] + tToT*(par[1] +tToT*(par[2] +tToT*(par[3] +tToT*(par[4] +tToT*par[5])))); 
572       AliDebug(2,Form(" The Time slewing + delay (ns): %f: ",timeCorr));
573       // add global time shift
574       //convert to ps
575       timeCorr*=1E3;
576       Double_t timeMis = (Double_t)(dig->GetTdc())*(Double_t)AliTOFGeometry::TdcBinWidth();
577       timeMis = timeMis+timeCorr;
578       AliDebug(2,Form(" The Miscalibrated time (ps): %e: ",timeMis));
579
580       // now update the digit info
581  
582       Int_t tdcCorr= (Int_t)(timeMis/AliTOFGeometry::TdcBinWidth());
583       AliDebug(2,Form(" Final Time (counts): %i: ",tdcCorr));
584       // Setting Decalibrated Time signal (TDC counts)    
585       dig->SetTdc(tdcCorr);   
586       // Setting realistic ToT signal (TDC counts) 
587       tToT*=1E3; //back to ps  
588       Int_t tot=(Int_t)(tToT/AliTOFGeometry::ToTBinWidth());//(factor 1E3 as input ToT is in ns)
589       dig->SetToT(tot); 
590       AliDebug(2,Form(" Final Time and ToT (counts): %d: , %d:",dig->GetTdc(),dig->GetToT()));
591       if(tdcCorr<0){
592         AliWarning (Form(" The bad Slewed Time(TDC counts)= %d ", tdcCorr)); 
593         AliWarning(Form(" The bad ToT (TDC counts)= %d ", tot)); 
594       }
595     }
596     else{
597     // For Data with no Miscalibration, set ToT signal == Adc
598       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)
599     }
600   }
601
602   if(!isToTSimulated)
603     AliDebug(1,"Standard Production, no miscalibrated digits");
604   else
605     if(!misCalibPars)
606       AliDebug(1,"Standard Production, no miscalibrated digits");
607     else
608       AliDebug(1,"Simulating miscalibrated digits");
609
610   return;
611 }
612 */