]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - TPC/AliTPCcalibDB.cxx
adding the libPHOSUtils dependency for the libAliHLTGlobal
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibDB.cxx
index 49b72e74e497aaa6dcae993558f08639b3bb1b9c..4e737b7c29db6da0868d7cf394352a3cb7ac6222 100644 (file)
@@ -87,6 +87,7 @@
 #include <AliSplineFit.h>
 
 #include "AliTPCcalibDB.h"
+#include "AliTPCcalibDButil.h"
 #include "AliTPCAltroMapping.h"
 #include "AliTPCExB.h"
 
@@ -103,10 +104,12 @@ class AliTPCCalDet;
 
 #include "TFile.h"
 #include "TKey.h"
+#include "TGraphErrors.h"
 
 #include "TObjArray.h"
 #include "TObjString.h"
 #include "TString.h"
+#include "TDirectory.h"
 #include "AliTPCCalPad.h"
 #include "AliTPCCalibPulser.h"
 #include "AliTPCCalibPedestal.h"
@@ -115,6 +118,7 @@ class AliTPCCalDet;
 #include "AliTPCTempMap.h"
 #include "AliTPCCalibVdrift.h"
 #include "AliTPCCalibRaw.h"
+#include "AliTPCParam.h"
 
 #include "AliTPCPreprocessorOnline.h"
 
@@ -188,8 +192,8 @@ AliTPCcalibDB::AliTPCcalibDB():
   fTemperatureArray(100000),    //! array of temperature sensors - per run - Just for calibration studies
   fVdriftArray(100000),                 //! array of v drift interfaces
   fDriftCorrectionArray(100000),  //! array of drift correction
-  fRunList(100000)              //! run list - indicates try to get the run param 
-
+  fRunList(100000),              //! run list - indicates try to get the run param 
+  fDButil(0)
 {
   //
   // constructor
@@ -225,7 +229,8 @@ AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& ):
   fTemperatureArray(0),   //! array of temperature sensors - per run - Just for calibration studies
   fVdriftArray(0),         //! array of v drift interfaces
   fDriftCorrectionArray(0),         //! array of v drift interfaces
