]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
calibration updates (Marian)
authormarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sun, 19 Oct 2008 15:08:27 +0000 (15:08 +0000)
committermarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sun, 19 Oct 2008 15:08:27 +0000 (15:08 +0000)
AliTPCcalibDB

Combine P and T measurement - provide correction factor
+  static Double_t GetPTRelative(UInt_t timeSec, Int_t run,  Int_t side);

AliTPCcalibTime, AliTPCcalibLaser
Adding the correction factors to debug streams - (see calibDB)

AliTPCCalibVdrift
Adding interface for pressure sensors

AliTPCcalibCosmic
Changed default name of the debug streamer

TPC/AliTPCCalibVdrift.cxx
TPC/AliTPCCalibVdrift.h
TPC/AliTPCcalibCosmic.cxx
TPC/AliTPCcalibDB.cxx
TPC/AliTPCcalibDB.h
TPC/AliTPCcalibTime.cxx
TPC/AliTPCcalibTime.h

index 5d5879af8168bbaf69b6ecee1925b967861bd645..b130494ad92ddace3120a0f7553057e0aacbc6f7 100644 (file)
@@ -31,7 +31,7 @@
 
 ClassImp(AliTPCCalibVdrift)
 
-  namespace paramDefinitions {
+namespace paramDefinitions {
     
     // Standard Conditions used as origin in the Magbolz simulations
     // Dimesions E [kV/cm], T [K], P [TORR], Cco2 [%], Cn2 [%]
@@ -57,10 +57,11 @@ ClassImp(AliTPCCalibVdrift)
 
 using namespace paramDefinitions;
 
-AliTPCCalibVdrift::AliTPCCalibVdrift(AliTPCSensorTempArray *SensTemp, TObject *SensPres, TObject *SensGasComp):
+AliTPCCalibVdrift::AliTPCCalibVdrift(AliTPCSensorTempArray *SensTemp, AliDCSSensor *SensPres, TObject *SensGasComp):
   TNamed(),
   fSensTemp(0),
   fSensPres(0),
+  fTempMap(0),
   fSensGasComp(0)
 {
   //
@@ -69,6 +70,7 @@ AliTPCCalibVdrift::AliTPCCalibVdrift(AliTPCSensorTempArray *SensTemp, TObject *S
 
   fSensTemp = SensTemp;
   fSensPres = SensPres;
+  fTempMap  = new AliTPCTempMap(fSensTemp);
   fSensGasComp = SensGasComp;
 }
 
@@ -76,6 +78,7 @@ AliTPCCalibVdrift::AliTPCCalibVdrift(const AliTPCCalibVdrift& source) :
   TNamed(source),
   fSensTemp(source.fSensTemp),
   fSensPres(source.fSensPres),
+  fTempMap(source.fTempMap),
   fSensGasComp(source.fSensGasComp)
 {
   //
@@ -100,9 +103,33 @@ AliTPCCalibVdrift::~AliTPCCalibVdrift()
 {
   //
   // AliTPCCalibVdrift destructor
+  // 
+
+}
+
+Double_t AliTPCCalibVdrift::GetPTRelative(UInt_t timeSec, Int_t side){
   //
+  // Get Relative difference of p/T for given time stamp
+  // timeSec - absolute time
+  // side    - 0 - A side 1 -C side
+  //
+  TTimeStamp tstamp(timeSec);
+  if (!fSensPres) return 0;
+  Double_t pressure    = fSensPres->GetValue(tstamp);
+  TLinearFitter * fitter = fTempMap->GetLinearFitter(3,side,tstamp);
+  if (!fitter) return 0;
+  TVectorD vec;
+  fitter->GetParameters(vec);
+  delete fitter;
+  if (vec[0]<10) return 0;
+  Double_t  temperature = vec[0]+273.15;
+  Double_t povertMeas = pressure/temperature;
+  const Double_t torrTokPascal = 0.75006;
+  Double_t povertNom =  kstdP/(torrTokPascal*kstdT);
+  return povertMeas/povertNom;
 }
 
+
 //_____________________________________________________________________________
 Double_t AliTPCCalibVdrift::VdriftLinearHyperplaneApprox(Double_t dE, Double_t dT, Double_t dP, Double_t dCco2, Double_t dCn2) 
 {
@@ -137,7 +164,7 @@ Double_t AliTPCCalibVdrift::GetVdriftChange(Double_t x, Double_t y, Double_t z,
   Double_t dE = 0; //FIXME: eventually include Field-Inhomogenities
 
   // Get Temperature Value ----------------------  
-  AliTPCTempMap *tempMap = new AliTPCTempMap(fSensTemp);
+  AliTPCTempMap *tempMap = fTempMap;
   Double_t tempValue = tempMap->GetTemperature(x, y, z, timeSec);
   Double_t dT = tempValue+273.15 - kstdT;
   
@@ -159,9 +186,7 @@ Double_t AliTPCCalibVdrift::GetVdriftChange(Double_t x, Double_t y, Double_t z,
   // Calculate change in drift velocity in terms of Vdrift_nominal
   Double_t vdrift = VdriftLinearHyperplaneApprox(dE, dT, dP, dCco2, dCn2); 
   
-  return vdrift;
-    
-  tempMap->~AliTPCTempMap();
+  return vdrift;    
 }
 
 //_____________________________________________________________________________
index f505bb017dfb69c635c43e8ad75de6b276eb9905..90aa3c76b840df3a48ee8718757ae4be47bb8972 100644 (file)
@@ -15,11 +15,18 @@ class TGraph;
 class AliTPCCalibVdrift : public TNamed {
 
 public:
-  AliTPCCalibVdrift(AliTPCSensorTempArray *SensTemp, TObject *SensPres, TObject *SensGasComp);
+  AliTPCCalibVdrift(AliTPCSensorTempArray *SensTemp, AliDCSSensor *SensPres, TObject *SensGasComp);
   AliTPCCalibVdrift(const AliTPCCalibVdrift& source);
   virtual ~AliTPCCalibVdrift();
   AliTPCCalibVdrift& operator=(const AliTPCCalibVdrift& source);
-
+  //
+  // Interface for the reconstruction
+  //
+  Double_t GetPTRelative(UInt_t timeSec, Int_t side);
+
+  //
+  // Stefan interfaces - for v drift study
+  //
   Double_t VdriftLinearHyperplaneApprox(Double_t dE, Double_t dT, Double_t dP, Double_t dCco2, Double_t dCn2);
   
   Double_t GetVdriftNominal();
@@ -32,7 +39,8 @@ public:
 protected:
 
   AliTPCSensorTempArray *fSensTemp;   // Temperature sensors 
-  TObject *fSensPres;         // Placeholder for Pressure sensors
+  AliDCSSensor          *fSensPres;   // pressure sensors
+  AliTPCTempMap         *fTempMap;    // Temerature sensor map
   TObject *fSensGasComp;      // placeholder for GasConzentration infos  
   
   ClassDef(AliTPCCalibVdrift,1);
index 35f442437eee9a0fae04fa255c7352cb820cdf73..ab3d7b5f7bebdcf2eb6d44a7617ec27e9e26ef6d 100644 (file)
@@ -545,46 +545,46 @@ gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
 gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
 AliXRDPROOFtoolkit tool; 
 TChain * chain = tool.MakeChain("cosmic.txt","Track0",0,1000000);
-chain->Lookup();
+chainCosmic->Lookup();
 
 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.01");  // OK
 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<2");     // OK
 TCut cutP1("cutP1","abs(Tr0.fP[1]-Tr1.fP[1])<10");   // OK
 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<0.1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
 TCut cutN("cutN","min(Orig0.fTPCncls,Orig1.fTPCncls)>100");
-TCut cutA=cutT+cutD+cutPt+cutN+cutP1;
+TCut cutA=cutT+cutD+cutPt+cutN+cutP1+"trigger!=16";
 
 TCut cutS("cutS","Orig0.fIp.fP[1]*Orig1.fIp.fP[1]>0");
 
 //
-chain->Draw(">>listEL",cutA,"entryList");
+chainCosmic->Draw(">>listEL",cutA,"entryList");
 TEntryList *elist = (TEntryList*)gDirectory->Get("listEL");
-chain->SetEntryList(elist);
+chainCosmic->SetEntryList(elist);
 //
-chain->Draw(">>listV40Z100","abs(d0)<40&&abs(v01)<100","entryList");
+chainCosmic->Draw(">>listV40Z100","abs(d0)<40&&abs(v01)<100","entryList");
 TEntryList *elistV40Z100 = (TEntryList*)gDirectory->Get("listV40Z100");
-chain->SetEntryList(elistV40Z100);
+chainCosmic->SetEntryList(elistV40Z100);
 
 //
 // Aliases
 //
-chain->SetAlias("side","(-1+(Tr0.fP[1]>0)*2)");
-chain->SetAlias("hpt","abs(Tr0.fP[4])<0.2");
-chain->SetAlias("signy","(-1+(Tr0.fP[0]>0)*2)");
+chainCosmic->SetAlias("side","(-1+(Tr0.fP[1]>0)*2)");
+chainCosmic->SetAlias("hpt","abs(Tr0.fP[4])<0.2");
+chainCosmic->SetAlias("signy","(-1+(Tr0.fP[0]>0)*2)");
 
-chain->SetAlias("dy","Tr0.fP[0]+Tr1.fP[0]");
-chain->SetAlias("dz","Tr0.fP[1]-Tr1.fP[1]");
-chain->SetAlias("d1pt","Tr0.fP[4]+Tr1.fP[4]");
-chain->SetAlias("dtheta","(Tr0.fP[3]+Tr1.fP[3])");
-chain->SetAlias("dphi","(Tr0.fAlpha-Tr1.fAlpha-pi)");
+chainCosmic->SetAlias("dy","Tr0.fP[0]+Tr1.fP[0]");
+chainCosmic->SetAlias("dz","Tr0.fP[1]-Tr1.fP[1]");
+chainCosmic->SetAlias("d1pt","Tr0.fP[4]+Tr1.fP[4]");
+chainCosmic->SetAlias("dtheta","(Tr0.fP[3]+Tr1.fP[3])");
+chainCosmic->SetAlias("dphi","(Tr0.fAlpha-Tr1.fAlpha-pi)");
 
-chain->SetAlias("mtheta","(Tr0.fP[3]-Tr1.fP[3])*0.5")
-chain->SetAlias("sa","(sin(Tr0.fAlpha+0.))");
-chain->SetAlias("ca","(cos(Tr0.fAlpha+0.))");
+chainCosmic->SetAlias("mtheta","(Tr0.fP[3]-Tr1.fP[3])*0.5")
+chainCosmic->SetAlias("sa","(sin(Tr0.fAlpha+0.))");
+chainCosmic->SetAlias("ca","(cos(Tr0.fAlpha+0.))");
 
 
 
-chain->Draw("dy:sqrt(abs(Tr0.fP[4]))>>hisdyA(5,0,1,50,-1,1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]>0");
+chainCosmic->Draw("dy:sqrt(abs(Tr0.fP[4]))>>hisdyA(5,0,1,50,-1,1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]>0");
 hisdyA->FitSlicesY();
 hisdyA_2->SetXTitle("#sqrt{1/p_{t}}");
 hisdyA_2->SetYTitle("#sigma_{y}(cm)");
@@ -592,7 +592,7 @@ hisdyA_2->SetTitle("Cosmic - Y matching");
 hisdyA_2->SetMaximum(0.5);
 
 
-chain->Draw("dy:sqrt(abs(Tr0.fP[4]))>>hisdyC(5,0,1,50,-1,1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]<0");
+chainCosmic->Draw("dy:sqrt(abs(Tr0.fP[4]))>>hisdyC(5,0,1,50,-1,1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]<0");
 hisdyC->FitSlicesY();
 hisdyC_2->SetXTitle("#sqrt{1/p_{t}}");
 hisdyC_2->SetYTitle("#sigma_{y}(cm)");
@@ -602,19 +602,19 @@ hisdyC_2->SetMinimum(0);
 hisdyC_2->SetMarkerStyle(22);
 hisdyA_2->SetMarkerStyle(21);
 hisdyC_2->SetMarkerSize(1.5);
-hisdzA_2->SetMarkerSize(1.5);
+hisdyA_2->SetMarkerSize(1.5);
 hisdyC_2->Draw();
 hisdyA_2->Draw("same");
 gPad->SaveAs("~/Calibration/Cosmic/pic/ymatching.gif")
 
-chain->Draw("dz:sqrt(abs(Tr0.fP[4]))>>hisdzA(5,0,1,50,-1,1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]>0");
+chainCosmic->Draw("dz:sqrt(abs(Tr0.fP[4]))>>hisdzA(5,0,1,50,-1,1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]>0");
 hisdzA->FitSlicesY();
 hisdzA_2->SetXTitle("#sqrt{1/p_{t}}");
 hisdzA_2->SetYTitle("#sigma_{z}(cm)");
 hisdzA_2->SetTitle("Cosmic - Z matching - A side ");
 hisdzA_2->SetMaximum(0.5);
 
-chain->Draw("dz:sqrt(abs(Tr0.fP[4]))>>hisdzC(5,0,1,50,-1,1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]<0");
+chainCosmic->Draw("dz:sqrt(abs(Tr0.fP[4]))>>hisdzC(5,0,1,50,-1,1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]<0");
 hisdzC->FitSlicesY();
 hisdzC_2->SetXTitle("#sqrt{1/p_{t}}");
 hisdzC_2->SetYTitle("#sigma_{z}(cm)");
@@ -632,14 +632,14 @@ hisdzA_2->Draw("same");
 //
 // PICTURE 1/pt
 //
-chain->Draw("d1pt:sqrt(abs(Tr0.fP[4]))>>hisd1ptA(5,0,1,30,-0.1,0.1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]>0");
+chainCosmic->Draw("d1pt:sqrt(abs(Tr0.fP[4]))>>hisd1ptA(5,0,1,30,-0.1,0.1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]>0");
 hisd1ptA->FitSlicesY();
 hisd1ptA_2->SetXTitle("#sqrt{1/p_{t}}");
 hisd1ptA_2->SetYTitle("#sigma_{z}(cm)");
 hisd1ptA_2->SetTitle("Cosmic - Z matching - A side ");
 hisd1ptA_2->SetMaximum(0.5);
 
-chain->Draw("d1pt:sqrt(abs(Tr0.fP[4]))>>hisd1ptC(5,0,1,30,-0.1,0.1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]<0");
+chainCosmic->Draw("d1pt:sqrt(abs(Tr0.fP[4]))>>hisd1ptC(5,0,1,30,-0.1,0.1)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]<0");
 hisd1ptC->FitSlicesY();
 hisd1ptC_2->SetXTitle("#sqrt{1/p_{t}}");
 hisd1ptC_2->SetYTitle("#sigma_{1/pt}(1/GeV)");
@@ -657,14 +657,14 @@ gPad->SaveAs("~/Calibration/Cosmic/pic/1ptmatching.gif")
 //
 // Theta
 //
-chain->Draw("dtheta:sqrt(abs(Tr0.fP[4]))>>hisdthetaA(5,0,1,30,-0.01,0.01)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]>0");
+chainCosmic->Draw("dtheta:sqrt(abs(Tr0.fP[4]))>>hisdthetaA(5,0,1,30,-0.01,0.01)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]>0");
 hisdthetaA->FitSlicesY();
 hisdthetaA_2->SetXTitle("#sqrt{1/p_{t}}");
 hisdthetaA_2->SetYTitle("#sigma_{#theta}(cm)");
 hisdthetaA_2->SetTitle("Cosmic - Z matching - A side ");
 hisdthetaA_2->SetMaximum(0.5);
 
-chain->Draw("dtheta:sqrt(abs(Tr0.fP[4]))>>hisdthetaC(5,0,1,30,-0.01,0.01)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]<0");
+chainCosmic->Draw("dtheta:sqrt(abs(Tr0.fP[4]))>>hisdthetaC(5,0,1,30,-0.01,0.01)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]<0");
 hisdthetaC->FitSlicesY();
 hisdthetaC_2->SetXTitle("#sqrt{1/p_{t}}");
 hisdthetaC_2->SetYTitle("#sigma_{#theta}(rad)");
@@ -682,14 +682,14 @@ gPad->SaveAs("~/Calibration/Cosmic/pic/thetamatching.gif")
 //
 // Phi
 //
-chain->Draw("dphi:sqrt(abs(Tr0.fP[4]))>>hisdphiA(5,0,1,30,-0.01,0.01)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]>0");
+chainCosmic->Draw("dphi:sqrt(abs(Tr0.fP[4]))>>hisdphiA(5,0,1,30,-0.01,0.01)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]>0");
 hisdphiA->FitSlicesY();
 hisdphiA_2->SetXTitle("#sqrt{1/p_{t}}");
 hisdphiA_2->SetYTitle("#sigma_{#phi}(rad)");
 hisdphiA_2->SetTitle("Cosmic - Z matching - A side ");
 hisdphiA_2->SetMaximum(0.5);
 
-chain->Draw("dphi:sqrt(abs(Tr0.fP[4]))>>hisdphiC(5,0,1,30,-0.01,0.01)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]<0");
+chainCosmic->Draw("dphi:sqrt(abs(Tr0.fP[4]))>>hisdphiC(5,0,1,30,-0.01,0.01)","!crossO&&!crossI&&abs(d0)<40&&Tr0.fP[1]<0");
 hisdphiC->FitSlicesY();
 hisdphiC_2->SetXTitle("#sqrt{1/p_{t}}");
 hisdphiC_2->SetYTitle("#sigma_{#phi}(rad)");
@@ -743,7 +743,7 @@ fstring+="side*sa*mtheta++";
 
 
 TString *strTheta0 = toolkit.FitPlane(chain,"dtheta",fstring->Data(), "hpt&&!crossI&&!crossO", chi2,npoints,fitParamA0,covMatrix,0.8);
-chain->SetAlias("dtheta0",strTheta0.Data());
+chainCosmic->SetAlias("dtheta0",strTheta0.Data());
 strTheta0->Tokenize("+")->Print();
 
 
@@ -774,9 +774,9 @@ strTheta0->Tokenize("+")->Print();
  TMatrixD covMatrix;
  
  //Theta
-chain->SetAlias("dthe","(Tr0.fP[3]+Tr1.fP[3])");
-chain->SetAlias("sign","(-1+(Tr0.fP[1]>0)*2)");
-chain->SetAlias("di","(1-abs(Tr0.fP[1])/250)");
+chainCosmic->SetAlias("dthe","(Tr0.fP[3]+Tr1.fP[3])");
+chainCosmic->SetAlias("sign","(-1+(Tr0.fP[1]>0)*2)");
+chainCosmic->SetAlias("di","(1-abs(Tr0.fP[1])/250)");
 //
 //
 TString strFit="";
@@ -797,8 +797,8 @@ strFit+="cos(Tr0.fAlpha)*Tr0.fP[3]++";         //8
  //                                        
  TString * thetaParam = toolkit.FitPlane(chain,"dthe", strFit.Data(),"1", chi2,npoints,fitParam,covMatrix,0.8,0,10000)
  
- chain->SetAlias("corTheta",thetaParam->Data());
- chain->Draw("dthe:Tr0.fP[1]","","",50000);
+ chainCosmic->SetAlias("corTheta",thetaParam->Data());
+ chainCosmic->Draw("dthe:Tr0.fP[1]","","",50000);
 
 
 
@@ -823,12 +823,12 @@ void AliTPCcalibCosmic::dEdxCorrection(){
   gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
   gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
  AliXRDPROOFtoolkit tool; 
-  TChain * chain = tool.MakeChain("cosmic.txt","Track0",0,1000000);
-  chain->Lookup();
+  TChain * chainCosmic = tool.MakeChain("cosmic.txt","Track0",0,1000000);
+  chainCosmic->Lookup();
   
-  chain->Draw(">>listEL",cutA,"entryList");
+  chainCosmic->Draw(">>listEL",cutA,"entryList");
   TEntryList *elist = (TEntryList*)gDirectory->Get("listEL");
-  chain->SetEntryList(elist);
+  chainCosmic->SetEntryList(elist);
 
   .x ~/rootlogon.C
    gSystem->Load("libSTAT.so");
@@ -838,7 +838,7 @@ void AliTPCcalibCosmic::dEdxCorrection(){
   TVectorD fitParam;
   TMatrixD covMatrix;
   
-  chain->Draw("Tr0.fP[4]+Tr1.fP[4]","OK"+cutA);
+  chainCosmic->Draw("Tr0.fP[4]+Tr1.fP[4]","OK"+cutA);
   
   TString strFit;
   strFit+="(Tr0.fP[1]/250)++";
@@ -872,15 +872,15 @@ TCut cutA=cutT+cutD+cutPt+cutN;
 
 
 
-TTree * chain = Track0;
+TTree * chainCosmic = Track0;
 
 
-chain->SetAlias("norm","signalTot0.fElements[3]/signalTot1.fElements[3]");
+chainCosmic->SetAlias("norm","signalTot0.fElements[3]/signalTot1.fElements[3]");
 //
-chain->SetAlias("dr1","(signalTot0.fElements[1]/signalTot0.fElements[3])");
-chain->SetAlias("dr2","(signalTot0.fElements[2]/signalTot0.fElements[3])");
-chain->SetAlias("dr4","(signalTot0.fElements[4]/signalTot0.fElements[3])");
-chain->SetAlias("dr5","(signalTot0.fElements[5]/signalTot0.fElements[3])");
+chainCosmic->SetAlias("dr1","(signalTot0.fElements[1]/signalTot0.fElements[3])");
+chainCosmic->SetAlias("dr2","(signalTot0.fElements[2]/signalTot0.fElements[3])");
+chainCosmic->SetAlias("dr4","(signalTot0.fElements[4]/signalTot0.fElements[3])");
+chainCosmic->SetAlias("dr5","(signalTot0.fElements[5]/signalTot0.fElements[3])");
 
 TString fstring="";  
 fstring+="dr1++";
@@ -899,17 +899,17 @@ fstring+="dr4*dr5++";
 
 TString *strqdedx = toolkit.FitPlane(chain,"norm",fstring->Data(), cutA, chi2,npoints,fitParam,covMatrix,-1,0,200000);
   
-chain->SetAlias("corQT",strqdedx->Data());
+chainCosmic->SetAlias("corQT",strqdedx->Data());
 
 */
 
 
 /*
-  chain->SetProof(kTRUE);
-  chain->Draw("Seed0.CookdEdxNorm(0,0.6,1,0,159,0,kTRUE,kTRUE):Seed0.CookdEdxNorm(0,0.6,1,0,159,0,kFALSE,kTRUE)",""+cutA,"",1000);
+  chainCosmic->SetProof(kTRUE);
+  chainCosmic->Draw("Seed0.CookdEdxNorm(0,0.6,1,0,159,0,kTRUE,kTRUE):Seed0.CookdEdxNorm(0,0.6,1,0,159,0,kFALSE,kTRUE)",""+cutA,"",1000);
 
 
-chain->Draw("Seed0.CookdEdxNorm(0,0.6,1,0,159,0,kTRUE,kTRUE)/Seed1.CookdEdxNorm(0,0.6,1,0,159,0,kTRUE,kTRUE)>>his(100,0.5,1.5)","min(Orig0.fTPCncls,Orig1.fTPCncls)>130"+cutA,"",50000);
+chainCosmic->Draw("Seed0.CookdEdxNorm(0,0.6,1,0,159,0,kTRUE,kTRUE)/Seed1.CookdEdxNorm(0,0.6,1,0,159,0,kTRUE,kTRUE)>>his(100,0.5,1.5)","min(Orig0.fTPCncls,Orig1.fTPCncls)>130"+cutA,"",50000);
 
 */
 
index 1edee3cc8ece5f3e6ed2eef2601bc57c30c9bd62..3facb6b59088ed2f70d06260306adc0bae15edf0 100644 (file)
@@ -113,6 +113,7 @@ class AliTPCCalDet;
 #include "AliTPCCalibCE.h"
 #include "AliTPCExBFirst.h"
 #include "AliTPCTempMap.h"
+#include "AliTPCCalibVdrift.h"
 
 
 
@@ -173,9 +174,10 @@ AliTPCcalibDB::AliTPCcalibDB():
   fMapping(0),
   fParam(0),
   fClusterParam(0),  
-  fGRPArray(100000),          //! array of GRPs  -  per run  - JUST for calibration studies
-  fGoofieArray(100000),        //! array of GOOFIE values -per run - Just for calibration studies
-  fTemperatureArray(100000),   //! array of temperature sensors - per run - Just for calibration studies
+  fGRPArray(100000),            //! array of GRPs  -  per run  - JUST for calibration studies
+  fGoofieArray(100000),         //! array of GOOFIE values -per run - Just for calibration studies
+  fTemperatureArray(100000),    //! array of temperature sensors - per run - Just for calibration studies
+  fVdriftArray(100000),                 //! array of v drift interfaces
   fRunList(100000)              //! run list - indicates try to get the run param 
 
 {
@@ -203,6 +205,7 @@ AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& ):
   fGRPArray(0),          //! array of GRPs  -  per run  - JUST for calibration studies
   fGoofieArray(0),        //! array of GOOFIE values -per run - Just for calibration studies
   fTemperatureArray(0),   //! array of temperature sensors - per run - Just for calibration studies
+  fVdriftArray(0),         //! array of v drift interfaces
   fRunList(0)              //! run list - indicates try to get the run param 
 {
   //
@@ -726,7 +729,10 @@ void AliTPCcalibDB::GetRunInformations( Int_t run){
   AliCDBEntry * entry = 0;
   if (run>= fRunList.GetSize()){
     fRunList.Set(run*2+1);
-    fGRPArray.Expand(run*2+1);fGoofieArray.Expand(run*2+1); fTemperatureArray.Expand(run*2+1);
+    fGRPArray.Expand(run*2+1);
+    fGoofieArray.Expand(run*2+1); 
+    fTemperatureArray.Expand(run*2+1);
+    fVdriftArray.Expand(run*2+1);
   }
   if (fRunList[run]>0) return;
   entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
@@ -736,6 +742,13 @@ void AliTPCcalibDB::GetRunInformations( Int_t run){
   entry = AliCDBManager::Instance()->Get("TPC/Calib/Temperature",run);
   if (entry)  fTemperatureArray.AddAt(entry->GetObject(),run);
   fRunList[run]=1;  // sign as used
+
+  AliDCSSensor * press = GetPressureSensor(run);
+  AliTPCSensorTempArray * temp = GetTemperatureSensor(run);
+  if (press && temp){
+    AliTPCCalibVdrift * vdrift = new AliTPCCalibVdrift(temp, press,0);
+    fVdriftArray.AddAt(vdrift,run);
+  }
 }
 
 
@@ -747,7 +760,7 @@ Float_t AliTPCcalibDB::GetGain(Int_t sector, Int_t row, Int_t pad){
   return calPad->GetCalROC(sector)->GetValue(row,pad);
 }
 
-AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run){
+AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run, Int_t type){
   //
   //
   AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
@@ -757,8 +770,8 @@ AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run){
     if (!grpRun) return 0; 
   }
   AliDCSSensor * sensor = grpRun->GetCavernAtmosPressure();
-  return sensor;
-  
+  if (type==1) sensor = grpRun->GetSurfaceAtmosPressure();
+  return sensor; 
 }
 
 AliTPCSensorTempArray * AliTPCcalibDB::GetTemperatureSensor(Int_t run){
@@ -785,15 +798,27 @@ AliDCSSensorArray * AliTPCcalibDB::GetGoofieSensors(Int_t run){
   return goofieArray;
 }
 
+AliTPCCalibVdrift *     AliTPCcalibDB::GetVdrift(Int_t run){
+  //
+  // Get the interface to the the vdrift 
+  //
+  AliTPCCalibVdrift  * vdrift = (AliTPCCalibVdrift*)fVdriftArray.At(run);
+  if (!vdrift) {
+    GetRunInformations(run);
+    vdrift= (AliTPCCalibVdrift*)fVdriftArray.At(run);
+  }
+  return vdrift;
+}
+
 
 
 
-Float_t AliTPCcalibDB::GetPressure(Int_t timeStamp, Int_t run){
+Float_t AliTPCcalibDB::GetPressure(Int_t timeStamp, Int_t run, Int_t type){
   //
   // GetPressure for given time stamp and runt
   //
   TTimeStamp stamp(timeStamp);
-  AliDCSSensor * sensor = Instance()->GetPressureSensor(run);
+  AliDCSSensor * sensor = Instance()->GetPressureSensor(run,type);
   if (!sensor) return 0;
   if (!sensor->GetFit()) return 0;
   return sensor->GetValue(stamp);
@@ -834,6 +859,17 @@ Float_t AliTPCcalibDB::GetTemperature(Int_t timeStamp, Int_t run, Int_t side){
 }
 
 
+Double_t AliTPCcalibDB::GetPTRelative(UInt_t timeSec, Int_t run, Int_t side){
+  //
+  // Get relative P/T 
+  // time - absolute time
+  // run  - run number
+  // side - 0 - A side   1-C side
+  AliTPCCalibVdrift * vdrift =  Instance()->GetVdrift(run);
+  if (!vdrift) return 0;
+  return vdrift->GetPTRelative(timeSec,side);
+}
+
 
 void AliTPCcalibDB::ProcessEnv(const char * runList){
   //
index 329fd8627c31381bf5baf5e3c48a51071c26c74f..ee79d6bf0f091178d403b9a80b6b9279249c6303 100644 (file)
@@ -26,6 +26,7 @@ class AliTPCAltroMapping;
 class AliTPCClusterParam;
 class AliDCSSensor;
 class AliDCSSensorArray;
+class AliTPCCalibVdrift;
 //class AliCDBStorage;
 
 class AliTPCcalibDB : public TObject
@@ -51,12 +52,15 @@ class AliTPCcalibDB : public TObject
   AliTPCAltroMapping ** GetMapping(){ return fMapping;}
   AliTPCClusterParam *GetClusterParam(){ return fClusterParam;}
   //
-  static Float_t GetPressure(Int_t timeStamp, Int_t run);
+  //
+  static Float_t GetPressure(Int_t timeStamp, Int_t run, Int_t type=0);
   static Bool_t  GetTemperatureFit(Int_t timeStamp, Int_t run, Int_t side,TVectorD& fit);
   static Float_t GetTemperature(Int_t timeStamp, Int_t run, Int_t side);
-  AliDCSSensor * GetPressureSensor(Int_t run);
+  static Double_t GetPTRelative(UInt_t timeSec, Int_t run,  Int_t side);
+  AliDCSSensor * GetPressureSensor(Int_t run, Int_t type=0);
   AliTPCSensorTempArray * GetTemperatureSensor(Int_t run);
   AliDCSSensorArray *     GetGoofieSensors(Int_t run);
+  AliTPCCalibVdrift *     GetVdrift(Int_t run);
   static Float_t GetGain(Int_t sector, Int_t row, Int_t pad);
   //
   static void     CreateObjectList(const Char_t *filename, TObjArray *calibObjects);
@@ -96,6 +100,7 @@ protected:
   TObjArray      fGRPArray;           //! array of GRPs  -  per run
   TObjArray      fGoofieArray;        //! array of GOOFIE values -per run
   TObjArray      fTemperatureArray;   //! array of temperature sensors - per run
+  TObjArray      fVdriftArray;        //! array of v drift interfaces
   TArrayI        fRunList;            //! run list - indicates try to get the run param
   //
   static AliTPCcalibDB* fgInstance;  // singleton control
index 73c139670451c1a095f51d8670f85ce9842a4385..3e65acdb841d74d01a653e1811ebe6be7095c733 100644 (file)
@@ -43,6 +43,13 @@ cal->GetHistVdrift()->Projection(1,0)->Draw()
 
     4. Analysis using debug streamers.    
 
+    gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
+    gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
+    AliXRDPROOFtoolkit tool;
+    TChain * chainTime = tool.MakeChain("time.txt","timeInfo",0,10200);
+    chainTime->Lookup();
+
+
 */
 
 
@@ -81,18 +88,28 @@ cal->GetHistVdrift()->Projection(1,0)->Draw()
 
 #include "TTreeStream.h"
 #include "AliTPCTracklet.h"
+#include "TTimeStamp.h"
+#include "AliTPCcalibDB.h"
+#include "AliTPCcalibLaser.h"
 
 ClassImp(AliTPCcalibTime)
 
 
 AliTPCcalibTime::AliTPCcalibTime() 
-  :AliTPCcalibBase(),
+  :AliTPCcalibBase(), 
+   fTriggerMask(0),
    fHistDeDxTgl(0),
    fHistDeDx(0),
    fHistVdrift(0),
    fIntegrationTimeDeDx(0),
-   fIntegrationTimeVdrift(0),
+   fIntegrationTimeVdrift(0),  
+   fLaser(0),       // pointer to laser calibration
+   fDz(0),          // current delta z
+   fdEdx(0),        // current dEdx
+   fdEdxRatio(0),   // current dEdx ratio
+   fTl(0),          // current tan(lambda)
    fCutMaxD(5),        // maximal distance in rfi ditection
+   fCutMaxDz(20),        // maximal distance in rfi ditection
    fCutTheta(0.03),    // maximal distan theta
    fCutMinDir(-0.99)   // direction vector products
 
@@ -108,8 +125,14 @@ AliTPCcalibTime::AliTPCcalibTime(const Text_t *name, const Text_t *title, ULong6
    fHistDeDx(0),
    fHistVdrift(0),
    fIntegrationTimeDeDx(0),
-   fIntegrationTimeVdrift(0),
+   fIntegrationTimeVdrift(0),  
+   fLaser(0),       // pointer to laser calibration
+   fDz(0),          // current delta z
+   fdEdx(0),        // current dEdx
+   fdEdxRatio(0),   // current dEdx ratio
+   fTl(0),          // current tan(lambda)
    fCutMaxD(5),        // maximal distance in rfi ditection
+   fCutMaxDz(20),        // maximal distance in rfi ditection
    fCutTheta(0.03),    // maximal distan theta
    fCutMinDir(-0.99)   // direction vector products
 {
@@ -150,14 +173,78 @@ AliTPCcalibTime::~AliTPCcalibTime(){
   //
   //
 }
+void AliTPCcalibTime::ResetCurrent(){
+  //
+  // reset current values
+  //
+  fDz=0;          // current delta z
+  fdEdx=0;        // current dEdx
+  fdEdxRatio=0;   // current dEdx ratio
+  fTl=0;          // current tan(lambda)
+
+}
+
 
 void AliTPCcalibTime::Process(AliESDEvent *event) {
   //
   //
   //
-  
+  Int_t ntracks=event->GetNumberOfTracks();
+  if (ntracks<2) return;
+  ResetCurrent();
+  //
   ProcessCosmic(event);
-
+  if (fTrigger==16){
+    if (!fLaser) fLaser =  new AliTPCcalibLaser("laserTPC","laserTPC",kFALSE);
+    fLaser->Process(event);
+  }
+  //  
+  // fill debug streamer
+  if (fStreamLevel>0 && fDz!=0){
+    TTreeSRedirector *cstream = GetDebugStreamer();
+    if (cstream){
+      TTimeStamp tstamp(fTime);
+      Float_t valuePressure0  = AliTPCcalibDB::GetPressure(tstamp,fRun,0);
+      Float_t valuePressure1 = AliTPCcalibDB::GetPressure(tstamp,fRun,1);
+      Double_t ptrelative0   = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
+      Double_t ptrelative1   = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
+      Double_t temp0         = AliTPCcalibDB::GetTemperature(tstamp,fRun,0);
+      Double_t temp1         = AliTPCcalibDB::GetTemperature(tstamp,fRun,1);
+      TVectorD vdriftA, vdriftC,vdriftAC;
+      if (fLaser && fTrigger==16) {
+       if (fLaser->fFitAside)  vdriftA=*(fLaser->fFitAside);
+       if (fLaser->fFitCside)  vdriftC=*(fLaser->fFitCside);
+       if (fLaser->fFitACside) vdriftAC=*(fLaser->fFitACside);
+      }
+      (*cstream)<<"timeInfo"<<
+       "run="<<fRun<<              //  run number
+       "event="<<fEvent<<          //  event number
+       "time="<<fTime<<            //  time stamp of event
+       "trigger="<<fTrigger<<      //  trigger
+       "mag="<<fMagF<<             //  magnetic field
+       // Environment values
+       "press0="<<valuePressure0<<
+       "press1="<<valuePressure1<<
+       "pt0="<<ptrelative0<<
+       "pt1="<<ptrelative1<<
+       "temp0="<<temp0<<
+       "temp1="<<temp1<<
+       //
+       // accumulated values
+       //
+       "fDz="<<fDz<<          //! current delta z
+       "fdEdx="<<fdEdx<<        //! current dEdx
+       "fdEdxRatio="<<fdEdxRatio<<   //! current dEdx ratio
+       "fTl="<<fTl<<          //! current tan(lambda)
+       //
+       //laser
+       //
+       "laserA.="<<&vdriftA<<
+       "laserC.="<<&vdriftC<<
+       "laserAC.="<<&vdriftAC<<
+       "\n";
+    }
+  }
 }
 
 
@@ -218,7 +305,10 @@ void AliTPCcalibTime::ProcessCosmic(AliESDEvent *event) {
        Double_t meanP = 0.5*(trackIn->GetP() + trackOut->GetP());
        Double_t TPCsignal = seed->CookdEdxNorm(0.0,0.6,1,0,159,0x0,kTRUE,kTRUE);
        Double_t vecDeDx[2] = {time, TPCsignal};
-       if (meanP > 12) fHistDeDx->Fill(vecDeDx);
+       if (meanP > 12) {
+        fdEdx = TPCsignal;
+        fHistDeDx->Fill(vecDeDx);
+       }
      }
    }
   }
@@ -302,12 +392,18 @@ void AliTPCcalibTime::ProcessCosmic(AliESDEvent *event) {
        if (isPair && z0>0 && z0inner>0 && z1<0 && z1inner<0) {
         Double_t vecVdrift[2] = {time, param0.GetZ() - param1.GetZ()};
         
-        if (track0->GetTPCNcls() > 80) fHistVdrift->Fill(vecVdrift);
+        if (track0->GetTPCNcls() > 80) {
+          fHistVdrift->Fill(vecVdrift);
+          fDz = param0.GetZ() - param1.GetZ();
+        }
        }
        if (isPair && z0 > 0 && z1 > 0) {
         if (track1->GetTPCNcls()> 110 && track0->GetTPCNcls()> 110 && seed1->CookdEdxNorm(0,0.6,1,0,159,0,kFALSE,kTRUE) > 0) {
+          Float_t dedxratio = seed0->CookdEdxNorm(0,0.6,1,0,159,0,kFALSE,kTRUE)/seed1->CookdEdxNorm(0,0.6,1,0,159,0,kFALSE,kTRUE);
           Double_t vecDeDxTgl[3] = {time, track0->GetTgl(), seed0->CookdEdxNorm(0,0.6,1,0,159,0,kFALSE,kTRUE)/seed1->CookdEdxNorm(0,0.6,1,0,159,0,kFALSE,kTRUE)};
           fHistDeDxTgl->Fill(vecDeDxTgl);
+          fdEdxRatio = dedxratio;
+          fTl = track0->GetTgl() ;
         }
        }
        
@@ -371,6 +467,7 @@ Bool_t  AliTPCcalibTime::IsPair(AliExternalTrackParam *tr0, AliExternalTrackPara
   const Double_t *p1 = tr1->GetParameter();
   if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
   if (TMath::Abs(p0[0]+p1[0])>fCutMaxD)  return kFALSE;
+  if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz)  return kFALSE;
   Double_t d0[3], d1[3];
   tr0->GetDirection(d0);    
   tr1->GetDirection(d1);       
index 4a5396781a6aff987c00c92b09bc878e25c8104e..41a34d5f1b4cb092eedf6cba7404b079d768e4cc 100644 (file)
@@ -17,6 +17,7 @@ class THnSparse;
 class TList;
 class AliESDEvent;
 class AliESDtrack;
+class AliTPCcalibLaser;
 
 #include "TTreeStream.h"
 
@@ -41,7 +42,7 @@ public:
   
 
 private:
-
+  void ResetCurrent();                  // reset current values
   ULong64_t fTriggerMask;               // select certain trigger within one run
 
   THnSparse * fHistDeDxTgl;             // dEdx vs. dip angle vs time histogram
@@ -51,9 +52,19 @@ private:
   Float_t fIntegrationTimeDeDx;         // required statistics for each dEdx time bin
   Float_t fIntegrationTimeVdrift;       // required statistics for each Vdrift time bin
 
+  AliTPCcalibLaser * fLaser;            //! laser calibration
+  //
+  // current information
+  //
+  Float_t fDz;          //! current delta z
+  Float_t fdEdx;        //! current dEdx
+  Float_t fdEdxRatio;   //! current dEdx ratio
+  Float_t fTl;          //! current tan(lambda)
+  
   // cuts
   //
   Float_t fCutMaxD;     // maximal distance in rfi ditection
+  Float_t fCutMaxDz;     // maximal distance in z ditection
   Float_t fCutTheta;    // maximal distance in theta ditection
   Float_t fCutMinDir;   // direction vector products