]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
New VZERO tender supply which can be used in order to correct for the LHC clock phase...
authorcvetan <cvetan@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 6 Sep 2010 13:58:55 +0000 (13:58 +0000)
committercvetan <cvetan@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 6 Sep 2010 13:58:55 +0000 (13:58 +0000)
ANALYSIS/TenderSupplies/AliVZEROTenderSupply.cxx [new file with mode: 0644]
ANALYSIS/TenderSupplies/AliVZEROTenderSupply.h [new file with mode: 0644]
ANALYSIS/TenderSuppliesLinkDef.h
ANALYSIS/libTENDERSupplies.pkg

diff --git a/ANALYSIS/TenderSupplies/AliVZEROTenderSupply.cxx b/ANALYSIS/TenderSupplies/AliVZEROTenderSupply.cxx
new file mode 100644 (file)
index 0000000..f84f3b9
--- /dev/null
@@ -0,0 +1,223 @@
+/**************************************************************************
+ * 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.                  *
+ **************************************************************************/
+
+
+///////////////////////////////////////////////////////////////////////////////
+//                                                                           //
+//  Recalculate VZERO timing and decision using the tender                   //
+//  (in case the LHC phase drift is updated in OCDB)                         //
+//                                                                           //
+///////////////////////////////////////////////////////////////////////////////
+
+#include <TList.h>
+#include <TObjString.h>
+#include <TChain.h>
+#include <TF1.h>
+
+#include <AliLog.h>
+#include <AliESDEvent.h>
+#include <AliESDVZERO.h>
+#include <AliCDBId.h>
+#include <AliCDBManager.h>
+#include <AliCDBEntry.h>
+#include <AliTender.h>
+#include <AliLHCClockPhase.h>
+#include <AliVZEROCalibData.h>
+#include <AliVZEROTriggerMask.h>
+
+#include "AliVZEROTenderSupply.h"
+
+ClassImp(AliVZEROTenderSupply)
+
+AliVZEROTenderSupply::AliVZEROTenderSupply() :
+  AliTenderSupply(),
+  fCalibData(NULL),
+  fTimeSlewing(NULL),
+  fLHCClockPhase(0),
+  fDebug(kFALSE)
+{
+  //
+  // default ctor
+  //
+}
+
+//_____________________________________________________
+AliVZEROTenderSupply::AliVZEROTenderSupply(const char *name, const AliTender *tender) :
+  AliTenderSupply(name,tender),
+  fCalibData(NULL),
+  fTimeSlewing(NULL),
+  fLHCClockPhase(0),
+  fDebug(kFALSE)
+{
+  //
+  // named ctor
+  //
+}
+
+//_____________________________________________________
+void AliVZEROTenderSupply::Init()
+{
+  //
+  // Initialise VZERO tender
+  //
+}
+
+//_____________________________________________________
+void AliVZEROTenderSupply::ProcessEvent()
+{
+  //
+  // Reapply the LHC-clock phase drift
+  //
+
+  AliESDEvent *event=fTender->GetEvent();
+  if (!event) return;
+  
+  //load gain correction if run has changed
+  if (fTender->RunChanged()){
+    if (fDebug) printf("AliVZEROTenderSupply::ProcessEvent - Run Changed (%d)\n",fTender->GetRun());
+    GetPhaseCorrection();
+
+    AliCDBEntry *entryCal = fTender->GetCDBManager()->Get("VZERO/Calib/Data",fTender->GetRun());
+    if (!entryCal) {
+      AliError("No VZERO calibration entry is found");
+      fCalibData = NULL;
+      return;
+    } else {
+      fCalibData = (AliVZEROCalibData*)entryCal->GetObject();
+      if (fDebug) printf("AliVZEROTenderSupply::Used VZERO calibration entry: %s\n",entryCal->GetId().ToString().Data());
+    }
+
+    AliCDBEntry *entrySlew = fTender->GetCDBManager()->Get("VZERO/Calib/TimeSlewing",fTender->GetRun());
+    if (!entrySlew) {
+      AliError("VZERO time slewing function is not found in OCDB !");
+      fTimeSlewing = NULL;
+      return;
+    } else {
+      fTimeSlewing = (TF1*)entrySlew->GetObject();
+      if (fDebug) printf("AliVZEROTenderSupply::Used VZERO time slewing entry: %s\n",entrySlew->GetId().ToString().Data());
+    }
+  }
+
+  if (!fCalibData || !fTimeSlewing) return;
+
+  //
+  // correct VZERO time signals and decision
+  //
+  AliESDVZERO *esdVZERO = event->GetVZEROData();
+  if (!esdVZERO) {
+    AliError("No VZERO object is found inside ESD!");
+    return;
+  }
+  if (!esdVZERO->TestBit(AliESDVZERO::kDecisionFilled)) {
+    AliWarning("VZERO offline trigger decisions were not filled in ESD, the tender supply is disabled");
+    return;
+  }
+
+  Float_t time[64];
+  for(Int_t i = 0; i < 64; ++i) {
+    time[i] = esdVZERO->GetTime(i);
+    time[i] += fLHCClockPhase;
+  }
+  esdVZERO->SetTime(time);
+
+  {
+    AliVZEROTriggerMask triggerMask;
+    triggerMask.FillMasks(esdVZERO, fCalibData, fTimeSlewing);
+  }
+
+}
+
+//_____________________________________________________
+void AliVZEROTenderSupply::GetPhaseCorrection()
+{
+  //
+  // Get Gain splines from OCDB
+  //
+
+  AliInfo("Get LHC-clock phase correction");
+
+  //
+  //find previous entry from the UserInfo
+  //
+  TTree *tree=((TChain*)fTender->GetInputData(0))->GetTree();
+  if (!tree) {
+    AliError("Tree not found in ESDhandler");
+    return;
+  }
+  
+  TList *userInfo=(TList*)tree->GetUserInfo();
+  if (!userInfo) {
+    AliError("No UserInfo found in tree");
+    return;
+  }
+
+  TList *cdbList=(TList*)userInfo->FindObject("cdbList");
+  if (!cdbList) {
+    AliError("No cdbList found in UserInfo");
+    if (AliLog::GetGlobalLogLevel()>=AliLog::kError) userInfo->Print();
+    return;
+  }
+
+  Float_t oldPhase = 0;
+
+  TIter nextCDB(cdbList);
+  TObjString *os=0x0;
+  while ( (os=(TObjString*)nextCDB()) ){
+    if (!(os->GetString().Contains("GRP/Calib/LHCClockPhase"))) continue;
+    AliCDBId *id=AliCDBId::MakeFromString(os->GetString());
+    
+    AliCDBEntry *entry=fTender->GetCDBManager()->Get(*id);
+    if (!entry) {
+      AliError("The previous LHC-clock phase entry is not found");
+      delete id;
+      return;
+    }
+    
+    if (fDebug) printf("AliVZEROTenderSupply::Used old LHC-clock phase entry: %s\n",entry->GetId().ToString().Data());
+    
+    AliLHCClockPhase *phase = (AliLHCClockPhase*)entry->GetObject();
+    if (!phase) {
+      AliError("Phase object is not found in the calibration entry");
+      delete id;
+      return;
+    }
+
+    oldPhase = phase->GetMeanPhase();
+
+    delete id;
+    break;
+  }
+
+  //
+  //new LHC-clock phase entry
+  //
+  Float_t newPhase = 0;
+  AliCDBEntry *entryNew=fTender->GetCDBManager()->Get("GRP/Calib/LHCClockPhase",fTender->GetRun());
+  if (!entryNew) {
+    AliError("No new LHC-clock phase calibration entry is found");
+    return;
+  }
+  if (fDebug) printf("AliVZEROTenderSupply::Used new LHC-clock phase entry: %s\n",entryNew->GetId().ToString().Data());
+  
+  AliLHCClockPhase *phase2 = (AliLHCClockPhase*)entryNew->GetObject();
+  if (!phase2) {
+    AliError("Phase object is not found in the calibration entry");
+    return;
+  }
+
+  newPhase = phase2->GetMeanPhase();
+
+  fLHCClockPhase = newPhase - oldPhase;
+}
diff --git a/ANALYSIS/TenderSupplies/AliVZEROTenderSupply.h b/ANALYSIS/TenderSupplies/AliVZEROTenderSupply.h
new file mode 100644 (file)
index 0000000..c1ccf56
--- /dev/null
@@ -0,0 +1,50 @@
+#ifndef ALIVZEROTENDERSUPPLY_H
+#define ALIVZEROTENDERSUPPLY_H
+
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+////////////////////////////////////////////////////////////////////////
+//                                                                    //
+//  Recalculate VZERO timing and decision using the tender            //
+//  (in case the LHC phase drift is updated in OCDB)                  //
+//                                                                    //
+////////////////////////////////////////////////////////////////////////
+
+
+
+#include <AliTenderSupply.h>
+
+class TF1;
+class AliVZEROCalibData;
+
+class AliVZEROTenderSupply: public AliTenderSupply {
+  
+public:
+  AliVZEROTenderSupply();
+  AliVZEROTenderSupply(const char *name, const AliTender *tender=NULL);
+  
+  virtual ~AliVZEROTenderSupply(){;}
+
+  virtual void              Init();
+  virtual void              ProcessEvent();
+  
+  void GetPhaseCorrection();
+
+  void SetDebug(Bool_t flag) { fDebug = flag; }
+
+private:
+  AliVZEROCalibData* fCalibData;      //! calibration data
+  TF1*               fTimeSlewing;    //! Function for time slewing correction
+  Float_t            fLHCClockPhase;  //! the correction to the LHC-clock phase
+  Bool_t             fDebug;          //  debug on/off
+  
+  AliVZEROTenderSupply(const AliVZEROTenderSupply&c);
+  AliVZEROTenderSupply& operator= (const AliVZEROTenderSupply&c);
+  
+  ClassDef(AliVZEROTenderSupply, 1)  // VZERO tender task
+};
+
+
+#endif
+
index 9060554d2ba931bf040376fff3ee258403d33a14..387133e2794f94d24782b0d04e033ffe229c5bcf 100644 (file)
@@ -14,4 +14,5 @@
 #pragma link C++ class AliTPCTenderSupply+;
 #pragma link C++ class AliTRDTenderSupply+;
 #pragma link C++ class AliVtxTenderSupply+;