-  fRunList(0)              //! run list - indicates try to get the run param 
+  fRunList(0),              //! run list - indicates try to get the run param 
+  fDButil(0)
 {
   //
   // Copy constructor invalid -- singleton implementation
@@ -293,11 +298,10 @@ void AliTPCcalibDB::SetRun(Long64_t run)
 
 void AliTPCcalibDB::Update(){
        //
-       AliCDBEntry * entry=0;
-  
+  AliCDBEntry * entry=0;
   Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
   AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
-  
+  fDButil = new AliTPCcalibDButil;   
   //
   entry          = GetCDBEntry("TPC/Calib/PadGainFactor");
   if (entry){
@@ -410,9 +414,9 @@ void AliTPCcalibDB::Update(){
   //  fExB=dynamic_cast<AliTPCExB*>(entry->GetObject()->Clone());
   //}
   //
-  // ExB  - calculate during initialization 
-  //      - 
-  fExB =  GetExB(-5,kTRUE);
+  // ExB  - calculate during initialization - in simulation /reconstruction
+  //      - not invoked here anymore
+  //fExB =  GetExB(-5,kTRUE);
      //
   if (!fTransform) {
     fTransform=new AliTPCTransform(); 
@@ -421,7 +425,6 @@ void AliTPCcalibDB::Update(){
 
   //
   AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
-  
 }
 
 
@@ -740,7 +743,7 @@ void AliTPCcalibDB::RegisterExB(Int_t index, Float_t bz, Bool_t bdelete){
   Float_t factor =  bz/(5.);  // default b filed in Cheb with minus sign
                               // was chenged in the Revision ???? (Ruben can you add here number)
   
-  AliMagF*   bmap = new AliMagF("MapsExB","MapsExB", 2,factor,1., 10.,AliMagF::k5kG,"$(ALICE_ROOT)/data/maps/mfchebKGI_sym.root");
+  AliMagF*   bmap = new AliMagF("MapsExB","MapsExB", factor,TMath::Sign(1.f,factor),AliMagF::k5kG);
   
   AliTPCExBFirst *exb  = new  AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
   AliTPCExB::SetInstance(exb);
@@ -792,6 +795,8 @@ void AliTPCcalibDB::UpdateRunInformations( Int_t run, Bool_t force){
   //
   // - > Don't use it for reconstruction - Only for Calibration studies
   //
+  if (run<0) return;
+  if (fRunList[run]>0 &&force==kFALSE) return;
   AliCDBEntry * entry = 0;
   if (run>= fRunList.GetSize()){
     fRunList.Set(run*2+1);
@@ -804,7 +809,9 @@ void AliTPCcalibDB::UpdateRunInformations( Int_t run, Bool_t force){
     fDriftCorrectionArray.Expand(run*2+1);
     fTimeGainSplinesArray.Expand(run*2+1);
   }
-  if (fRunList[run]>0 &&force==kFALSE) return;
+
+  fRunList[run]=1;  // sign as used
+
   //
   entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
   if (entry)  {
@@ -845,14 +852,31 @@ void AliTPCcalibDB::UpdateRunInformations( Int_t run, Bool_t force){
   if (entry)  {
     fTemperatureArray.AddAt(entry->GetObject(),run);
   }
-  fRunList[run]=1;  // sign as used
+  //apply fDButil filters
+
+  fDButil->UpdateFromCalibDB();
+  if (fTemperature) fDButil->FilterTemperature(fTemperature);
 
   AliDCSSensor * press = GetPressureSensor(run,0);
   AliTPCSensorTempArray * temp = GetTemperatureSensor(run);
-  if (press && temp){
+  Bool_t accept=kTRUE;
+  if (temp) {
+    accept = fDButil->FilterTemperature(temp)>0.1;
+  }
+  if (press) {
+    const Double_t kMinP=950.;
+    const Double_t kMaxP=1050.;
+    const Double_t kMaxdP=10.;
+    const Double_t kSigmaCut=4.;
+    fDButil->FilterSensor(press,kMinP,kMaxP,kMaxdP,kSigmaCut);
+    if (press->GetFit()==0) accept=kFALSE;
+  }
+  if (press && temp &&accept){
     AliTPCCalibVdrift * vdrift = new AliTPCCalibVdrift(temp, press,0);
     fVdriftArray.AddAt(vdrift,run);
   }
+  fDButil->FilterCE(120., 3., 4.,0);
+  fDButil->FilterTracks(run, 10.,0);
 }
 
 
@@ -864,6 +888,30 @@ Float_t AliTPCcalibDB::GetGain(Int_t sector, Int_t row, Int_t pad){
   return calPad->GetCalROC(sector)->GetValue(row,pad);
 }
 
+AliSplineFit* AliTPCcalibDB::GetVdriftSplineFit(const char* name, Int_t run){
+  //
+  //
+  //
+  TObjArray *arr=GetTimeVdriftSplineRun(run);
+  if (!arr) return 0;
+  return dynamic_cast<AliSplineFit*>(arr->FindObject(name));
+}
+
+AliSplineFit* AliTPCcalibDB::CreateVdriftSplineFit(const char* graphName, Int_t run){
+  //
+  // create spline fit from the drift time graph in TimeDrift
+  //
+  TObjArray *arr=GetTimeVdriftSplineRun(run);
+  if (!arr) return 0;
+  TGraph *graph=dynamic_cast<TGraph*>(arr->FindObject(graphName));
+  if (!graph) return 0;
+  AliSplineFit *fit = new AliSplineFit();
+  fit->SetGraph(graph);
+  fit->SetMinPoints(graph->GetN()+1);
+  fit->InitKnots(graph,2,0,0.001);
+  fit->SplineFit(0);
+  return fit;
+}
 
 AliGRPObject *AliTPCcalibDB::GetGRP(Int_t run){
   //
@@ -950,6 +998,18 @@ TObjArray * AliTPCcalibDB::GetTimeGainSplinesRun(Int_t run){
   return gainSplines;
 }
 
+TObjArray * AliTPCcalibDB::GetTimeVdriftSplineRun(Int_t run){
+  //
+  // Get drift spline array
+  //
+  TObjArray * driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
+  if (!driftSplines) {
+    UpdateRunInformations(run);
+    driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
+  }
+  return driftSplines;
+}
+
 AliDCSSensorArray * AliTPCcalibDB::GetVoltageSensors(Int_t run){
   //
   // Get temperature sensor array
@@ -1059,7 +1119,7 @@ Float_t AliTPCcalibDB::GetDCSSensorValue(AliDCSSensorArray *arr, Int_t timeStamp
     for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
       Double_t x,y;
       gr->GetPoint(ipoint,x,y);
-      Int_t time=sensor->GetStartTime()+x*3600; //time in graph is hours
+      Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
       if (time<timeStamp) continue;
       val=y;
       break;
@@ -1072,14 +1132,14 @@ Float_t AliTPCcalibDB::GetDCSSensorValue(AliDCSSensorArray *arr, Int_t timeStamp
     if (val==0 ){
       Double_t x,y;
       gr->GetPoint(0,x,y);
-      Int_t time=sensor->GetStartTime()+x*3600; //time in graph is hours
+      Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
       if ((time-timeStamp)<5*60) val=y;
     }
     //last point
     if (val==0 ){
       Double_t x,y;
       gr->GetPoint(gr->GetN()-1,x,y);
-      Int_t time=sensor->GetStartTime()+x*3600; //time in graph is hours
+      Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
       if ((timeStamp-time)<5*60) val=y;
     }
   } else {
@@ -1323,14 +1383,21 @@ Char_t  AliTPCcalibDB::GetL3Polarity(Int_t run) {
   //
   // get l3 polarity from GRP
   //
-  return AliTPCcalibDB::GetGRP(run)->GetL3Polarity();
+  Char_t pol=-100;
+  AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
+  if (grp) pol=grp->GetL3Polarity();
+  return pol;
 }
 
 TString AliTPCcalibDB::GetRunType(Int_t run){
   //
   // return run type from grp
   //
-  return AliTPCcalibDB::GetGRP(run)->GetRunType();
+
+//   TString type("UNKNOWN");
+  AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
+  if (grp) return grp->GetRunType();
+  return "UNKNOWN";
 }
 
 Float_t AliTPCcalibDB::GetValueGoofie(Int_t timeStamp, Int_t run, Int_t type){
@@ -1467,3 +1534,142 @@ Bool_t AliTPCcalibDB::CreateGUITree(Int_t run, const char* filename)
   return kTRUE;
 }
 
+Bool_t AliTPCcalibDB::CreateRefFile(Int_t run, const char* filename)
+{
+  //
+  // Create a gui tree for run number 'run'
+  //
+  
+  if (!AliCDBManager::Instance()->GetDefaultStorage()){
+    AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
+                    MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
+    return kFALSE;
+  }
+  TString file(filename);
+  if (file.IsNull()) file=Form("RefCalPads_%d.root",run);
+  TDirectory *currDir=gDirectory;
+  //db instance
+  AliTPCcalibDB *db=AliTPCcalibDB::Instance();
+  // retrieve cal pad objects
+  db->SetRun(run);
+  //open file
+  TFile f(file.Data(),"recreate");
+  //noise and pedestals
+  db->GetPedestals()->Write("Pedestals");
+  db->GetPadNoise()->Write("PadNoise");
+  //pulser data
+  db->GetPulserTmean()->Write("PulserTmean");
+  db->GetPulserTrms()->Write("PulserTrms");
+  db->GetPulserQmean()->Write("PulserQmean");
+  //CE data
+  db->GetCETmean()->Write("CETmean");
+  db->GetCETrms()->Write("CETrms");
+  db->GetCEQmean()->Write("CEQmean");
+  //Altro data
+  db->GetALTROAcqStart() ->Write("ALTROAcqStart");
+  db->GetALTROZsThr()    ->Write("ALTROZsThr");
+  db->GetALTROFPED()     ->Write("ALTROFPED");
+  db->GetALTROAcqStop()  ->Write("ALTROAcqStop");
+  db->GetALTROMasked()   ->Write("ALTROMasked");
+  //
+  f.Close();
+  currDir->cd();
+  return kTRUE;
+}
+
+
+
+Double_t AliTPCcalibDB::GetVDriftCorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
+  //
+  // Get time dependent drift velocity correction
+  // multiplication factor        vd = vdnom *(1+vdriftcorr)
+  // Arguments:
+  // mode determines the algorith how to combine the Laser Track, LaserCE and physics tracks
+  // timestamp - timestamp
+  // run       - run number
+  // side      - the drift velocity per side (possible for laser and CE)
+  //
+  // Notice - Extrapolation outside of calibration range  - using constant function
+  //
+  Double_t result;
+  // mode 1  automatic mode - according to the distance to the valid calibration
+  //                        -  
+  Double_t deltaP=0,  driftP=0,  wP  = 0.;
+  Double_t deltaLT=0, driftLT=0, wLT = 0.;
+  Double_t deltaCE=0, driftCE=0, wCE = 0.;
+  driftP  = fDButil->GetVDriftTPC(deltaP,run,timeStamp); 
+  driftCE = fDButil->GetVDriftTPCCE(deltaCE, run,timeStamp,36000,2);
+  driftLT = fDButil->GetVDriftTPCLaserTracks(deltaLT,run,timeStamp,36000,2);
+  deltaP   = TMath::Abs(deltaP);
+  deltaLT  = TMath::Abs(deltaLT);
+  deltaCE  = TMath::Abs(deltaCE);
+  if (mode==1) {
+    const Double_t kEpsilon=0.0000000001;
+    Double_t meanDist= (deltaP+deltaLT+deltaCE)*0.3;
+    if (meanDist<1.) return driftLT;
+    wP  = meanDist/(deltaP +0.005*meanDist);
+    wLT = meanDist/(deltaLT+0.005*meanDist);
+    wCE = meanDist/(deltaCE+0.001*meanDist);
+    if (TMath::Abs(driftCE)<kEpsilon) wCE=0;  // invalid calibration
+    result = (driftP*wP+driftLT*wLT+driftCE*wCE)/(wP+wLT+wCE);
+  }
+
+  return result;
+}
+
+Double_t AliTPCcalibDB::GetTime0CorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
+  //
+  // Get time dependent time 0 (trigger delay in cm) correction
+  // additive correction        time0 = time0+ GetTime0CorrectionTime
+  // Value etracted combining the vdrift correction using laser tracks and CE and the physics track matchin
+  // Arguments:
+  // mode determines the algorith how to combine the Laser Track and physics tracks
+  // timestamp - timestamp
+  // run       - run number
+  // side      - the drift velocity per side (possible for laser and CE)
+  //
+  // Notice - Extrapolation outside of calibration range  - using constant function
+  //
+  Double_t result=0;
+  if (mode==1) result=fDButil->GetTriggerOffsetTPC(run,timeStamp);    
+  result  *=fParam->GetZLength();
+
+  return result;
+
+}
+
+
+
+
+Double_t AliTPCcalibDB::GetVDriftCorrectionGy(Int_t timeStamp, Int_t run, Int_t side, Int_t /*mode*/){
+  //
+  // Get global y correction drift velocity correction factor
+  // additive factor        vd = vdnom*(1+GetVDriftCorrectionGy *gy)
+  // Value etracted combining the vdrift correction using laser tracks and CE
+  // Arguments:
+  // mode determines the algorith how to combine the Laser Track, LaserCE
+  // timestamp - timestamp
+  // run       - run number
+  // side      - the drift velocity gy correction per side (CE and Laser tracks)
+  //
+  // Notice - Extrapolation outside of calibration range  - using constant function
+  // 
+  if (run<=0 && fTransform) run = fTransform->GetCurrentRunNumber();
+  UpdateRunInformations(run,kFALSE);
+  TObjArray *array =AliTPCcalibDB::Instance()->GetTimeVdriftSplineRun(run);
+  if (!array) return 0;
+  TGraphErrors *laserA= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
+  TGraphErrors *laserC= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
+  
+  Double_t result=0;
+  if (laserA && laserC){
+   result= (laserA->Eval(timeStamp)+laserC->Eval(timeStamp))*0.5;
+  }
+  if (laserA && side==0){
+    result = (laserA->Eval(timeStamp));
+  }
+  if (laserC &&side==1){
+    result = (laserC->Eval(timeStamp));
+  }
+  return -result/250.; //normalized before
+}