]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRawUtils.cxx
593cca6dd4f6be870bc12d39a5f586038c18e1d7
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALRawUtils.cxx
1 // -*- mode: c++ -*-
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /* $Id$ */
18
19
20 //_________________________________________________________________________
21 //  Utility Class for handling Raw data
22 //  Does all transitions from Digits to Raw and vice versa, 
23 //  for simu and reconstruction
24 //
25 //  Note: the current version is still simplified. Only 
26 //    one raw signal per digit is generated; either high-gain or low-gain
27 //    Need to add concurrent high and low-gain info in the future
28 //    No pedestal is added to the raw signal.
29 //*-- Author: Marco van Leeuwen (LBL)
30
31 #include "AliEMCALRawUtils.h"
32 #include <stdexcept>
33   
34 #include "TF1.h"
35 #include "TGraph.h"
36 #include <TRandom.h>
37 class TSystem;
38   
39 class AliLog;
40 #include "AliRun.h"
41 #include "AliRunLoader.h"
42 class AliCaloAltroMapping;
43 #include "AliAltroBuffer.h"
44 #include "AliRawReader.h"
45 #include "AliCaloRawStreamV3.h"
46 #include "AliDAQ.h"
47   
48 #include "AliEMCALRecParam.h"
49 #include "AliEMCALLoader.h"
50 #include "AliEMCALGeometry.h"
51 class AliEMCALDigitizer;
52 #include "AliEMCALDigit.h"
53 #include "AliEMCALRawDigit.h"
54 #include "AliEMCAL.h"
55 #include "AliCaloCalibPedestal.h"  
56 #include "AliCaloFastAltroFitv0.h"
57 #include "AliCaloNeuralFit.h"
58 #include "AliCaloBunchInfo.h"
59 #include "AliCaloFitResults.h"
60 #include "AliCaloRawAnalyzerFastFit.h"
61 #include "AliCaloRawAnalyzerNN.h"
62 #include "AliCaloRawAnalyzerLMS.h"
63 #include "AliCaloRawAnalyzerPeakFinder.h"
64 #include "AliCaloRawAnalyzerCrude.h"
65 #include "AliEMCALTriggerRawDigitMaker.h"
66 #include "AliEMCALTriggerSTURawStream.h"
67 #include "AliEMCALTriggerData.h"
68
69 ClassImp(AliEMCALRawUtils)
70   
71 // Signal shape parameters
72 Int_t    AliEMCALRawUtils::fgTimeBins     = 256;       // number of sampling bins of the raw RO signal (we typically use 15-50; theoretical max is 1k+) 
73 Double_t AliEMCALRawUtils::fgTimeBinWidth = 100E-9 ;   // each sample is 100 ns
74 Double_t AliEMCALRawUtils::fgTimeTrigger  = 600E-9 ;   // the time of the trigger as approximately seen in the data
75
76 // some digitization constants
77 Int_t    AliEMCALRawUtils::fgThreshold         = 1;
78 Int_t    AliEMCALRawUtils::fgDDLPerSuperModule = 2;  // 2 ddls per SuperModule
79 Int_t    AliEMCALRawUtils::fgPedestalValue     = 0;  // pedestal value for digits2raw, default generate ZS data
80 Double_t AliEMCALRawUtils::fgFEENoise          = 3.; // 3 ADC channels of noise (sampled)
81
82 AliEMCALRawUtils::AliEMCALRawUtils( Algo::fitAlgorithm fitAlgo)
83   : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
84     fNPedSamples(0), fGeom(0), fOption(""),
85     fRemoveBadChannels(kTRUE),fFittingAlgorithm(0),  
86     fTimeMin(-1.),fTimeMax(1.),
87     fUseFALTRO(kFALSE),fRawAnalyzer(0),
88     fTriggerRawDigitMaker(0x0)
89 {
90
91   //These are default parameters.  
92   //Can be re-set from without with setter functions
93   //Already set in the OCDB and passed via setter in the AliEMCALReconstructor
94   fHighLowGainFactor = 16. ;   // Adjusted for a low gain range of 82 GeV (10 bits) 
95   fOrder             = 2;      // Order of gamma fn
96   fTau               = 2.35;   // in units of timebin, from CERN 2007 testbeam
97   fNoiseThreshold    = 3;      // 3 ADC counts is approx. noise level
98   fNPedSamples       = 4;      // Less than this value => likely pedestal samples
99   fRemoveBadChannels = kFALSE; // Do not remove bad channels before fitting
100   fUseFALTRO         = kTRUE;  // Get the trigger FALTRO information and pass it to digits.
101   SetFittingAlgorithm(fitAlgo);
102
103   //Get Mapping RCU files from the AliEMCALRecParam                                 
104   const TObjArray* maps = AliEMCALRecParam::GetMappings();
105   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
106
107   for(Int_t i = 0; i < 4; i++) {
108     fMapping[i] = (AliAltroMapping*)maps->At(i);
109   }
110
111   //To make sure we match with the geometry in a simulation file,
112   //let's try to get it first.  If not, take the default geometry
113   AliRunLoader *rl = AliRunLoader::Instance();
114   if (rl && rl->GetAliRun()) {
115     AliEMCAL * emcal = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"));
116     if(emcal)fGeom = emcal->GetGeometry();
117     else {
118       AliDebug(1, Form("Using default geometry in raw reco"));
119       fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
120     }
121
122   } else {
123     AliDebug(1, Form("Using default geometry in raw reco"));
124     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
125   }
126
127   if(!fGeom) AliFatal(Form("Could not get geometry!"));
128         
129   fTriggerRawDigitMaker = new AliEMCALTriggerRawDigitMaker();
130
131 }
132
133 //____________________________________________________________________________
134 AliEMCALRawUtils::AliEMCALRawUtils(AliEMCALGeometry *pGeometry, Algo::fitAlgorithm fitAlgo)
135   : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
136     fNPedSamples(0), fGeom(pGeometry), fOption(""),
137     fRemoveBadChannels(kTRUE),fFittingAlgorithm(0),
138     fTimeMin(-1.),fTimeMax(1.),
139     fUseFALTRO(kFALSE),fRawAnalyzer(),
140     fTriggerRawDigitMaker(0x0)
141 {
142   //
143   // Initialize with the given geometry - constructor required by HLT
144   // HLT does not use/support AliRunLoader(s) instances
145   // This is a minimum intervention solution
146   // Comment by MPloskon@lbl.gov
147   //
148
149   //These are default parameters. 
150   //Can be re-set from without with setter functions 
151   //Already set in the OCDB and passed via setter in the AliEMCALReconstructor
152   fHighLowGainFactor = 16. ;   // adjusted for a low gain range of 82 GeV (10 bits)
153   fOrder             = 2;      // order of gamma fn
154   fTau               = 2.35;   // in units of timebin, from CERN 2007 testbeam
155   fNoiseThreshold    = 3;      // 3 ADC counts is approx. noise level
156   fNPedSamples       = 4;      // Less than this value => likely pedestal samples
157   fRemoveBadChannels = kFALSE; // Do not remove bad channels before fitting
158   fUseFALTRO         = kTRUE;  // Get the trigger FALTRO information and pass it to digits.
159   SetFittingAlgorithm(fitAlgo);
160
161   //Get Mapping RCU files from the AliEMCALRecParam
162   const TObjArray* maps = AliEMCALRecParam::GetMappings();
163   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
164
165   for(Int_t i = 0; i < 4; i++) {
166     fMapping[i] = (AliAltroMapping*)maps->At(i);
167   }
168
169   if(!fGeom) AliFatal(Form("Could not get geometry!"));
170         
171   fTriggerRawDigitMaker = new AliEMCALTriggerRawDigitMaker();   
172 }
173
174 //____________________________________________________________________________
175 AliEMCALRawUtils::AliEMCALRawUtils(const AliEMCALRawUtils& rawU)
176   : TObject(),
177     fHighLowGainFactor(rawU.fHighLowGainFactor), 
178     fOrder(rawU.fOrder),
179     fTau(rawU.fTau),
180     fNoiseThreshold(rawU.fNoiseThreshold),
181     fNPedSamples(rawU.fNPedSamples),
182     fGeom(rawU.fGeom), 
183     fOption(rawU.fOption),
184     fRemoveBadChannels(rawU.fRemoveBadChannels),
185     fFittingAlgorithm(rawU.fFittingAlgorithm),
186     fTimeMin(rawU.fTimeMin),fTimeMax(rawU.fTimeMax),
187     fUseFALTRO(rawU.fUseFALTRO),
188     fRawAnalyzer(rawU.fRawAnalyzer),
189     fTriggerRawDigitMaker(rawU.fTriggerRawDigitMaker)
190 {
191   //copy ctor
192   fMapping[0] = rawU.fMapping[0];
193   fMapping[1] = rawU.fMapping[1];
194   fMapping[2] = rawU.fMapping[2];
195   fMapping[3] = rawU.fMapping[3];
196 }
197
198 //____________________________________________________________________________
199 AliEMCALRawUtils& AliEMCALRawUtils::operator =(const AliEMCALRawUtils &rawU)
200 {
201   //assignment operator
202
203   if(this != &rawU) {
204     fHighLowGainFactor = rawU.fHighLowGainFactor;
205     fOrder             = rawU.fOrder;
206     fTau               = rawU.fTau;
207     fNoiseThreshold    = rawU.fNoiseThreshold;
208     fNPedSamples       = rawU.fNPedSamples;
209     fGeom              = rawU.fGeom;
210     fOption            = rawU.fOption;
211     fRemoveBadChannels = rawU.fRemoveBadChannels;
212     fFittingAlgorithm  = rawU.fFittingAlgorithm;
213     fTimeMin           = rawU.fTimeMin;
214     fTimeMax           = rawU.fTimeMax;
215     fUseFALTRO         = rawU.fUseFALTRO;
216     fRawAnalyzer       = rawU.fRawAnalyzer;
217     fMapping[0]        = rawU.fMapping[0];
218     fMapping[1]        = rawU.fMapping[1];
219     fMapping[2]        = rawU.fMapping[2];
220     fMapping[3]        = rawU.fMapping[3];
221     fTriggerRawDigitMaker = rawU.fTriggerRawDigitMaker;  
222   }
223
224   return *this;
225
226 }
227
228 //____________________________________________________________________________
229 AliEMCALRawUtils::~AliEMCALRawUtils() {
230   //dtor
231
232 }
233
234 //____________________________________________________________________________
235 void AliEMCALRawUtils::Digits2Raw()
236 {
237   // convert digits of the current event to raw data
238   
239   AliRunLoader *rl = AliRunLoader::Instance();
240   AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
241   
242   // get the digits
243   loader->LoadDigits("EMCAL");
244   loader->GetEvent();
245   TClonesArray* digits = loader->Digits() ;
246   
247   if (!digits) {
248     Warning("Digits2Raw", "no digits found !");
249     return;
250   }
251   
252   static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
253   AliAltroBuffer* buffers[nDDL];
254   for (Int_t i=0; i < nDDL; i++)
255     buffers[i] = 0;
256   
257   TArrayI adcValuesLow(fgTimeBins);
258   TArrayI adcValuesHigh(fgTimeBins);
259   
260   // loop over digits (assume ordered digits)
261   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
262     AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
263     if(!digit){
264       AliFatal("NULL Digit");
265     }
266     else{
267       if (digit->GetAmplitude() < fgThreshold) 
268         continue;
269       
270       //get cell indices
271       Int_t nSM = 0;
272       Int_t nIphi = 0;
273       Int_t nIeta = 0;
274       Int_t iphi = 0;
275       Int_t ieta = 0;
276       Int_t nModule = 0;
277       fGeom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
278       fGeom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
279       
280       //Check which is the RCU, 0 or 1, of the cell.
281       Int_t iRCU = -111;
282       //RCU0
283       if (0<=iphi&&iphi<8) iRCU=0; // first cable row
284       else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half; 
285       //second cable row
286       //RCU1
287       else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half; 
288       //second cable row
289       else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
290       
291       if (nSM%2==1) iRCU = 1 - iRCU; // swap for odd=C side, to allow us to cable both sides the same
292       
293       if (iRCU<0) 
294         Fatal("Digits2Raw()","Non-existent RCU number: %d", iRCU);
295       
296       //Which DDL?
297       Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
298       if (iDDL < 0 || iDDL >= nDDL){
299         Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
300       }
301       else{
302         if (buffers[iDDL] == 0) {      
303           // open new file and write dummy header
304           TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
305           //Select mapping file RCU0A, RCU0C, RCU1A, RCU1C
306           Int_t iRCUside=iRCU+(nSM%2)*2;
307           //iRCU=0 and even (0) SM -> RCU0A.data   0
308           //iRCU=1 and even (0) SM -> RCU1A.data   1
309           //iRCU=0 and odd  (1) SM -> RCU0C.data   2
310           //iRCU=1 and odd  (1) SM -> RCU1C.data   3
311           //cout<<" nSM "<<nSM<<"; iRCU "<<iRCU<<"; iRCUside "<<iRCUside<<endl;
312           buffers[iDDL] = new AliAltroBuffer(fileName.Data(),fMapping[iRCUside]);
313           buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
314         }
315         
316         // out of time range signal (?)
317         if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
318           AliInfo("Signal is out of time range.\n");
319           buffers[iDDL]->FillBuffer((Int_t)digit->GetAmplitude());
320           buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() );  // time bin
321           buffers[iDDL]->FillBuffer(3);          // bunch length      
322           buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM);  // trailer
323           // calculate the time response function
324         } else {
325           Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmplitude(), adcValuesHigh.GetArray(), adcValuesLow.GetArray()) ; 
326           if (lowgain) 
327             buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow.GetArray(), fgThreshold);
328           else 
329             buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh.GetArray(), fgThreshold);
330         }
331       }// iDDL under the limits
332     }//digit exists
333   }//Digit loop
334   
335   // write headers and close files
336   for (Int_t i=0; i < nDDL; i++) {
337     if (buffers[i]) {
338       buffers[i]->Flush();
339       buffers[i]->WriteDataHeader(kFALSE, kFALSE);
340       delete buffers[i];
341     }
342   }
343   
344   loader->UnloadDigits();
345 }
346
347 //____________________________________________________________________________
348 void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr, const AliCaloCalibPedestal* pedbadmap, TClonesArray *digitsTRG, AliEMCALTriggerData* trgData)
349 {
350   // convert raw data of the current event to digits                                                                                     
351   
352   if(digitsArr) digitsArr->Clear("C"); 
353   
354   if (!digitsArr) {
355     Error("Raw2Digits", "no digits found !");
356     return;
357   }
358   if (!reader) {
359     Error("Raw2Digits", "no raw reader found !");
360     return;
361   }
362         
363   AliEMCALTriggerSTURawStream inSTU(reader);
364         
365   AliCaloRawStreamV3 in(reader,"EMCAL",fMapping);       
366         
367   // Select EMCAL DDL's;
368   reader->Select("EMCAL",0,43); // 43 = AliEMCALGeoParams::fgkLastAltroDDL
369   
370   fTriggerRawDigitMaker->Reset();       
371   fTriggerRawDigitMaker->SetIO(reader, in, inSTU, digitsTRG, trgData);
372         
373   // fRawAnalyzer setup
374   fRawAnalyzer->SetNsampleCut(5); // requirement for fits to be done, for the new methods
375   fRawAnalyzer->SetOverflowCut(fgkOverflowCut);
376   fRawAnalyzer->SetAmpCut(fNoiseThreshold);
377   fRawAnalyzer->SetFitArrayCut(fNoiseThreshold);
378   fRawAnalyzer->SetIsZeroSuppressed(true); // TMP - should use stream->IsZeroSuppressed(), or altro cfg registers later
379   
380   // channel info parameters
381   Int_t lowGain  = 0;
382   Int_t caloFlag = 0; // low, high gain, or TRU, or LED ref.
383   
384   // start loop over input stream 
385   while (in.NextDDL()) {
386           
387     //    if ( in.GetDDLNumber() != 0 && in.GetDDLNumber() != 2 ) continue;
388     
389     while (in.NextChannel()) {
390       
391       //Check if the signal  is high or low gain and then do the fit, 
392       //if it  is from TRU or LEDMon do not fit
393       caloFlag = in.GetCaloFlag();
394       //                if (caloFlag != 0 && caloFlag != 1) continue; 
395       if (caloFlag > 2) continue; // Work with ALTRO and FALTRO 
396       
397       //Do not fit bad channels of ALTRO
398       if(caloFlag < 2 && fRemoveBadChannels && pedbadmap->IsBadChannel(in.GetModule(),in.GetColumn(),in.GetRow())) {
399         //printf("Tower from SM %d, column %d, row %d is BAD!!! Skip \n", in.GetModule(),in.GetColumn(),in.GetRow());
400         continue;
401       }  
402       
403       vector<AliCaloBunchInfo> bunchlist; 
404       while (in.NextBunch()) {
405         bunchlist.push_back( AliCaloBunchInfo(in.GetStartTimeBin(), in.GetBunchLength(), in.GetSignals() ) );
406       } // loop over bunches
407       if (bunchlist.size() == 0) continue;
408       if ( caloFlag < 2 )
409         { // ALTRO
410           Float_t time = 0; 
411           Float_t amp  = 0; 
412           short timeEstimate  = 0;
413           Float_t ampEstimate = 0;
414           Bool_t fitDone = kFALSE;
415           Float_t chi2 = 0;
416           Int_t ndf = 0;
417
418           AliDebug(1,Form("trunk lms? %d, %d\n", fFittingAlgorithm, Algo::kLMS));
419
420           if ( fFittingAlgorithm == Algo::kFastFit || fFittingAlgorithm == Algo::kNeuralNet || 
421                fFittingAlgorithm == Algo::kLMS || fFittingAlgorithm == Algo::kPeakFinder || 
422                fFittingAlgorithm == Algo::kCrude) {
423           // all functionality to determine amp and time etc is encapsulated inside the Evaluate call for these methods 
424           AliCaloFitResults fitResults = fRawAnalyzer->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2()); 
425           
426           amp          = fitResults.GetAmp();
427           time         = fitResults.GetTime();
428           timeEstimate = fitResults.GetMaxTimebin();
429           ampEstimate  = fitResults.GetMaxSig();
430           chi2 = fitResults.GetChi2();
431           ndf = fitResults.GetNdf();
432           if (fitResults.GetStatus() == Ret::kFitPar) {
433             fitDone = kTRUE;
434           }
435         }
436         else { // for the other methods we for now use the functionality of 
437           // AliCaloRawAnalyzer as well, to select samples and prepare for fits, 
438           // if it looks like there is something to fit
439           
440           // parameters init.
441           Float_t pedEstimate  = 0;
442           short maxADC = 0;
443           Int_t first = 0;
444           Int_t last = 0;
445           Int_t bunchIndex = 0;
446           //
447           // The PreFitEvaluateSamples + later call to FitRaw will hopefully 
448           // be replaced by a single Evaluate call or so soon, like for the other
449           // methods, but this should be good enough for evaluation of 
450           // the methods for now (Jan. 2010)
451           //
452           int nsamples = fRawAnalyzer->PreFitEvaluateSamples( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2(), bunchIndex, ampEstimate, maxADC, timeEstimate, pedEstimate, first, last); 
453           
454           if (ampEstimate >= fNoiseThreshold) { // something worth looking at
455             
456             time = timeEstimate; // maxrev in AliCaloRawAnalyzer speak; comes with an offset w.r.t. real timebin
457             Int_t timebinOffset = bunchlist.at(bunchIndex).GetStartBin() - (bunchlist.at(bunchIndex).GetLength()-1); 
458             amp = ampEstimate; 
459             
460             if ( nsamples > 1 && maxADC<fgkOverflowCut ) { // possibly something to fit
461               FitRaw(first, last, amp, time, chi2, fitDone);
462               time += timebinOffset;
463               timeEstimate += timebinOffset;
464               ndf = nsamples - 2;
465             }
466             
467           } // ampEstimate check
468         } // method selection
469
470         if ( fitDone ) { // brief sanity check of fit results       
471           Float_t ampAsymm = (amp - ampEstimate)/(amp + ampEstimate);
472           Float_t timeDiff = time - timeEstimate;
473           if ( (TMath::Abs(ampAsymm) > 0.1) || (TMath::Abs(timeDiff) > 2) ) {
474             // AliDebug(2,Form("Fit results amp %f time %f not consistent with expectations amp %f time %d", amp, time, ampEstimate, timeEstimate));
475             
476             // for now just overwrite the fit results with the simple/initial estimate
477             amp     = ampEstimate;
478             time    = timeEstimate; 
479             fitDone = kFALSE;
480           } 
481         } // fitDone
482         
483         if (amp >= fNoiseThreshold) { // something to be stored
484           if ( ! fitDone) { // smear ADC with +- 0.5 uniform (avoid discrete effects)
485             amp += (0.5 - gRandom->Rndm()); // Rndm generates a number in ]0,1]
486           }
487           
488           Int_t id = fGeom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
489           lowGain  = in.IsLowGain();
490           
491           // go from time-bin units to physical time fgtimetrigger
492           time = time * GetRawFormatTimeBinWidth(); // skip subtraction of fgTimeTrigger?
493           // subtract RCU L1 phase (L1Phase is in seconds) w.r.t. L0:
494           time -= in.GetL1Phase();
495           
496           AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
497           // printf("Added tower: SM %d, row %d, column %d, amp %3.2f\n",in.GetModule(), in.GetRow(), in.GetColumn(),amp);
498
499           AddDigit(digitsArr, id, lowGain, amp, time, chi2, ndf); 
500         }
501         
502       }//ALTRO
503       else if(fUseFALTRO)
504       {// Fake ALTRO
505         fTriggerRawDigitMaker->Add( bunchlist );
506       }//Fake ALTRO
507     } // end while over channel   
508   } //end while over DDL's, of input stream 
509   
510   fTriggerRawDigitMaker->PostProcess(); 
511         
512   TrimDigits(digitsArr);
513         
514   return ; 
515 }
516
517 //____________________________________________________________________________ 
518 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Float_t amp, Float_t time, Float_t chi2, Int_t ndf) {
519   //
520   // Add a new digit. 
521   // This routine checks whether a digit exists already for this tower 
522   // and then decides whether to use the high or low gain info
523   //
524   // Called by Raw2Digits
525   
526   AliEMCALDigit *digit = 0, *tmpdigit = 0;
527   TIter nextdigit(digitsArr);
528   while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
529     if (tmpdigit->GetId() == id) digit = tmpdigit;
530   }
531
532   if (!digit) { // no digit existed for this tower; create one
533                 Int_t type = AliEMCALDigit::kHG; // use enum in AliEMCALDigit
534                 if (lowGain) { 
535                         amp *= fHighLowGainFactor;
536                         type = AliEMCALDigit::kLGnoHG;
537                 } 
538                 Int_t idigit = digitsArr->GetEntries();
539                 new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, type, idigit, chi2, ndf); 
540                 AliDebug(2,Form("Add digit Id %d for the first time, type %d", id, type));
541   }//digit added first time
542   else { // a digit already exists, check range 
543                 // (use high gain if signal < cut value, otherwise low gain)
544                 if (lowGain) { // new digit is low gain
545                         if (digit->GetAmplitude() > fgkOverflowCut) {  // use if previously stored (HG) digit is out of range
546                                 digit->SetAmplitude(fHighLowGainFactor * amp);
547                                 digit->SetTime(time);
548                                 digit->SetType(AliEMCALDigit::kLG);
549                                 AliDebug(2,Form("Add LG digit ID %d for the second time, type %d", digit->GetId(), digit->GetType()));
550                         }
551                 }//new low gain digit
552                 else { // new digit is high gain 
553                         if (amp < fgkOverflowCut) { // new digit is high gain; use if not out of range
554                                 digit->SetAmplitude(amp);
555                                 digit->SetTime(time);
556                                 digit->SetType(AliEMCALDigit::kHG);
557                                 AliDebug(2,Form("Add HG digit ID %d for the second time, type %d", digit->GetId(), digit->GetType()));
558                         }
559                         else { // HG out of range, just change flag value to show that HG did exist
560                                 digit->SetType(AliEMCALDigit::kLG);
561                                 AliDebug(2,Form("Change LG digit to HG, ID %d, type %d", digit->GetId(), digit->GetType()));
562                         }
563                 }//new high gain digit
564   }//digit existed replace it
565   
566 }
567
568 //____________________________________________________________________________ 
569 void AliEMCALRawUtils::TrimDigits(TClonesArray *digitsArr) 
570 {
571   // Remove digits with only low gain and large time
572   
573   AliEMCALDigit *digit = 0;
574   Int_t n = 0;
575   Int_t nDigits = digitsArr->GetEntriesFast();
576   TIter nextdigit(digitsArr);
577   while ((digit = (AliEMCALDigit*) nextdigit())) {
578
579     //Check if only LG existed, remove if so
580     if (digit->GetType() == AliEMCALDigit::kLGnoHG) {
581       AliDebug(10,Form("Remove digit with id %d, LGnoHG",digit->GetId()));
582       digitsArr->Remove(digit);
583     }
584     //Check if time is too large or too small, remove if so
585     else if(fTimeMin > digit->GetTime() || fTimeMax < digit->GetTime()) {
586       digitsArr->Remove(digit);
587       AliDebug(10,Form("Remove digit with id %d, Bad Time %e",digit->GetId(), digit->GetTime()));
588     }
589     // Check if Chi2 is undefined
590     else if (0 > digit->GetChi2()) {
591       digitsArr->Remove(digit);
592       AliDebug(10,Form("Remove digit with id %d, Bad Chi2 %e",digit->GetId(), digit->GetChi2()));
593     }
594     //Good digit, just reassign the index of the digit in case there was a previous removal
595     else {
596       digit->SetIndexInList(n); 
597       n++;
598     }    
599   }//while
600   
601   digitsArr->Compress();
602   AliDebug(1,Form("N Digits before trimming : %d; after array compression %d",nDigits,digitsArr->GetEntriesFast()));
603            
604 }
605         
606 //____________________________________________________________________________ 
607 void AliEMCALRawUtils::FitRaw(const Int_t firstTimeBin, const Int_t lastTimeBin, Float_t & amp, Float_t & time, Float_t & chi2, Bool_t & fitDone) const 
608 { // Fits the raw signal time distribution
609   
610   //--------------------------------------------------
611   //Do the fit, different fitting algorithms available
612   //--------------------------------------------------
613   int nsamples = lastTimeBin - firstTimeBin + 1;
614   fitDone = kFALSE;
615   
616   switch(fFittingAlgorithm) {
617   case Algo::kStandard:
618     {
619       if (nsamples < 3) { return; } // nothing much to fit
620       //printf("Standard fitter \n");
621
622       // Create Graph to hold data we will fit 
623       TGraph *gSig =  new TGraph( nsamples); 
624       for (int i=0; i<nsamples; i++) {
625         Int_t timebin = firstTimeBin + i;    
626         gSig->SetPoint(i, timebin, fRawAnalyzer->GetReversed(timebin)); 
627       }
628       
629       TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
630       signalF->SetParameters(10.,5.,fTau,fOrder,0.); //set all defaults once, just to be safe
631       signalF->SetParNames("amp","t0","tau","N","ped");
632       signalF->FixParameter(2,fTau); // tau in units of time bin
633       signalF->FixParameter(3,fOrder); // order
634       signalF->FixParameter(4, 0); // pedestal should be subtracted when we get here 
635       signalF->SetParameter(1, time);
636       signalF->SetParameter(0, amp);
637       // set rather loose parameter limits
638       signalF->SetParLimits(0, 0.5*amp, 2*amp );
639       signalF->SetParLimits(1, time - 4, time + 4); 
640       
641       try {                     
642         gSig->Fit(signalF, "QROW"); // Note option 'W': equal errors on all points
643         // assign fit results
644         amp  = signalF->GetParameter(0); 
645         time = signalF->GetParameter(1);
646         chi2 = signalF->GetChisquare();
647         fitDone = kTRUE;
648       }
649       catch (const std::exception & e) {
650         AliError( Form("TGraph Fit exception %s", e.what()) ); 
651         // stay with default amp and time in case of exception, i.e. no special action required
652         fitDone = kFALSE;
653       }
654       delete signalF;
655       
656       //printf("Std   : Amp %f, time %g\n",amp, time);
657       delete gSig; // delete TGraph
658       
659       break;
660     }//kStandard Fitter
661       //----------------------------
662   case Algo::kLogFit:
663     {
664       if (nsamples < 3) { return; } // nothing much to fit
665       //printf("LogFit \n");
666       
667       // Create Graph to hold data we will fit 
668       TGraph *gSigLog =  new TGraph( nsamples); 
669       for (int i=0; i<nsamples; i++) {
670         Int_t timebin = firstTimeBin + i;    
671         gSigLog->SetPoint(timebin, timebin, TMath::Log(fRawAnalyzer->GetReversed(timebin) ) ); 
672       }
673       
674       TF1 * signalFLog = new TF1("signalLog", RawResponseFunctionLog, 0, GetRawFormatTimeBins(), 5);
675       signalFLog->SetParameters(2.3, 5.,fTau,fOrder,0.); //set all defaults once, just to be safe
676       signalFLog->SetParNames("amplog","t0","tau","N","ped");
677       signalFLog->FixParameter(2,fTau); // tau in units of time bin
678       signalFLog->FixParameter(3,fOrder); // order
679       signalFLog->FixParameter(4, 0); // pedestal should be subtracted when we get here 
680       signalFLog->SetParameter(1, time);
681       if (amp>=1) {
682         signalFLog->SetParameter(0, TMath::Log(amp));
683       }
684       
685       gSigLog->Fit(signalFLog, "QROW"); // Note option 'W': equal errors on all points
686       
687       // assign fit results
688       Double_t amplog = signalFLog->GetParameter(0); //Not Amp, but Log of Amp
689       amp = TMath::Exp(amplog);
690       time = signalFLog->GetParameter(1);
691       fitDone = kTRUE;
692       
693       delete signalFLog;
694       //printf("LogFit: Amp %f, time %g\n",amp, time);
695       delete gSigLog; 
696       break;
697     } //kLogFit 
698       //----------------------------    
699       
700       //----------------------------
701   }//switch fitting algorithms
702   
703   return;
704 }
705
706 //__________________________________________________________________
707 void AliEMCALRawUtils::FitParabola(const TGraph *gSig, Float_t & amp) const 
708 {
709   //BEG YS alternative methods to calculate the amplitude
710   Double_t * ymx = gSig->GetX() ; 
711   Double_t * ymy = gSig->GetY() ; 
712   const Int_t kN = 3 ; 
713   Double_t ymMaxX[kN] = {0., 0., 0.} ; 
714   Double_t ymMaxY[kN] = {0., 0., 0.} ; 
715   Double_t ymax = 0. ; 
716   // find the maximum amplitude
717   Int_t ymiMax = 0 ;  
718   for (Int_t ymi = 0; ymi < gSig->GetN(); ymi++) {
719     if (ymy[ymi] > ymMaxY[0] ) {
720       ymMaxY[0] = ymy[ymi] ; //<========== This is the maximum amplitude
721       ymMaxX[0] = ymx[ymi] ;
722       ymiMax = ymi ; 
723     }
724   }
725   // find the maximum by fitting a parabola through the max and the two adjacent samples
726   if ( ymiMax < gSig->GetN()-1 && ymiMax > 0) {
727     ymMaxY[1] = ymy[ymiMax+1] ;
728     ymMaxY[2] = ymy[ymiMax-1] ; 
729     ymMaxX[1] = ymx[ymiMax+1] ;
730     ymMaxX[2] = ymx[ymiMax-1] ; 
731     if (ymMaxY[0]*ymMaxY[1]*ymMaxY[2] > 0) {
732       //fit a parabola through the 3 points y= a+bx+x*x*x
733       Double_t sy = 0 ; 
734       Double_t sx = 0 ; 
735       Double_t sx2 = 0 ; 
736       Double_t sx3 = 0 ; 
737       Double_t sx4 = 0 ; 
738       Double_t sxy = 0 ; 
739       Double_t sx2y = 0 ; 
740       for (Int_t i = 0; i < kN ; i++) {
741         sy += ymMaxY[i] ; 
742         sx += ymMaxX[i] ;               
743         sx2 += ymMaxX[i]*ymMaxX[i] ; 
744         sx3 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ; 
745         sx4 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ; 
746         sxy += ymMaxX[i]*ymMaxY[i] ; 
747         sx2y += ymMaxX[i]*ymMaxX[i]*ymMaxY[i] ; 
748       }
749       Double_t cN = (sx2y*kN-sy*sx2)*(sx3*sx-sx2*sx2)-(sx2y*sx-sxy*sx2)*(sx3*kN-sx*sx2); 
750       Double_t cD = (sx4*kN-sx2*sx2)*(sx3*sx-sx2*sx2)-(sx4*sx-sx3*sx2)*(sx3*kN-sx*sx2) ;
751       Double_t c  = cN / cD ; 
752       Double_t b  = ((sx2y*kN-sy*sx2)-c*(sx4*kN-sx2*sx2))/(sx3*kN-sx*sx2) ;
753       Double_t a  = (sy-b*sx-c*sx2)/kN  ;
754       Double_t xmax = -b/(2*c) ; 
755       ymax = a + b*xmax + c*xmax*xmax ;//<========== This is the maximum amplitude
756       amp = ymax;
757     }
758   }
759   
760   Double_t diff = TMath::Abs(1-ymMaxY[0]/amp) ; 
761   if (diff > 0.1) 
762     amp = ymMaxY[0] ; 
763   //printf("Yves   : Amp %f, time %g\n",amp, time);
764   //END YS
765   return;
766 }
767
768 //__________________________________________________________________
769 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
770 {
771   // Matches version used in 2007 beam test
772   //
773   // Shape of the electronics raw reponse:
774   // It is a semi-gaussian, 2nd order Gamma function of the general form
775   //
776   // xx = (t - t0 + tau) / tau  [xx is just a convenient help variable]
777   // F = A * (xx**N * exp( N * ( 1 - xx) )   for xx >= 0
778   // F = 0                                   for xx < 0 
779   //
780   // parameters:
781   // A:   par[0]   // Amplitude = peak value
782   // t0:  par[1]
783   // tau: par[2]
784   // N:   par[3]
785   // ped: par[4]
786   //
787   Double_t signal = 0.;
788   Double_t tau    = par[2];
789   Double_t n      = par[3];
790   Double_t ped    = par[4];
791   Double_t xx     = ( x[0] - par[1] + tau ) / tau ;
792
793   if (xx <= 0) 
794     signal = ped ;  
795   else {  
796     signal = ped + par[0] * TMath::Power(xx , n) * TMath::Exp(n * (1 - xx )) ; 
797   }
798   return signal ;  
799 }
800
801 //__________________________________________________________________
802 Double_t AliEMCALRawUtils::RawResponseFunctionLog(Double_t *x, Double_t *par)
803 {
804   // Matches version used in 2007 beam test
805   //
806   // Shape of the electronics raw reponse:
807   // It is a semi-gaussian, 2nd order Gamma function of the general form
808   //
809   // xx = (t - t0 + tau) / tau  [xx is just a convenient help variable]
810   // F = A * (xx**N * exp( N * ( 1 - xx) )   for xx >= 0
811   // F = 0                                   for xx < 0 
812   //
813   // parameters:
814   // Log[A]:   par[0]   // Amplitude = peak value
815   // t0:  par[1]
816   // tau: par[2]
817   // N:   par[3]
818   // ped: par[4]
819   //
820   Double_t signal = 0. ;
821   Double_t tau    = par[2];
822   Double_t n      = par[3];
823   //Double_t ped = par[4]; // not used
824   Double_t xx     = ( x[0] - par[1] + tau ) / tau ;
825
826   if (xx < 0) 
827     signal = par[0] - n*TMath::Log(TMath::Abs(xx)) + n * (1 - xx ) ;  
828   else {  
829     signal = par[0] + n*TMath::Log(xx) + n * (1 - xx ) ; 
830   }
831   return signal ;  
832 }
833
834 //__________________________________________________________________
835 Bool_t AliEMCALRawUtils::RawSampledResponse(const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL, const Int_t keyErr) const 
836 {
837   // for a start time dtime and an amplitude damp given by digit, 
838   // calculates the raw sampled response AliEMCAL::RawResponseFunction
839   Bool_t lowGain = kFALSE ; 
840   
841   // A:   par[0]   // Amplitude = peak value
842   // t0:  par[1]                            
843   // tau: par[2]                            
844   // N:   par[3]                            
845   // ped: par[4]
846   
847   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
848   signalF.SetParameter(0, damp) ; 
849   signalF.SetParameter(1, (dtime + fgTimeTrigger)/fgTimeBinWidth) ; 
850   signalF.SetParameter(2, fTau) ; 
851   signalF.SetParameter(3, fOrder);
852   signalF.SetParameter(4, fgPedestalValue);
853         
854   Double_t signal=0.0, noise=0.0;
855   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
856     signal = signalF.Eval(iTime) ;  
857     // Next lines commeted for the moment but in principle it is not necessary to add
858     // extra noise since noise already added at the digits level.       
859     
860     //According to Terry Awes, 13-Apr-2008
861     //add gaussian noise in quadrature to each sample
862     //Double_t noise = gRandom->Gaus(0.,fgFEENoise);
863     //signal = sqrt(signal*signal + noise*noise);
864     
865     // March 17,09 for fast fit simulations by Alexei Pavlinov.
866     // Get from PHOS analysis. In some sense it is open questions.
867     if(keyErr>0) {
868       noise = gRandom->Gaus(0.,fgFEENoise);
869       signal += noise; 
870     }
871           
872     adcH[iTime] =  static_cast<Int_t>(signal + 0.5) ;
873     if ( adcH[iTime] > fgkRawSignalOverflow ){  // larger than 10 bits 
874       adcH[iTime] = fgkRawSignalOverflow ;
875       lowGain = kTRUE ; 
876     }
877     
878     signal /= fHighLowGainFactor;
879     
880     adcL[iTime] =  static_cast<Int_t>(signal + 0.5) ;
881     if ( adcL[iTime] > fgkRawSignalOverflow)  // larger than 10 bits 
882       adcL[iTime] = fgkRawSignalOverflow ;
883         
884   }
885   
886   return lowGain ; 
887 }
888
889 //__________________________________________________________________
890 void AliEMCALRawUtils::CalculateChi2(const Double_t* t, const Double_t* y, const Int_t nPoints, 
891 const Double_t sig, const Double_t tau, const Double_t amp, const Double_t t0, Double_t &chi2)
892 {
893   //   Input:
894   //   t[]   - array of time bins
895   //   y[]   - array of amplitudes after pedestal subtractions;
896   //   nPoints  - number of points 
897   //   sig   - error of amplitude measurement (one value for all channels)
898   //           if sig<0 that mean sig=1.
899   //   tau   - filter time response (in timebin units)
900   //   amp   - amplitude at t0;
901   //   t0    - time of max amplitude; 
902   // Output:
903   //   chi2 - chi2
904   //   ndf = nPoints - 2 when tau fixed 
905   //   ndf = nPoints - 3 when tau free
906   static Double_t par[5]={0.0, 0.0, 0.0, 2.0, 0.0};
907
908   par[0] = amp;
909   par[1] = t0;
910   par[2] = tau;
911   // par[3]=n=2.; par[4]=ped=0.0
912
913   Double_t dy = 0.0, x = 0.0, f=0.0;
914   for(Int_t i=0; i<nPoints; i++){
915     x     = t[i];
916     f     = RawResponseFunction(&x, par);
917     dy    = y[i] - f;
918     chi2 += dy*dy;
919     //printf(" AliEMCALRawUtils::CalculateChi2 : %i : y %f -> f %f : dy %f \n", i, y[i], f, dy); 
920   }
921   if(sig>0.0) chi2 /= (sig*sig);
922 }
923
924 //__________________________________________________________________
925 void AliEMCALRawUtils::SetFittingAlgorithm(Int_t fitAlgo)              
926 {
927         //Set fitting algorithm and initialize it if this same algorithm was not set before.
928         //printf("**** Set Algorithm , number %d ****\n",fitAlgo);
929
930         if(fitAlgo == fFittingAlgorithm && fRawAnalyzer) {
931                 //Do nothing, this same algorithm already set before.
932                 //printf("**** Algorithm already set before, number %d, %s ****\n",fitAlgo, fRawAnalyzer->GetName());
933                 return;
934         }
935         //Initialize the requested algorithm
936           AliDebug(1,Form("trunk INIT lms? %d, %d\n", fFittingAlgorithm, Algo::kLMS));
937
938         if(fitAlgo != fFittingAlgorithm || !fRawAnalyzer) {
939                 //printf("**** Init Algorithm , number %d ****\n",fitAlgo);
940                 
941                 fFittingAlgorithm = fitAlgo; 
942                 if (fRawAnalyzer) delete fRawAnalyzer;  // delete prev. analyzer if existed.
943                 
944                 if (fitAlgo == Algo::kFastFit) {
945                         fRawAnalyzer = new AliCaloRawAnalyzerFastFit();
946                 }
947                 else if (fitAlgo == Algo::kNeuralNet) {
948                         fRawAnalyzer = new AliCaloRawAnalyzerNN();
949                 }
950                 else if (fitAlgo == Algo::kLMS) {
951                         fRawAnalyzer = new AliCaloRawAnalyzerLMS();
952                 }
953                 else if (fitAlgo == Algo::kPeakFinder) {
954                         fRawAnalyzer = new AliCaloRawAnalyzerPeakFinder();
955                 }
956                 else if (fitAlgo == Algo::kCrude) {
957                         fRawAnalyzer = new AliCaloRawAnalyzerCrude();
958                 }
959                 else {
960                   fRawAnalyzer = new AliCaloRawAnalyzer();
961                   //fRawAnalyzer = 0;
962                 }
963         }
964         
965 }
966
967