]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - EMCAL/AliEMCALQADataMakerRec.cxx
Move setting of cluster in array after the cluster has been modified with distance...
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALQADataMakerRec.cxx
index b56eb3a533e27a4a8ec4eaaa0cb70aab6429e97b..5080fa4750fbfa4d041a7f98714ae9f200e79f5d 100644 (file)
@@ -37,6 +37,7 @@ Also calculate the ratio of amplitude from LED Monitor system (current/Reference
 #include <TLine.h>
 #include <TText.h>
 #include <TProfile.h> 
+#include <TProfile2D.h> 
 #include <TStyle.h>
 // --- Standard library ---
 
@@ -69,6 +70,8 @@ Also calculate the ratio of amplitude from LED Monitor system (current/Reference
 #include "AliCaloRawAnalyzerKStandard.h"
 #include "AliCaloRawAnalyzerPeakFinder.h"
 #include "AliCaloRawAnalyzerCrude.h"
+#include "AliEMCALGeometry.h"
+#include "AliEMCALTriggerSTURawStream.h"
 
 #include "AliCaloRawAnalyzerFactory.h"
 
@@ -82,6 +85,7 @@ AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(fitAlgorithm fitAlgo) :
   fFittingAlgorithm(0),
   fRawAnalyzer(0),
   fRawAnalyzerTRU(0),
+       fGeom(0),
   fSuperModules(10), // FIXME!!! number of SuperModules; 10 for 2011; update default for later runs 
   fFirstPedestalSample(0),
   fLastPedestalSample(3),
@@ -115,6 +119,8 @@ AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(fitAlgorithm fitAlgo) :
   
   fRawAnalyzerTRU->SetFixTau(kTRUE); 
   fRawAnalyzerTRU->SetTau(2.5); // default for TRU shaper
+
+       fGeom = new AliEMCALGeometry("EMCAL_COMPLETEV1", "EMCAL");
 //  for (Int_t sm = 0 ; sm < fSuperModules ; sm++){
 //    fTextSM[sm] = NULL ;
 //  }
@@ -126,6 +132,7 @@ AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qad
   fFittingAlgorithm(0),
   fRawAnalyzer(0),
   fRawAnalyzerTRU(0),
+       fGeom(0),
   fSuperModules(qadm.GetSuperModules()), 
   fFirstPedestalSample(qadm.GetFirstPedestalSample()), 
   fLastPedestalSample(qadm.GetLastPedestalSample()),  
@@ -324,7 +331,9 @@ void AliEMCALQADataMakerRec::InitRaws()
   Int_t nSMSectors = fSuperModules / 2; // 2 SMs per sector
   Int_t nbinsZ = 2*AliEMCALGeoParams::fgkEMCALCols;
   Int_t nbinsPhi = nSMSectors * AliEMCALGeoParams::fgkEMCALRows;
-   
+       
+  Int_t nTRUCols = 2*AliEMCALGeoParams::fgkEMCALTRUCols; //total TRU columns for 2D TRU histos
+  Int_t nTRURows = nSMSectors*AliEMCALGeoParams::fgkEMCALTRUsPerSM*AliEMCALGeoParams::fgkEMCALTRURows; //total TRU rows for 2D TRU histos
    // counter info: number of channels per event (bins are SM index)
   TProfile * h0 = new TProfile("hLowEmcalSupermodules", "Low Gain EMC: # of towers vs SuperMod;SM Id;# of towers",
                               fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
@@ -344,10 +353,10 @@ void AliEMCALQADataMakerRec::InitRaws()
   // how much above pedestal was the max sample?  (bins are towers)
   TProfile * h4 = new TProfile("hLowEmcalRawMaxMinusMin", "Low Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]", 
                               nTot, -0.5, nTot-0.5, profileOption) ;
-  Add2RawsList(h4, kSigLG, expert, image, !saveCorr) ;
+  Add2RawsList(h4, kSigLG, expert, !image, !saveCorr) ;
   TProfile * h5 = new TProfile("hHighEmcalRawMaxMinusMin", "High Gain EMC: Max - Min vs towerId;Tower Id;Max-Min [ADC counts]",
                               nTot, -0.5, nTot-0.5, profileOption) ;
-  Add2RawsList(h5, kSigHG, expert, image, !saveCorr) ;
+  Add2RawsList(h5, kSigHG, expert, !image, !saveCorr) ;
 
   // total counter: channels per event
   TH1I * h6 = new TH1I("hLowNtot", "Low Gain EMC: Total Number of found towers;# of Towers;Counts", 200, 0, nTot) ;
@@ -365,32 +374,6 @@ void AliEMCALQADataMakerRec::InitRaws()
                               nTot, -0.5, nTot-0.5, profileOption) ;
   Add2RawsList(h9, kPedHG, expert, !image, !saveCorr) ;
        
-  //temp 2D amplitude histogram for the current run
-  fHighEmcHistoH2F = new TH2F("h2DHighEC2", "High Gain EMC:Max - Min [ADC counts]", nbinsZ, -0.5 , nbinsZ-0.5, nbinsPhi, -0.5, nbinsPhi-0.5);
-   fHighEmcHistoH2F->SetDirectory(0) ; // this histo must be memory resident
-  //add ratio histograms: to comapre the current run with the reference data 
-  TH2F * h15 = new TH2F("h2DRatioAmp", "High Gain Ratio to Reference:Amplitude_{current run}/Amplitude_{reference run}", nbinsZ, -0.5 , nbinsZ-0.5, 
-                        nbinsPhi, -0.5, nbinsPhi-0.5);
-  //settings for display in amore
-  h15->SetTitle("Amplitude_{current run}/Amplitude_{reference run}"); 
-  h15->SetMaximum(2.0);
-  h15->SetMinimum(0.1);
-  h15->SetOption("COLZ");
-  gStyle->SetOptStat(0);
-  Int_t color[] = {4,3,2} ;
-  gStyle->SetPalette(3,color);
-  h15->GetZaxis()->SetNdivisions(3);
-  h15->UseCurrentStyle();
-  h15->SetDirectory(0);
-  Add2RawsList(h15, k2DRatioAmp, !expert, image, !saveCorr) ;
-
-  TH1F * h16 = new TH1F("hRatioDist", "Amplitude_{current run}/Amplitude_{reference run} ratio distribution", nTot, 0., 2.);
-  h16->SetMinimum(0.1); 
-  h16->SetMaximum(100.);
-  gStyle->SetOptStat(0);
-  h16->UseCurrentStyle();
-  h16->SetDirectory(0);
-  Add2RawsList(h16, kRatioDist, !expert, image, !saveCorr) ;
  
   // now repeat the same for TRU and LEDMon data
   Int_t nTot2x2 = fSuperModules * AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // max number of TRU channels for all SuperModules
@@ -400,34 +383,35 @@ void AliEMCALQADataMakerRec::InitRaws()
                                fSuperModules, -0.5, fSuperModules-0.5, profileOption) ;
   Add2RawsList(hT0, kNsmodTRU, expert, !image, !saveCorr) ;
 
-  // where did max sample occur? (bins are TRU channels)
-  TProfile * hT1 = new TProfile("hTRUEmcalRawtime", "TRU EMC: Time at Max vs 2x2Id;2x2 Id;Time [ticks]", 
-                               nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
-  Add2RawsList(hT1, kTimeTRU, expert, !image, !saveCorr) ;
-
   // how much above pedestal was the max sample?  (bins are TRU channels)
-  TProfile * hT2 = new TProfile("hTRUEmcalRawMaxMinusMin", "TRU EMC: Max - Min vs 2x2Id;2x2 Id;Max-Min [ADC counts]", 
+  TProfile * hT1 = new TProfile("hTRUEmcalRawMaxMinusMin", "TRU EMC: Max - Min vs 2x2Id;2x2 Id;Max-Min [ADC counts]", 
                                nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
-  Add2RawsList(hT2, kSigTRU, expert, !image, !saveCorr) ;
+  Add2RawsList(hT1, kSigTRU, expert, !image, !saveCorr) ;
 
   // total counter: channels per event
-  TH1I * hT3 = new TH1I("hTRUNtot", "TRU EMC: Total Number of found TRU channels;# of TRU Channels;Counts", 200, 0, nTot2x2) ;
-  hT3->Sumw2() ;
-  Add2RawsList(hT3, kNtotTRU, expert, !image, !saveCorr) ;
-
-  // pedestal (bins are TRU channels)
-  TProfile * hT4 = new TProfile("hTRUEmcalRawPed", "TRU EMC: Pedestal vs 2x2Id;2x2 Id;Pedestal [ADC counts]", 
-                               nTot2x2, -0.5, nTot2x2-0.5, profileOption) ;
-  Add2RawsList(hT4, kPedTRU, expert, !image, !saveCorr) ;
+  TH1I * hT2 = new TH1I("hTRUNtot", "TRU EMC: Total Number of found TRU channels;# of TRU Channels;Counts", 200, 0, nTot2x2) ;
+  hT2->Sumw2() ;
+  Add2RawsList(hT2, kNtotTRU, expert, !image, !saveCorr) ;
 
   // L0 trigger hits: # of hits (bins are TRU channels)
-  TH1I * hT5 = new TH1I("hTRUEmcalL0hits", "L0 trigger hits: Total number of 2x2 L0 generated", nTot2x2, -0.5, nTot2x2);
-  hT5->Sumw2();
-  Add2RawsList(hT5, kNL0TRU, expert, !image, !saveCorr);
+  TH2I * hT3 = new TH2I("hTRUEmcalL0hits", "L0 trigger hits: Total number of 2x2 L0 generated",  nTRUCols, -0.5, nTRUCols - 0.5, nTRURows, -0.5, nTRURows-0.5);
+  hT3->SetOption("COLZ");
+  //hT3->Sumw2();
+  Add2RawsList(hT3, kNL0TRU, !expert, image, !saveCorr);
 
   // L0 trigger hits: average time (bins are TRU channels)
-  TProfile * hT6 = new TProfile("hTRUEmcalL0hitsAvgTime", "L0 trigger hits: average time bin", nTot2x2, -0.5, nTot2x2, profileOption); 
-  Add2RawsList(hT6, kTimeL0TRU, expert, !image, !saveCorr);
+  TProfile2D * hT4 = new TProfile2D("hTRUEmcalL0hitsAvgTime", "L0 trigger hits: average time bin", nTRUCols, -0.5, nTRUCols - 0.5, nTRURows, -0.5, nTRURows-0.5, profileOption);
+  hT4->SetOption("COLZ");
+  Add2RawsList(hT4, kTimeL0TRU, !expert, image, !saveCorr);
+
+  // L0 trigger hits: first in the event (bins are TRU channels)
+  TH1I * hT5 = new TH1I("hTRUEmcalL0hitsFirst", "L0 trigger hits: First hit in the event", nTot2x2, -0.5, nTot2x2);
+  hT5->Sumw2();
+  Add2RawsList(hT5, kNL0FirstTRU, expert, !image, !saveCorr);
+       
+  // L0 trigger hits: average time of first hit in the event (bins are TRU channels)
+  TProfile * hT6 = new TProfile("hTRUEmcalL0hitsFirstAvgTime", "L0 trigger hits: average time of first hit", nTot2x2, -0.5, nTot2x2, profileOption); 
+  Add2RawsList(hT6, kTimeL0FirstTRU, expert, !image, !saveCorr);
 
   // and also LED Mon..
   // LEDMon has both high and low gain channels, just as regular FEE/towers
@@ -473,6 +457,33 @@ void AliEMCALQADataMakerRec::InitRaws()
                               nTotLEDMon, -0.5, nTotLEDMon-0.5, profileOption) ;
   Add2RawsList(hL9, kPedHGLEDMon, expert, !image, !saveCorr) ;
   
+  //temp 2D amplitude histogram for the current run
+  fHighEmcHistoH2F = new TH2F("h2DHighEC2", "High Gain EMC:Max - Min [ADC counts]", nbinsZ, -0.5 , nbinsZ-0.5, nbinsPhi, -0.5, nbinsPhi-0.5);
+   fHighEmcHistoH2F->SetDirectory(0) ; // this histo must be memory resident
+  //add ratio histograms: to comapre the current run with the reference data 
+  TH2F * h15 = new TH2F("h2DRatioAmp", "High Gain Ratio to Reference:Amplitude_{current run}/Amplitude_{reference run}", nbinsZ, -0.5 , nbinsZ-0.5, 
+                        nbinsPhi, -0.5, nbinsPhi-0.5);
+  //settings for display in amore
+  h15->SetTitle("Amplitude_{current run}/Amplitude_{reference run}"); 
+  h15->SetMaximum(2.0);
+  h15->SetMinimum(0.1);
+  h15->SetOption("COLZ");
+  gStyle->SetOptStat(0);
+  Int_t color[] = {4,3,2} ;
+  gStyle->SetPalette(3,color);
+  h15->GetZaxis()->SetNdivisions(3);
+  h15->UseCurrentStyle();
+  h15->SetDirectory(0);
+  Add2RawsList(h15, k2DRatioAmp, !expert, image, !saveCorr) ;
+
+  TH1F * h16 = new TH1F("hRatioDist", "Amplitude_{current run}/Amplitude_{reference run} ratio distribution", nTot, 0., 2.);
+  h16->SetMinimum(0.1); 
+  h16->SetMaximum(100.);
+  gStyle->SetOptStat(0);
+  h16->UseCurrentStyle();
+  h16->SetDirectory(0);
+  Add2RawsList(h16, kRatioDist, !expert, image, !saveCorr) ;
+
   //add two histograms for shifter from the LED monitor system: comapre LED monitor with the reference run
   //to be used for decision whether we need to change reference data
   TH1F * hL10 = new TH1F("hMaxMinusMinLEDMonRatio", "LEDMon amplitude, Ratio to reference run", nTotLEDMon, -0.5, nTotLEDMon-0.5) ;
@@ -494,6 +505,36 @@ void AliEMCALQADataMakerRec::InitRaws()
   Add2RawsList(hL11, kLEDMonRatioDist, !expert, image, !saveCorr) ;
   
   GetCalibRefFromOCDB();   
+
+
+       //STU histgrams
+
+ //histos
+ Int_t nSTUCols = AliEMCALGeoParams::fgkEMCALSTUCols;
+ Int_t nSTURows = AliEMCALGeoParams::fgkEMCALSTURows;
+//             kAmpL1, kGL1, kJL1,
+//             kGL1V0, kJL1V0, kSTUTRU  
+       
+ TProfile2D *hS0 = new TProfile2D("hL1Amp", "Mean STU signal per Row and Column", nSTUCols, -0.5, nSTUCols-0.5, nSTURows, -0.5, nSTURows-0.5);
+ Add2RawsList(hS0, kAmpL1, expert, !image, !saveCorr) ;
+       
+ TH2F *hS1 = new TH2F("hL1Gamma", "L1 Gamma patch position (FastOR top-left)", nSTUCols, -0.50, nSTUCols-0.5, nSTURows, -0.5, nSTURows-0.5);
+ Add2RawsList(hS1, kGL1, expert, image, !saveCorr) ;
+       
+ TH2F *hS2 = new TH2F("hL1Jet", "L1 Jet patch position (FastOR top-left)", 12, -0.5, nSTUCols-0.5, 16, 0, nSTURows-0.5);
+ Add2RawsList(hS2, kJL1, expert, image, !saveCorr) ;
+       
+ TH2I *hS3 = new TH2I("hL1GV0", "L1 Gamma patch amplitude versus V0 signal", 500, 0, 50000, 1500, 0, 1500);
+ Add2RawsList(hS3, kGL1V0, expert, !image, !saveCorr) ;
+       
+ TH2I *hS4 = new TH2I("hL1JV0", "L1 Jet patch amplitude versus V0 signal", 500, 0, 50000, 1000, 0, 1000);
+ Add2RawsList(hS4, kJL1V0, expert, !image, !saveCorr) ;
+
+ TH1I *hS5 = new TH1I("hFrameR","Link between TRU and STU", 32, 0, 32);
+ Add2RawsList(hS5, kSTUTRU, expert, !image, !saveCorr) ;
+
+
+
   //
   ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
 }
@@ -547,14 +588,15 @@ void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
 
   AliRecoParam::EventSpecie_t saveSpecie = fEventSpecie ;
   if (rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent) { 
-    SetEventSpecie(AliRecoParam::kCalib) ;
+    SetEventSpecie(AliRecoParam::kCalib) ;     
   }
-  
+
   const Int_t nTowersPerSM = AliEMCALGeoParams::fgkEMCALRows * AliEMCALGeoParams::fgkEMCALCols; // number of towers in a SuperModule; 24x48
   const Int_t nRows        = AliEMCALGeoParams::fgkEMCALRows; // number of rows per SuperModule
   const Int_t nStripsPerSM = AliEMCALGeoParams::fgkEMCALLEDRefs; // number of strips per SuperModule
   const Int_t n2x2PerSM    = AliEMCALGeoParams::fgkEMCALTRUsPerSM * AliEMCALGeoParams::fgkEMCAL2x2PerTRU; // number of TRU 2x2's per SuperModule
   const Int_t n2x2PerTRU   = AliEMCALGeoParams::fgkEMCAL2x2PerTRU;
+  const Int_t nTot2x2     = fSuperModules * n2x2PerSM; // total TRU channel
 
   // SM counters; decl. should be safe, assuming we don't get more than expected SuperModules..
   Int_t nTotalSMLG[AliEMCALGeoParams::fgkEMCALModules]       = {0};
@@ -564,13 +606,20 @@ void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
   Int_t nTotalSMHGLEDMon[AliEMCALGeoParams::fgkEMCALModules] = {0};
 
   const Int_t nTRUL0ChannelBits = 10; // used for L0 trigger bits checks
+       int firstL0TimeBin = 999;
+  int triggers[nTot2x2][24]; //auxiliary array for L0 trigger - TODO remove hardcoded 24
+       memset(triggers, 0, sizeof(int) * 24 * nTot2x2);
+
   Int_t iSM = 0; // SuperModule index 
   // start loop over input stream  
   while (in.NextDDL()) {
     Int_t iRCU = in.GetDDLNumber() % 2; // RCU0 or RCU1, within SuperModule
+    Int_t iDDL = in.GetDDLNumber();
     fRawAnalyzer->SetIsZeroSuppressed( in.GetZeroSupp() ); 
-
+    
     while (in.NextChannel()) {
+      Int_t iBranch = in.GetBranch();
+      
       iSM = in.GetModule(); // SuperModule
       //prInt_tf("iSM %d DDL %d", iSM, in.GetDDLNumber()); 
       if (iSM>=0 && iSM<fSuperModules) { // valid module reading
@@ -603,7 +652,7 @@ void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
            time = fitResults.GetTof(); 
            firstPedSample = fFirstPedestalSampleTRU;
            lastPedSample  = fLastPedestalSampleTRU;
-           if (in.GetColumn() > n2x2PerTRU) {
+           if (in.GetColumn() >= n2x2PerTRU) {
              isTRUL0IdData = true;
            }
          }
@@ -650,15 +699,21 @@ void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
                  if(iTRUIdInSM < n2x2PerTRU) {
                    Int_t iTRUAbsId = iTRUIdInSM + n2x2PerTRU * iTRUId;
                    // Fill the histograms
-                   FillRawsData(kNL0TRU,iTRUAbsId);
-                   FillRawsData(kTimeL0TRU,iTRUAbsId, startBin);
+                   Int_t globTRUCol, globTRURow;
+                   GetTruChannelPosition(globTRURow, globTRUCol, iSM, iDDL, iBranch, iTRUIdInSM );
+                   
+                   FillRawsData(kNL0TRU, globTRUCol, globTRURow);
+                   FillRawsData(kTimeL0TRU, globTRUCol, globTRURow, startBin);
+                   triggers[iTRUAbsId][startBin] = 1;
+                   
+                   if((int)startBin < firstL0TimeBin) firstL0TimeBin = startBin;
                  }
                }
              }
              startBin--;
            } // i      
          } // TRU L0 Id data                   
-
+         
          // fill histograms
          if ( in.IsLowGain() || in.IsHighGain() ) { // regular towers
            Int_t towerId = iSM*nTowersPerSM + in.GetColumn()*nRows + in.GetRow();
@@ -696,13 +751,13 @@ void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
            nTotalSMTRU[iSM]++; 
            if ( (amp > fMinSignalTRU) && (amp < fMaxSignalTRU) ) { 
              FillRawsData(kSigTRU,iTRU2x2Id, amp);
-             FillRawsData(kTimeTRU,iTRU2x2Id, time);
-           }
-           if (nPed > 0) {
-             for (Int_t i=0; i<nPed; i++) {
-               FillRawsData(kPedTRU,iTRU2x2Id, pedSamples[i]);
-             }
+             //FillRawsData(kTimeTRU,iTRU2x2Id, time);
            }
+           //if (nPed > 0) {
+             //for (Int_t i=0; i<nPed; i++) {
+               //FillRawsData(kPedTRU,iTRU2x2Id, pedSamples[i]);
+             //}
+           //}
          }
          // LED Mon
          else if ( in.IsLEDMonData() ) {
@@ -743,7 +798,16 @@ void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
     }// end while over channel 
    
   }//end while over DDL's, of input stream 