+#pragma link C++ class AliVZEROTenderSupply+;
 
index 99ed61ee01b75906b5a20b56a83811df50fea897..22c6b29b026bf457d8dbd7ad1bf883cb9637b249 100644 (file)
@@ -7,17 +7,18 @@ SRCS= TenderSupplies/AliTOFcalibESD.cxx \
       TenderSupplies/AliTOFTenderSupply.cxx \
       TenderSupplies/AliTPCTenderSupply.cxx \
       TenderSupplies/AliTRDTenderSupply.cxx \
-      TenderSupplies/AliVtxTenderSupply.cxx
+      TenderSupplies/AliVtxTenderSupply.cxx \
+      TenderSupplies/AliVZEROTenderSupply.cxx
 
 HDRS= $(SRCS:.cxx=.h)
 
 DHDR= TenderSuppliesLinkDef.h
 
-EINCLUDE:= ANALYSIS ANALYSIS/Tender STEER TOF ANALYSIS/TenderSupplies
+EINCLUDE:= ANALYSIS ANALYSIS/Tender STEER TOF VZERO ANALYSIS/TenderSupplies
 
 ifeq (win32gcc,$(ALICE_TARGET))
 PACKSOFLAGS:= $(SOFLAGS) -L$(ALICE_ROOT)/lib/tgt_$(ALICE_TARGET) -lSTEERBase \
--lESD -lSTEER -lANALYSISalice -lANALYSIS -lCORRFW -lTENDER -lTOFbase
+-lESD -lSTEER -lANALYSISalice -lANALYSIS -lCORRFW -lTENDER -lTOFbase -lVZEROBase -lVZERORec
 endif