]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
AliTPCcalibDB.cxx - Check the OCDB info consistency
authormarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 22 Oct 2009 08:36:04 +0000 (08:36 +0000)
committermarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 22 Oct 2009 08:36:04 +0000 (08:36 +0000)
AliTPCTransform.cxx  - cach current correction
CalibEnv.C           - exporting drift velocity for different type of post-proces algorithm

TPC/AliTPCTransform.cxx
TPC/AliTPCcalibDB.cxx
TPC/CalibMacros/CalibEnv.C

index 8877982cdfcc25351b65c679ec199bda90538e96..6a9f4268afca51169ffbbc964e6357b11be92cd8 100755 (executable)
@@ -246,35 +246,38 @@ void AliTPCTransform::Local2RotatedGlobal(Int_t sector, Double_t *x) const {
     }
   }
   //
-  //
-  Double_t vdcorrectionTime=1;
-  Double_t time0corrTime=0;
-  //
-  if(fCurrentRecoParam&&fCurrentRecoParam->GetUseDriftCorrectionTime()>0) {
-    vdcorrectionTime = (1+AliTPCcalibDB::Instance()->
-      GetVDriftCorrectionTime(fCurrentTimeStamp, 
-                             AliTPCcalibDB::Instance()->GetRun(),
-                             sector%36>=18,
-                             fCurrentRecoParam->GetUseDriftCorrectionTime()));                   
-    time0corrTime= AliTPCcalibDB::Instance()->
-      GetTime0CorrectionTime(fCurrentTimeStamp, 
-                            AliTPCcalibDB::Instance()->GetRun(),
-                            sector%36>=18,
-                            fCurrentRecoParam->GetUseDriftCorrectionTime());   
-  }
-  //
-  if(fCurrentRecoParam&&fCurrentRecoParam->GetUseDriftCorrectionGY()>0) {
-    Float_t xyzPad[3];
-    AliTPCROC::Instance()->GetPositionGlobal(sector, TMath::Nint(x[0]) ,TMath::Nint(x[1]), xyzPad);
-    
-    Double_t corrGy= (1+(xyzPad[1])*AliTPCcalibDB::Instance()->
-                     GetVDriftCorrectionGy(fCurrentTimeStamp, 
-                                           AliTPCcalibDB::Instance()->GetRun(),
-                                           sector%36>=18,
-                                           fCurrentRecoParam->GetUseDriftCorrectionGY()));
-    vdcorrectionTime *=corrGy;
+  // simple caching non thread save
+  static Double_t vdcorrectionTime=1;
+  static Double_t time0corrTime=0;
+  static Int_t    lastStampT=-1;
+  //
+  if (lastStampT!=fCurrentTimeStamp){
+    lastStampT=fCurrentTimeStamp;
+    if(fCurrentRecoParam&&fCurrentRecoParam->GetUseDriftCorrectionTime()>0) {
+      vdcorrectionTime = (1+AliTPCcalibDB::Instance()->
+                         GetVDriftCorrectionTime(fCurrentTimeStamp, 
+                                                 AliTPCcalibDB::Instance()->GetRun(),
+                                                 sector%36>=18,
+                                                 fCurrentRecoParam->GetUseDriftCorrectionTime()));                       
+      time0corrTime= AliTPCcalibDB::Instance()->
+       GetTime0CorrectionTime(fCurrentTimeStamp, 
+                              AliTPCcalibDB::Instance()->GetRun(),
+                              sector%36>=18,
+                              fCurrentRecoParam->GetUseDriftCorrectionTime()); 
+    }
+    //
+    if(fCurrentRecoParam&&fCurrentRecoParam->GetUseDriftCorrectionGY()>0) {
+      Float_t xyzPad[3];
+      AliTPCROC::Instance()->GetPositionGlobal(sector, TMath::Nint(x[0]) ,TMath::Nint(x[1]), xyzPad);
+      
+      Double_t corrGy= (1+(xyzPad[1])*AliTPCcalibDB::Instance()->
+                       GetVDriftCorrectionGy(fCurrentTimeStamp, 
+                                             AliTPCcalibDB::Instance()->GetRun(),
+                                             sector%36>=18,
+                                             fCurrentRecoParam->GetUseDriftCorrectionGY()));
+      vdcorrectionTime *=corrGy;
+    }
   }
-  
 
 
   if (!param){
index d02a64daecaa9827b6e37a6a2234ecc53b4bcc78..9b38f2a8df7daf057f554a1f4c7117b044f936ba 100644 (file)
@@ -296,8 +296,7 @@ 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
   
@@ -424,7 +423,6 @@ void AliTPCcalibDB::Update(){
 
   //
   AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
-  
 }
 
 
@@ -1570,6 +1568,8 @@ Double_t AliTPCcalibDB::GetVDriftCorrectionTime(Int_t timeStamp, Int_t run, Int_
   //
   // 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_DRIFT_LASER_ALL_A");
@@ -1601,6 +1601,8 @@ Double_t AliTPCcalibDB::GetTime0CorrectionTime(Int_t timeStamp, Int_t run, Int_t
   //
   // 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_DRIFT_LASER_ALL_A");
@@ -1641,7 +1643,9 @@ Double_t AliTPCcalibDB::GetVDriftCorrectionGy(Int_t timeStamp, Int_t run, Int_t
   // 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");
index 5704e51830c650e358b45724d4a7ff800e90a371..41ff7bad3d47852b99c9ca38fb71879d163cccbf 100644 (file)
@@ -7,7 +7,8 @@ gSystem->AddIncludePath("-I$ALICE_ROOT/TPC");
 CalibEnv("run.list");
 TFile f("dcsTime.root")
 */
+#include "TMVA/TSpline1.h"
+#include "TMVA/TSpline2.h"
 #include <iostream>
 #include <fstream>
 #include <stdio.h>
@@ -50,6 +51,8 @@ AliTPCcalibDButil *dbutil =0;
 TTree * dcsTree=0;
 TString refFile="dummy.root"; 
 TTreeSRedirector *pcstream =0;
+void GetNearest(TGraph *graph, Double_t x, Double_t &dx, Double_t &y);
+void GetInterpoloationSigma(TGraph *graph, Double_t x, Double_t &dx, Double_t &sy, Double_t deltaX);
 //
 //
 void ProcessRun(Int_t irun, Int_t startTime, Int_t endTime);
@@ -450,6 +453,7 @@ void ProcessDrift(Int_t run, Int_t timeStamp){
   static Double_t     vlaserA[3]={0,0,0};
   static Double_t     vlaserC[3]={0,0,0};
   static Double_t     vcosmicAll=0;
+  
 
   if (array){
     laserA[0]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DELAY_LASER_ALL_A");
@@ -475,4 +479,79 @@ void ProcessDrift(Int_t run, Int_t timeStamp){
     "vlaserC1="<<vlaserC[1]<<
     "vlaserC2="<<vlaserC[2]<<
     "vcosmicAll="<<vcosmicAll;
+
+  //
+  // define distance to measurement
+  //
+  static Double_t dlaserA=0; 
+  static Double_t dlaserC=0; 
+  static Double_t dcosmic=0; 
+  static Double_t slaserA=0; 
+  static Double_t slaserC=0; 
+  static Double_t scosmic=0; 
+  static Double_t  vclaserA[3]={0,0,0};
+  static Double_t  vclaserC[3]={0,0,0};
+  static Double_t  vccosmicAll=0;
+  Double_t dummy;
+  Double_t deltaX=1800; // +-0.5 hour
+  for (Int_t i=0;i<3;i++){
+    if (laserA[i]) GetNearest(laserA[i],timeStamp,dlaserA,vclaserA[i]);
+    if (laserC[i]) GetNearest(laserC[i],timeStamp,dlaserC,vclaserC[i]);
+  }  
+  if (laserA[1]) GetInterpoloationSigma(laserA[1],timeStamp,dummy,slaserA,deltaX);
+  if (laserC[1]) GetInterpoloationSigma(laserC[1],timeStamp,dummy,slaserC,deltaX);
+
+  if (cosmicAll) GetNearest(cosmicAll,timeStamp,dcosmic,vccosmicAll);
+  if (cosmicAll) GetInterpoloationSigma(cosmicAll,timeStamp,dummy,scosmic,deltaX);
+  (*pcstream)<<"dcs"<<
+    "vclaserA0="<<vclaserA[0]<<
+    "vclaserA1="<<vclaserA[1]<<
+    "vclaserA2="<<vclaserA[2]<<
+    "vclaserC0="<<vclaserC[0]<<
+    "vclaserC1="<<vclaserC[1]<<
+    "vclaserC2="<<vclaserC[2]<<
+    "vccosmicAll="<<vccosmicAll<<
+    "dlaserA="<<dlaserA<<
+    "dlaserC="<<dlaserC<<
+    "dcosmic="<<dcosmic<<
+    "slaserA="<<slaserA<<
+    "slaserC="<<slaserC<<
+    "scosmic="<<scosmic;
+}
+
+void GetNearest(TGraph *graph, Double_t xref, Double_t &dx, Double_t &y){
+  //
+  // find the closest point to xref  in x  direction
+  // return dx and value 
+  Int_t index=0;
+  index = TMath::BinarySearch(graph->GetN(), graph->GetX(),xref);
+  if (index<0) index=0;
+  if (index>=graph->GetN()-1) index=graph->GetN()-2;
+  if (xref-graph->GetX()[index]>graph->GetX()[index]-xref) index++;
+  dx = xref-graph->GetX()[index];
+  y  = graph->GetY()[index];
+}
+
+
+void GetInterpoloationSigma(TGraph *graph, Double_t x, Double_t &/*dx*/, Double_t &sy,Double_t deltaX){
+  //
+  //
+  //
+  TMVA::TSpline1 * spline1= new TMVA::TSpline1("spline1",graph);
+  Double_t deltas[1000];
+  Int_t counter=0;
+  Int_t index=0;
+  index = TMath::BinarySearch(graph->GetN(), graph->GetX(),x);
+  if (index<1) index=1;
+  if (index>=graph->GetN()-1) index=graph->GetN()-2;
+  
+  for (Int_t idelta=-10; idelta<=10; idelta++){
+    Double_t y0=0,dummy=0;
+    Double_t lx=x+idelta*deltaX/20.;
+    deltas[counter]=spline1->Eval(lx)-graph->GetY()[index];
+    counter++;
+    deltas[counter]=spline1->Eval(lx)-graph->GetY()[index+1];
+    counter++;
+  }
+  sy=TMath::RMS(counter,deltas);
 }