-
+  //filling some L0 trigger histos
+  if( firstL0TimeBin < 999 ){
+    for(Int_t i = 0; i < nTot2x2; i++) {       
+      if( triggers[i][firstL0TimeBin] > 0 ) {
+       //histo->Fill(i,j);
+       FillRawsData(kNL0FirstTRU, i);
+       FillRawsData(kTimeL0FirstTRU, i, firstL0TimeBin);
+      }
+    }
+  }
   
   //calculate the ratio of the amplitude and fill the histograms, only if the events type is Calib
   // RS: operation on the group of histos kSigHG,k2DRatioAmp,kRatioDist,kLEDMonRatio,kLEDMonRatio,kSigLGLEDMon
@@ -837,9 +901,11 @@ void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
   IncEvCountCycleESDs();
   IncEvCountTotalESDs();
   SetEventSpecie(saveSpecie) ; 
+  
+       MakeRawsSTU(rawReader);
+
   // just in case the next rawreader consumer forgets to reset; let's do it here again..
   rawReader->Reset() ;
-  
   return;
 }
 
@@ -1009,3 +1075,191 @@ void AliEMCALQADataMakerRec::ConvertProfile2H(TProfile * p, TH2 * histo)
     histo->SetBinContent(col2d+1, row2d+1, binContent);
   }
 } 
