]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
silvermy@ornl.gov - adding AliCaloCalibSignal, from Josh H (UT) and David S (ORNL)
authordsilverm <dsilverm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 30 May 2008 13:42:34 +0000 (13:42 +0000)
committerdsilverm <dsilverm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 30 May 2008 13:42:34 +0000 (13:42 +0000)
EMCAL/AliCaloCalibSignal.cxx [new file with mode: 0644]
EMCAL/AliCaloCalibSignal.h [new file with mode: 0644]
EMCAL/EMCALLinkDefbase.h
EMCAL/libEMCALbase.pkg

diff --git a/EMCAL/AliCaloCalibSignal.cxx b/EMCAL/AliCaloCalibSignal.cxx
new file mode 100644 (file)
index 0000000..ab0a1b0
--- /dev/null
@@ -0,0 +1,520 @@
+/**************************************************************************
+ * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+/* $Id: AliCaloCalibSignal.cxx $ */
+
+//________________________________________________________________________
+//
+// A help class for monitoring and calibration tools: MOOD, AMORE etc.,
+// It can be created and used a la (ctor):
+/*
+  //Create the object for making the histograms
+  fSignals = new AliCaloCalibSignal( fDetType );
+  // AliCaloCalibSignal knows how many modules we have for PHOS or EMCAL
+  fNumModules = fSignals->GetModules();
+*/
+// fed an event:
+//  fSignals->ProcessEvent(fCaloRawStream,fRawEventHeaderBase);
+// asked to draw graphs or profiles:
+//  fSignals->GetGraphAmpVsTimeHighGain(module,column,row)->Draw("ap");
+// or
+//  fSignals->GetProfAmpVsTimeHighGain(module,column,row)->Draw();
+// etc.
+//________________________________________________________________________
+
+#include "TFile.h"
+
+#include "AliRawEventHeaderBase.h"
+#include "AliCaloRawStream.h"
+
+//The include file
+#include "AliCaloCalibSignal.h"
+
+ClassImp(AliCaloCalibSignal)
+
+using namespace std;
+
+// ctor; initialize everything in order to avoid compiler warnings
+// put some reasonable defaults
+AliCaloCalibSignal::AliCaloCalibSignal(kDetType detectorType) :  
+  TObject(),
+  fDetType(kNone),
+  fColumns(0),
+  fRows(0),
+  fModules(0),
+  fRunNumber(-1),
+  fStartTime(0),
+  fAmpCut(50),
+  fReqFractionAboveAmpCutVal(0.8),
+  fReqFractionAboveAmp(kTRUE),
+  fHour(0),
+  fLatestHour(0),
+  fUseAverage(kTRUE),
+  fSecInAverage(1800), 
+  fNEvents(0),
+  fNAcceptedEvents(0)
+{
+  //Default constructor. First we set the detector-type related constants.
+  if (detectorType == kPhos) {
+    fColumns = fgkPhosCols;
+    fRows = fgkPhosRows;
+    fModules = fgkPhosModules;
+  } 
+  else {
+    //We'll just trust the enum to keep everything in line, so that if detectorType
+    //isn't kPhos then it is kEmCal. Note, however, that this is not necessarily the
+    //case, if someone intentionally gives another number
+    fColumns = fgkEmCalCols;
+    fRows = fgkEmCalRows;
+    fModules = fgkEmCalModules;
+  }
+
+  fDetType = detectorType;
+
+  // Set the number of points for each Amp vs. Time graph to 0
+  memset(fNHighGain, 0, sizeof(fNHighGain));
+  memset(fNLowGain, 0, sizeof(fNLowGain));
+
+  CreateGraphs(); // set up the TGraphs
+
+  // init TProfiles to NULL=0 also
+  memset(fProfAmpVsTimeHighGain, 0, sizeof(fProfAmpVsTimeHighGain));
+  memset(fProfAmpVsTimeLowGain, 0, sizeof(fProfAmpVsTimeLowGain));
+} 
+
+// dtor
+//_____________________________________________________________________
+AliCaloCalibSignal::~AliCaloCalibSignal()
+{
+  ClearObjects();
+}
+
+//_____________________________________________________________________
+void AliCaloCalibSignal::ClearObjects()
+{
+  // delete what was created in the ctor (TGraphs), and possible later (TProfiles)
+  for (int i=0; i<fgkMaxTowers; i++) {
+    if ( fGraphAmpVsTimeHighGain[i] ) { delete fGraphAmpVsTimeHighGain[i]; }
+    if ( fGraphAmpVsTimeLowGain[i] ) { delete fGraphAmpVsTimeLowGain[i]; }
+    if ( fProfAmpVsTimeHighGain[i] ) { delete fProfAmpVsTimeHighGain[i]; }
+    if ( fProfAmpVsTimeLowGain[i] ) { delete fProfAmpVsTimeLowGain[i]; }
+  }
+  // set pointers
+  memset(fGraphAmpVsTimeHighGain, 0, sizeof(fGraphAmpVsTimeHighGain));
+  memset(fGraphAmpVsTimeLowGain, 0, sizeof(fGraphAmpVsTimeLowGain));
+  memset(fProfAmpVsTimeHighGain, 0, sizeof(fProfAmpVsTimeHighGain));
+  memset(fProfAmpVsTimeLowGain, 0, sizeof(fProfAmpVsTimeLowGain));
+
+  return;
+}
+
+// copy ctor
+//_____________________________________________________________________
+AliCaloCalibSignal::AliCaloCalibSignal(const AliCaloCalibSignal &sig) :
+  TObject(sig),
+  fDetType(sig.GetDetectorType()),
+  fColumns(sig.GetColumns()),
+  fRows(sig.GetRows()),
+  fModules(sig.GetModules()),
+  fRunNumber(sig.GetRunNumber()),
+  fStartTime(sig.GetStartTime()),
+  fAmpCut(sig.GetAmpCut()),
+  fReqFractionAboveAmpCutVal(sig.GetReqFractionAboveAmpCutVal()),
+  fReqFractionAboveAmp(sig.GetReqFractionAboveAmp()),
+  fHour(sig.GetHour()),
+  fLatestHour(sig.GetLatestHour()),
+  fUseAverage(sig.GetUseAverage()),
+  fSecInAverage(sig.GetSecInAverage()),
+  fNEvents(sig.GetNEvents()),
+  fNAcceptedEvents(sig.GetNAcceptedEvents())
+{
+  // also the TGraph contents
+  AddInfo(&sig);
+}
+
+// assignment operator; use copy ctor to make life easy..
+//_____________________________________________________________________
+AliCaloCalibSignal& AliCaloCalibSignal::operator = (const AliCaloCalibSignal &source)
+{
+  // assignment operator; use copy ctor
+  if (&source == this) return *this;
+
+  new (this) AliCaloCalibSignal(source);
+  return *this;
+}
+
+//_____________________________________________________________________
+void AliCaloCalibSignal::CreateGraphs()
+{
+  //Then, loop for the requested number of modules
+  TString title, name;
+  for (int i = 0; i < fModules; i++) {
+    
+    // Amplitude vs. Time graph for each channel
+    for(int ic=0;ic < fColumns;ic++){
+      for(int ir=0;ir < fRows;ir++){
+       
+       int id = GetTowerNum(i, ic, ir);
+       
+       // high gain graph
+       name = "fGraphAmpVsTimeHighGain_"; name += i;
+       name += "_"; name += ic;
+       name += "_"; name += ir;
+       title = "Amp vs. Time High Gain Mod "; title += i;
+       title += " Col "; title += ic;
+       title += " Row "; title += ir;
+       
+       fGraphAmpVsTimeHighGain[id] = new TGraph();
+       fGraphAmpVsTimeHighGain[id]->SetName(name);
+       fGraphAmpVsTimeHighGain[id]->SetTitle(title);
+       fGraphAmpVsTimeHighGain[id]->SetMarkerStyle(20);
+       
+       // Low Gain
+       name = "fGraphAmpVsTimeLowGain_"; name += i;
+       name += "_"; name += ic;
+       name += "_"; name += ir;
+       title = "Amp vs. Time Low Gain Mod "; title += i;
+       title += " Col "; title += ic;
+       title += " Row "; title += ir;
+       
+       fGraphAmpVsTimeLowGain[id] = new TGraph();
+       fGraphAmpVsTimeLowGain[id]->SetName(name);
+       fGraphAmpVsTimeLowGain[id]->SetTitle(title);
+       fGraphAmpVsTimeLowGain[id]->SetMarkerStyle(20);
+       
+      }
+    }
+
+  }//end for nModules 
+}
+
+//_____________________________________________________________________
+void AliCaloCalibSignal::Reset()
+{
+  Zero(); // set all counters to 0
+  ClearObjects(); // delete previous TGraphs and TProfiles
+  CreateGraphs(); // and create some new ones
+  return;
+}
+
+//_____________________________________________________________________
+void AliCaloCalibSignal::Zero()
+{
+  // set all counters to 0; not cuts etc.though
+  fHour = 0;
+  fLatestHour = 0;
+  fNEvents = 0;
+  fNAcceptedEvents = 0;
+  return;
+}
+
+//_____________________________________________________________________
+Bool_t AliCaloCalibSignal::CheckFractionAboveAmp(int *AmpVal, int nTotChan)
+{
+  int nAbove = 0;
+    
+  int TowerNum = 0;
+  for (int i = 0; i<fModules; i++) {
+    for (int j = 0; j<fColumns; j++) {
+      for (int k = 0; k<fRows; k++) {
+       TowerNum = GetTowerNum(i,j,k);
+       if (AmpVal[TowerNum] > fAmpCut) { 
+         nAbove++;
+       }
+      }
+    }
+  }
+  
+  double fraction = (1.0*nAbove) / nTotChan;
+  
+  if (fraction > fReqFractionAboveAmpCutVal) {  
+    return true;
+  }
+  else return false;
+}
+
+//_____________________________________________________________________
+Bool_t AliCaloCalibSignal::AddInfo(const AliCaloCalibSignal *sig)
+{
+  // just do this for the basic graphs/profiles that get filled in ProcessEvent
+  // may not have data for all channels, but let's just Add everything..
+  for (int i = 0; i < fModules; i++) {
+    for (int j = 0; j < fColumns; j++) {
+      for (int k = 0; k < fRows; k++) {
+       
+       int id = GetTowerNum(i,j,k);
+       if(fUseAverage){
+         GetProfAmpVsTimeHighGain(id)->Add(sig->GetProfAmpVsTimeHighGain(id));
+         GetProfAmpVsTimeLowGain(id)->Add(sig->GetProfAmpVsTimeLowGain(id));
+       }
+       else{     
+         //DS
+//       sig->GetGraphAmpVsTimeHighGain(i,j,k);
+//       sig->GetGraphAmpVsTimeLowGain(i,j,k);
+       }
+
+      }//end for nModules 
+    }//end for nColumns
+  }//end for nRows
+
+  return kTRUE;//We succesfully added info from the supplied object
+}
+
+
+//_____________________________________________________________________
+Bool_t AliCaloCalibSignal::ProcessEvent(AliCaloRawStream *in, AliRawEventHeaderBase *aliHeader)
+{ 
+  // Method to process=analyze one event in the data stream
+  if (!in) return kFALSE; //Return right away if there's a null pointer
+  
+  fNEvents++; // one more event
+
+  // PHOS has more towers than EMCAL, so use PHOS numbers to set array sizes
+  int AmpValHighGain[fgkMaxTowers];
+  int AmpValLowGain[fgkMaxTowers];
+
+  memset(AmpValHighGain, 0, sizeof(AmpValHighGain));
+  memset(AmpValLowGain, 0, sizeof(AmpValLowGain));
+
+  int sample, i = 0; //The sample temp, and the sample number in current event.
+  int max = fgkSampleMin, min = fgkSampleMax;//Use these for picking the signal
+  int gain = 0;
+  
+  // Number of Low and High gain channels for this event:
+  int nLowChan = 0; 
+  int nHighChan = 0; 
+
+  int TowerNum = 0; // array index for TGraphs etc.
+
+  // loop first to get the fraction of channels with amplitudes above cut
+  while (in->Next()) {
+    sample = in->GetSignal(); //Get the adc signal
+    if (sample < min) min = sample;
+    if (sample > max) max = sample;
+    i++;
+    if ( i >= in->GetTimeLength()) {
+      //If we're here then we're done with this tower
+      gain = 1 - in->IsLowGain();
+      
+      int arrayPos = in->GetModule(); //The modules are numbered starting from 0
+      if (arrayPos >= fModules) {
+       //TODO: return an error message, if appopriate (perhaps if debug>0?)
+       return kFALSE;
+      } 
+      
+      //Debug
+      if (arrayPos < 0 || arrayPos >= fModules) {
+       printf("Oh no: arrayPos = %i.\n", arrayPos); 
+      }
+
+      // get tower number for AmpVal array
+      TowerNum = GetTowerNum(arrayPos, in->GetColumn(), in->GetRow()); 
+
+      if (gain == 0) {
+       // fill amplitude into the array           
+        AmpValLowGain[TowerNum] = max - min;
+       nLowChan++;
+      } 
+      else if (gain==1) {//fill the high gain ones
+       // fill amplitude into the array
+       AmpValHighGain[TowerNum] = max - min;
+       nHighChan++;
+      }//end if gain
+
+      
+      max = fgkSampleMin; min = fgkSampleMax;
+      i = 0;
+      
+    }//End if end of tower
+   
+  }//end while, of stream
+  
+  // now check if it was a led event, only use high gain (that should be sufficient)
+  if (fReqFractionAboveAmp) {
+    bool ok = CheckFractionAboveAmp(AmpValHighGain, nHighChan);
+    if (!ok) return false; // skip event
+  }
+
+  fNAcceptedEvents++; // one more event accepted
+
+  if (fStartTime == 0) { // if start-timestamp wasn't set,we'll pick it up from the first event we encounter
+    fStartTime = aliHeader->Get("Timestamp");
+  }
+
+  fHour = (aliHeader->Get("Timestamp")-fStartTime)/(double)fgkNumSecInHr;
+  if (fLatestHour < fHour) {
+    fLatestHour = fHour; 
+  }
+  
+  // it is a led event, now fill graphs (maybe profiles later)
+  for(int i=0;i<fModules;i++){
+    for(int j=0;j<fColumns;j++){
+      for(int k=0;k<fRows;k++){
+       
+       TowerNum = GetTowerNum(i, j, k); 
+
+       if(AmpValHighGain[TowerNum]) {
+         fGraphAmpVsTimeHighGain[TowerNum]->SetPoint(fNHighGain[TowerNum]++,fHour,AmpValHighGain[TowerNum]);
+       }
+       if(AmpValLowGain[TowerNum]) {
+         fGraphAmpVsTimeLowGain[TowerNum]->SetPoint(fNLowGain[TowerNum]++,fHour,AmpValLowGain[TowerNum]);
+       }
+      }
+    }
+  }
+  
+  return kTRUE;
+}
+
+//_____________________________________________________________________
+void AliCaloCalibSignal::CreateProfile(int imod, int ic, int ir, int towerId, int gain,
+                                      int nbins, double min, double max)
+{ //! create/setup a TProfile
+  TString title, name;   
+  if (gain == 0) { 
+    name = "fProfAmpVsTimeLowGain_";   
+    title = "Amp vs. Time Low Gain Mod "; 
+  } 
+  else if (gain == 1) { 
+    name = "fProfAmpVsTimeHighGain_"; 
+    title = "Amp vs. Time High Gain Mod "; 
+  } 
+  name += imod;
+  name += "_"; name += ic;
+  name += "_"; name += ir;
+  title += imod;
+  title += " Col "; title += ic;
+  title += " Row "; title += ir;
+           
+  // use "s" option for RMS
+  if (gain==0) { 
+    fProfAmpVsTimeLowGain[towerId] = new TProfile(name,title, nbins, min, max,"s");
+  }
+  else if (gain==1) {
+    fProfAmpVsTimeHighGain[towerId] = new TProfile(name,title, nbins, min, max,"s");
+  }
+
+  return;
+}
+//_____________________________________________________________________
+Bool_t AliCaloCalibSignal::Save(TString fileName, Bool_t saveEmptyGraphs)
+{
+  //Saves all the histograms (or profiles, to be accurate) to the designated file
+  
+  TFile destFile(fileName, "recreate");
+  
+  if (destFile.IsZombie()) {
+    return kFALSE;
+  }
+  
+  destFile.cd();
+
+  // setup variables for the TProfile plot
+  int numProfBins = 0;
+  double timeMin = 0;
+  double timeMax = 0;
+  if (fUseAverage) {
+    if (fSecInAverage > 0) {
+      numProfBins = (int)( (fLatestHour*fgkNumSecInHr)/fSecInAverage + 1 ); // round-off
+    }
+    numProfBins += 2; // add extra buffer : first and last
+    double binSize = 1.0*fSecInAverage / fgkNumSecInHr;
+    timeMin = - binSize;
+    timeMax = timeMin + numProfBins*binSize;
+  }
+
+  int numGraphPoints= 0;
+  int TowerNum = 0;    
+  for (int i = 0; i < fModules; i++) {
+    
+    for(int ic=0;ic < fColumns;ic++){
+      for(int ir=0;ir < fRows;ir++){
+
+       TowerNum = GetTowerNum(i, ic, ir); 
+
+       // 1st: high gain
+       numGraphPoints= fGraphAmpVsTimeHighGain[TowerNum]->GetN();
+       if( numGraphPoints>0 || saveEmptyGraphs) {
+         
+         // average the graphs points over time if requested and put them in a profile plot
+         if(fUseAverage && numGraphPoints>0) {
+           
+           // get the values
+           double *graphX = fGraphAmpVsTimeHighGain[TowerNum]->GetX();
+           double *graphY = fGraphAmpVsTimeHighGain[TowerNum]->GetY();
+
+           // create the TProfile: 1 is for High gain              
+           CreateProfile(i, ic, ir, TowerNum, 1,
+                         numProfBins, timeMin, timeMax);
+
+           // loop over graph points and fill profile
+           for(int ip=0; ip < numGraphPoints; ip++){
+             fProfAmpVsTimeHighGain[TowerNum]->Fill(graphX[ip],graphY[ip]);
+           }
+           
+           fProfAmpVsTimeHighGain[TowerNum]->GetXaxis()->SetTitle("Hours");
+           fProfAmpVsTimeHighGain[TowerNum]->GetYaxis()->SetTitle("MaxAmplitude - Pedestal");
+           fProfAmpVsTimeHighGain[TowerNum]->Write();
+
+         }
+          else{
+            //otherwise, just save the graphs and forget the profiling
+            fGraphAmpVsTimeHighGain[TowerNum]->GetXaxis()->SetTitle("Hours");
+            fGraphAmpVsTimeHighGain[TowerNum]->GetYaxis()->SetTitle("MaxAmplitude - Pedestal");
+            fGraphAmpVsTimeHighGain[TowerNum]->Write();
+          }
+         
+       } // low gain graph info should be saved in one form or another
+       
+       // 2nd: now go to the low gain case
+       numGraphPoints= fGraphAmpVsTimeLowGain[TowerNum]->GetN();
+       if( numGraphPoints>0 || saveEmptyGraphs) {
+         
+         // average the graphs points over time if requested and put them in a profile plot
+         if(fUseAverage && numGraphPoints>0) {
+           
+           double *graphX = fGraphAmpVsTimeLowGain[TowerNum]->GetX();
+           double *graphY = fGraphAmpVsTimeLowGain[TowerNum]->GetY();
+           
+           // create the TProfile: 0 is for Low gain       
+           CreateProfile(i, ic, ir, TowerNum, 0,
+                         numProfBins, timeMin, timeMax);
+
+           // loop over graph points and fill profile
+           for(int ip=0; ip < numGraphPoints; ip++){
+             fProfAmpVsTimeLowGain[TowerNum]->Fill(graphX[ip],graphY[ip]);
+           }
+           
+           fProfAmpVsTimeLowGain[TowerNum]->GetXaxis()->SetTitle("Hours");
+           fProfAmpVsTimeLowGain[TowerNum]->GetYaxis()->SetTitle("MaxAmplitude - Pedestal");
+           fProfAmpVsTimeLowGain[TowerNum]->Write();
+
+         }
+         else{
+            //otherwise, just save the graphs and forget the profiling
+           fGraphAmpVsTimeLowGain[TowerNum]->GetXaxis()->SetTitle("Hours");
+           fGraphAmpVsTimeLowGain[TowerNum]->GetYaxis()->SetTitle("MaxAmplitude - Pedestal");
+           fGraphAmpVsTimeLowGain[TowerNum]->Write();
+         }
+         
+       } // low gain graph info should be saved in one form or another
+
+      } // fRows
+    } // fColumns
+
+  } // fModules
+  destFile.Close();
+  
+  return kTRUE;
+}
diff --git a/EMCAL/AliCaloCalibSignal.h b/EMCAL/AliCaloCalibSignal.h
new file mode 100644 (file)
index 0000000..578ffda
--- /dev/null
@@ -0,0 +1,188 @@
+#ifndef ALICALOCALIBSIGNAL_H
+#define ALICALOCALIBSIGNAL_H
+
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/* $Id: AliCaloCalibSignal.h  $ */
+
+// \file AliCaloCalibSignal.h
+//   \brief Description:
+//   A help class for monitoring and calibration tools: MOOD, AMORE etc.,
+//   that can process events from a standard AliCaloRawStream,
+//   most usually from LED/pulser runs. It stores signal info as
+//   typical (highest) amplitude vs time in TGraphs (one per channel)
+//   or TProfiles if we decide to just store the averages (and not all points) 
+//   for the detectors (EMCAL and PHOS).
+
+//   \author: Josh Hamblen (UTenn), original version. 
+//   [Consultant: D. Silvermyr (ORNL)]
+//   Partly based on AliCaloCalibPedestal.
+//   
+//   \version $Revision:  $
+//   \date $Date: $
+
+#include "TGraph.h"
+#include "TProfile.h"
+class AliCaloRawStream;
+class AliRawEventHeaderBase;
+
+class AliCaloCalibSignal : public TObject {
+  
+ public:
+
+  enum kDetType {kPhos, kEmCal, kNone};//The detector types
+  
+  AliCaloCalibSignal(kDetType detectorType = kPhos); //ctor
+  virtual ~AliCaloCalibSignal(); //dtor
+
+  // copy ctor, and '=' operator, are not fully tested/debugged yet
+  AliCaloCalibSignal(const AliCaloCalibSignal &sig); // copy ctor
+  AliCaloCalibSignal& operator = (const  AliCaloCalibSignal &source); //!
+  
+  Bool_t ProcessEvent(AliCaloRawStream *in, AliRawEventHeaderBase *aliHeader); // added header for time info
+  Bool_t CheckFractionAboveAmp(int *AmpVal, int nTotChan); // check fraction of signals to check for LED events
+
+  ////////////////////////////
+  //Simple getters
+  // need public access to the TGraphs.. 
+  TGraph * GetGraphAmpVsTimeHighGain(int imod, int icol, int irow) const // Return a pointer to the high gain graph
+    { int towId = GetTowerNum(imod, icol, irow); return fGraphAmpVsTimeHighGain[towId];}; //!
+  TGraph * GetGraphAmpVsTimeLowGain(int imod, int icol, int irow) const // Return a pointer to the low gain graph
+    {int towId = GetTowerNum(imod, icol, irow); return fGraphAmpVsTimeLowGain[towId];};        
+  TGraph * GetGraphAmpVsTimeHighGain(int towId) const // Return a pointer to the high gain graph
+    { return fGraphAmpVsTimeHighGain[towId];}; //!
+  TGraph * GetGraphAmpVsTimeLowGain(int towId) const // Return a pointer to the low gain graph
+    { return fGraphAmpVsTimeLowGain[towId];}; //!      
+
+  // and similarly for the TProfiles
+  TProfile * GetProfAmpVsTimeHighGain(int imod, int icol, int irow) const // Return a pointer to the high gain profile
+    { int towId = GetTowerNum(imod, icol, irow); return fProfAmpVsTimeHighGain[towId];}; //!   
+  TProfile * GetProfAmpVsTimeLowGain(int imod, int icol, int irow) const // Return a pointer to the low gain profile
+    { int towId = GetTowerNum(imod, icol, irow); return fProfAmpVsTimeLowGain[towId];}; //!    
+  TProfile * GetProfAmpVsTimeHighGain(int towId) const // Return a pointer to the high gain profile
+    { return fProfAmpVsTimeHighGain[towId];}; //!      
+  TProfile * GetProfAmpVsTimeLowGain(int towId) const // Return a pointer to the low gain profile
+    { return fProfAmpVsTimeLowGain[towId];}; //!       
+
+  // how many points do we have in each TGraph
+  int GetNHighGain(int imod, int icol, int irow) const //!
+    { int towId = GetTowerNum(imod, icol, irow); return fNHighGain[towId];};   //!
+  int GetNLowGain(int imod, int icol, int irow) const //!
+    { int towId = GetTowerNum(imod, icol, irow); return fNLowGain[towId];};    //!
+  int GetNHighGain(int towId) const { return fNHighGain[towId];};      //!
+  int GetNLowGain(int towId) const { return fNLowGain[towId];};        //!
+
+  // Basic info: getters  
+  kDetType GetDetectorType() const {return fDetType;};//Returns if this is a PHOS or EMCAL object
+  
+  int GetColumns() const {return fColumns;}; //The number of columns per module
+  int GetRows() const {return fRows;}; //The number of rows per module
+  int GetModules() const {return fModules;}; //The number of modules
+  int GetTowerNum(int imod, int icol, int irow) const { return imod*fColumns*fRows + icol*fRows + irow;}; // help index
+
+  // Basic Counters
+  int GetNEvents() const {return fNEvents;};
+  int GetNAcceptedEvents() const {return fNAcceptedEvents;};
+
+  ///////////////////////////////
+  //  Get and Set Cuts
+  // Section for if we should help with the event selection of what is likely LED events
+  void SetAmpCut(double d) { fAmpCut = d; } //!
+  double GetAmpCut() const { return fAmpCut; }; //!
+  void SetReqFractionAboveAmpCutVal(double d) { fReqFractionAboveAmpCutVal = d; } //!
+  double GetReqFractionAboveAmpCutVal() const { return fReqFractionAboveAmpCutVal; }; //!
+  void SetReqFractionAboveAmp(bool b) { fReqFractionAboveAmp = b; } //!
+  double GetReqFractionAboveAmp() const { return fReqFractionAboveAmp; }; //!
+
+  // We may select to only use the averaged info in the TProfiles rather than the
+  // the full in the TGraphs
+  void SetUseAverage(bool b) { fUseAverage = b; } //!
+  double GetUseAverage() const { return fUseAverage; }; //!
+  void SetSecInAverage(int secInAverage) {fSecInAverage = secInAverage;}; // length of the interval that should be used for the average calculation (determines number of bins in TProfile)
+  int GetSecInAverage() const {return fSecInAverage;}; //!
+
+  // Info on time since start of run
+  double GetHour() const { return fHour; }; // time info for current event
+  double GetCurrentHour() const { return fHour; }; // time info for current event (same as GetHour(), just more explicitly named)
+  double GetLatestHour() const { return fLatestHour; }; // the latest time encountered
+  // These times are typically the same, but not necessarily if the events do not come in order 
+  void SetLatestHour(double d) { fLatestHour = d; }; // could be useful when we know the length of the run (i.e. after it is over), e.g. for PreProcessor
+
+  // RunNumbers : setters and getters
+  void SetRunNumber(int runNo) {fRunNumber = runNo;}; //!
+  int GetRunNumber() const {return fRunNumber;};  //!
+  
+  // Start-of-run timestamp : set and get
+  void SetStartTime(int startTime) {fStartTime = startTime;}; //!
+  int GetStartTime() const {return fStartTime;}; //!
+
+  /////////////////////////////
+  //Analysis functions
+  void Reset();//Resets the whole class.
+  Bool_t AddInfo(const AliCaloCalibSignal *sig);//picks up new info from supplied argument  
+
+  //Saving functions
+  Bool_t Save(TString fileName, Bool_t saveEmptyGraphs = kFALSE); //Saves the TGraphs to a .root file
+
+ private:
+
+  void ClearObjects(); // delete old objects and set pointers
+  void Zero(); // set all counters to 0
+  void CreateGraphs(); //! create/setup the TGraphs
+  void CreateProfile(int imod, int icol, int irow, int towerId, int gain,
+                    int nbins, double min, double max); //! create/setup a TProfile
+    
+ private:
+
+  kDetType fDetType; //The detector type for this object
+  int fColumns;        //The number of columns per module
+  int fRows;   //The number of rows per module
+  int fModules;        //The number of modules
+  int fRunNumber; //The run number. Needs to be set by the user.
+  int fStartTime;  // Time of first event
+
+  double fAmpCut; // amplitude cut value
+  double fReqFractionAboveAmpCutVal; // required fraction that should be above cut
+  bool fReqFractionAboveAmp; // flag to select if we should do some event selection based on amplitudes
+
+  double fHour; // fraction of hour since beginning of run, for amp vs. time graphs, for current event
+  double fLatestHour; // largest fraction of hour since beginning of run, for amp vs. time graphs
+  bool fUseAverage; // flag to average graph points into over a time interval
+  int fSecInAverage; // time interval for the graph averaging
+
+  // status counters
+  int fNEvents; // # events processed
+  int fNAcceptedEvents; // # events accepted
+
+  //Constants needed by the class
+  static const int fgkSampleMax = 1023; // highest possible sample value (10-bit = 0x3ff)
+  static const int fgkSampleMin = 0; // lowest possible sample value 
+  
+  static const int fgkPhosRows = 64; // number of rows per module for PHOS
+  static const int fgkPhosCols = 56; // number of columns per module for PHOS
+  static const int fgkPhosModules = 5; // number of modules for PHOS
+  
+  static const int fgkEmCalRows = 24; // number of rows per module for EMCAL
+  static const int fgkEmCalCols = 48; // number of columns per module for EMCAL
+  static const int fgkEmCalModules = 12; // number of modules for EMCAL
+
+  // From numbers above: PHOS has more possible towers (17920) than EMCAL (13824) 
+  // so use PHOS numbers to set max. array sizes
+  static const int fgkMaxTowers = 17920; // fgkPhosModules * fgkPhosCols * fgkPhosRows; 
+  
+  static const int fgkNumSecInHr = 3600;  // number of seconds in an hour, for the fractional hour conversion on the time graph
+  
+  TGraph *fGraphAmpVsTimeHighGain[fgkMaxTowers]; // Amplitude vs. Time Graph for each high gain channel
+  TGraph *fGraphAmpVsTimeLowGain[fgkMaxTowers]; // Amplitude vs. Time Graph for each low gain channel
+  TProfile *fProfAmpVsTimeHighGain[fgkMaxTowers]; // Amplitude vs. Time Profile for each high gain channel
+  TProfile *fProfAmpVsTimeLowGain[fgkMaxTowers]; // Amplitude vs. Time Profile for each low gain channel
+  
+  int fNHighGain[fgkMaxTowers]; // Number of points for each Amp. vs. Time graph
+  int fNLowGain[fgkMaxTowers]; // Number of points for each Amp. vs. Time graph
+  
+  ClassDef(AliCaloCalibSignal,1)
+    
+};
+    
+#endif
index 76f2954391bc62f0585fdfecc1d5dd1228039665..9bbaf7d60bbc7429bae2823515ec59e056571b45 100644 (file)
@@ -21,6 +21,7 @@
 #pragma link C++ class AliEMCALSensorTemp+;
 #pragma link C++ class AliEMCALSensorTempArray+;
 #pragma link C++ class AliCaloCalibPedestal+;
+#pragma link C++ class AliCaloCalibSignal+;
 #pragma link C++ class AliEMCALSurvey+;
 #pragma link C++ class AliEMCALRecParam+;
 #pragma link C++ class AliEMCALQAChecker+;
index c7f3eab1f24a7048000c308a04fc3b13aa520518..6a0b272fed69c705ff50758e321275e9a21df8b5 100644 (file)
@@ -18,6 +18,7 @@ AliEMCALPreprocessor.cxx \
 AliEMCALSensorTemp.cxx \
 AliEMCALSensorTempArray.cxx \
 AliCaloCalibPedestal.cxx \
+AliCaloCalibSignal.cxx \
 AliEMCALSurvey.cxx \
 AliEMCALRecParam.cxx \
 AliEMCALQAChecker.cxx