+//____________________________________________________________________________ 
+void AliEMCALQADataMakerRec::GetTruChannelPosition( Int_t &globRow, Int_t &globColumn, Int_t module, Int_t ddl, Int_t branch, Int_t column )
+{
+  Int_t mrow;
+  Int_t mcol;
+  Int_t trow;
+  Int_t tcol;
+  Int_t drow;
+  Int_t rcu;
+  // RCU 0 or 1
+  rcu = ddl % 2;
+
+  // 12 rows of 2x2s in a module (3 TRUs by 4 rows)
+  mrow = (module/2) * 12;
+  // 24 columns per module, odd module numbers increased by 24
+  mcol = (module%2) * 24;
+
+  // position within TRU coordinates
+  tcol = column / 4;
+  trow = column % 4;
+
+  //.combine
+  if( module%2 == 0 ){   // A side
+    // mirror rows
+    trow = 3 - trow;
+
+    // TRU in module row addition
+    drow = (rcu*branch+rcu) * 4;
+
+  }
+  else{   // C side
+    // mirror columns
+    tcol = 23 - tcol;
+
+    // TRU in module row addition
+    drow = (2 - (rcu*branch+rcu)) * 4;
+  }
+
+  // output global row/collumn position (0,0 = SMA0, phi = 0, |eta| = max)
+  globRow = mrow + drow + trow;
+  globColumn = mcol + tcol;
+       return;
+
+}
+//____________________________________________________________________________ 
+void AliEMCALQADataMakerRec::MakeRawsSTU(AliRawReader* rawReader){
+
+       AliEMCALTriggerSTURawStream* inSTU = new AliEMCALTriggerSTURawStream(rawReader);
+       
+       rawReader->Reset();
+       rawReader->Select("EMCAL", 44);
+
+       
+       //L1 segmentation
+       Int_t sizeL1gsubr = 1;
+       Int_t sizeL1gpatch = 2; 
+       Int_t sizeL1jsubr = 4; 
+
+     Int_t EMCALtrig[AliEMCALGeoParams::fgkEMCALSTUCols][AliEMCALGeoParams::fgkEMCALSTURows];
+
+               memset(EMCALtrig, 0, sizeof(int) * AliEMCALGeoParams::fgkEMCALSTUCols * AliEMCALGeoParams::fgkEMCALSTURows);
+               
+
+       
+               
+     if (inSTU->ReadPayLoad()) 
+       {
+                       
+         //Fw version (use in case of change in L1 jet 
+         Int_t fw = inSTU->GetFwVersion();
+         Int_t sizeL1jpatch = 2+(fw >> 16);
+
+         //To check link
+       
+         Int_t mask = inSTU->GetFrameReceived();
+
+
+         for (int i = 0; i < 32; i++)
+           {
+              if ((mask >> i) & 0x1) FillRawsData(kSTUTRU, i);
+           }
+
+
+         //V0 signal in STU
+         Int_t V0Sig = inSTU->GetV0A()+inSTU->GetV0C();
+                       
+         //FastOR amplitude receive from TRU
+         for (Int_t i = 0; i < 32; i++)
+           {
+             UInt_t adc[96];
+             for (Int_t j = 0; j < 96; j++) adc[j] = 0;
+                               
+             inSTU->GetADC(i, adc);
+                               
+             Int_t iTRU = fGeom->GetTRUIndexFromSTUIndex(i);
+                               
+             for (Int_t j = 0; j < 96; j++)
+               {
+                 Int_t idx;
+                 fGeom->GetAbsFastORIndexFromTRU(iTRU, j, idx);
+                               
+                 Int_t px, py;
+                 fGeom->GetPositionInEMCALFromAbsFastORIndex(idx, px, py);
+                                       
+                 EMCALtrig[px][py] = adc[j];
+               }
+           }
+                       
+         //L1 Gamma patches
+         Int_t iTRU_STU, x, y;
+         for (Int_t i = 0; i < inSTU->GetNL1GammaPatch(); i++)
+           {
+             if (inSTU->GetL1GammaPatch(i, iTRU_STU, x, y)) // col (0..23), row (0..3)
+               {
+                 Int_t iTRU;
+                 iTRU = fGeom->GetTRUIndexFromSTUIndex(iTRU_STU);
+                                       
+                 Int_t etaG = 23-x, phiG = y + 4 * int(iTRU/2); //position in EMCal
+                                       
+                 if (iTRU%2) etaG += 24; //C-side
+                                       
+                 etaG = etaG - sizeL1gsubr * sizeL1gpatch + 1;
+                               
+                 //Position of patch L1G (bottom-left FastOR of the patch)
+                       FillRawsData(kGL1, etaG, phiG);
+                                       
+                 //loop to sum amplitude of FOR in the gamma patch
+                 Int_t L1G_PatchAmp = 0;
+                 for (Int_t L1Gx = 0; L1Gx < sizeL1gpatch; L1Gx ++)
+                   {
+                     for (Int_t L1Gy = 0; L1Gy < sizeL1gpatch; L1Gy ++)
+                       {
+                         if (etaG+L1Gx < 48 && phiG+L1Gy < 64) L1G_PatchAmp += EMCALtrig[etaG+L1Gx][phiG+L1Gy];
+                         //cout << EMCALtrig[etaG+L1Gx][phiG+L1Gy] << endl;
+                       }
+                   }
+                       
+                 //if (L1G_PatchAmp > 500) cout << "L1G amp =" << L1G_PatchAmp << endl;
+                       FillRawsData(kGL1V0, V0Sig, L1G_PatchAmp);
+                                       
+               }
+           }
+               
+                       
+         //L1 Jet patches
+         for (Int_t i = 0; i < inSTU->GetNL1JetPatch(); i++)
+           {
+             if (inSTU->GetL1JetPatch(i, x, y)) // col (0,15), row (0,11)
+               {
+                               
+                 Int_t etaJ = sizeL1jsubr * (11-y-sizeL1jpatch + 1);
+                 Int_t phiJ = sizeL1jsubr * (15-x-sizeL1jpatch + 1);
+                                       
+                 //position of patch L1J (FOR bottom-left)
+                       FillRawsData(kJL1, etaJ, phiJ);
+                                       
+                 //loop the sum aplitude of FOR in the jet patch
+                 Int_t L1J_PatchAmp = 0;
+                 for (Int_t L1Jx = 0; L1Jx < sizeL1jpatch*4; L1Jx ++)
+                   {
+                     for (Int_t L1Jy = 0; L1Jy < sizeL1jpatch*4; L1Jy ++)
+                       {
+                         if (etaJ+L1Jx < 48 && phiJ+L1Jy < 64) L1J_PatchAmp += EMCALtrig[etaJ+L1Jx][phiJ+L1Jy];
+                       }
+                   }
+               
+                 //cout << "L1J amp =" << L1J_PatchAmp << endl;
+                       FillRawsData(kJL1V0, V0Sig, L1J_PatchAmp);
+                                       
+               }
+           }
+
+       }
+               
+     //Fill FOR amplitude histo
+     for (Int_t i = 0; i < 48; i++)
+       {
+         for (Int_t j = 0; j < 60; j++)
+           {
+             if (EMCALtrig[i][j] != 0) FillRawsData(kAmpL1, i, j, EMCALtrig[i][j]);
+           }
+       }
+               
+       delete inSTU;
+       return;
+}
+
+