]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Updates for the TPC calibration using laser:
authormarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 2 Jul 2009 16:49:21 +0000 (16:49 +0000)
committermarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 2 Jul 2009 16:49:21 +0000 (16:49 +0000)
(Peter Christiansen)

I made some modification to the laser/TPC code to have
1) absolute residuals in AliTPCcalibLaser
2) be able to switch on/off the corrections I want.

The modified files are:
M      AliTPCcalibLaser.h    (for 1)
M      AliTPCcalibLaser.cxx  (for 1)
M      AliTPCcalibCalib.h    (for 2)
M      AliTPCcalibCalib.cxx  (for 2)
M      AliTPCRecoParam.h     (for 2)
M      AliTPCRecoParam.cxx   (for 2)
M      AliTPCTransform.h     (for 2)
M      AliTPCTransform.cxx   (for 2)

As can be seen most of the modifications are for 2 (but they are of
course small).
There is no dependence between the two so the files for 1 (most
important changes) can be committed without 2.

Here follows some more details about the modifications:
*******************************************************************
1) To have histograms with absolute residuals in the laser calibration
output:
AliTPCcalibLaser
Added 2 array of histograms and 1 method to calculate the linear
parametrisations based on the database.

1 problem was solved in a bit non-optimal way. In the laser database the
slopes in z are inverted, so I now as default invert them and has a flag
that can be set not to invert them for comparison with the simulated
data. Ideal solution is to update laser data base and remove this fix
(and regenerate the laser simulation).

I only added functionality. I did not remove any. (this also means that
I do not cut away rows where I know from the simulation that the beams
are overlapping as I have done in my previous presentations).
********************************************************************
2) Give the possibility to switch on and off corrections.
- AliTPCRecoParam.h
   Added the variable
   Int_t fUseExBCorrection;       // use ExB correction
   to be able to not apply ExB correction to the data
   default is 1 = on
NB!!! I do not know if AliTPCRecoParam is in the database as an object?
Could it cause problems that I added a variable? In principle it should
be set correctly to 1 since I bumped the version number....

In AliTPCTransform.cxx I made something similar to what was done for the
time correction:
+  if(fCurrentRecoParam&&fCurrentRecoParam->GetUseExBCorrection()) {
+
+    calib->GetExB()->Correct(x,xx);
+
+  } else {
+
+    xx[0] = x[0];
+    xx[1] = x[1];
+    xx[2] = x[2];
+  }

NB! This could cause trouble in the reconstruction in the case where no
fCurrentRecoParam is there because then no ExB correction will be done.
But as this is the same situation for the time correction I assumed that
the aim is that there will always be a fCurrentRecoParam set.

To control all the corrections done in AliTPCcalibCalib I made a flag
for each correction that can be set. For the ExB correction I had to add
a method to AliTPCTransform.h
 AliTPCRecoParam * GetCurrentRecoParamNonConst() const {return
fCurrentRecoParam;}
so that I can set the flag to not use ExB corrections.

This means that I can do what I want to do:
Run laser calibrations with my choice of calibrations applied
Get the absolute residuals in the output.

Here is a short presentation with the results I get:
http://www.hep.lu.se/staff/christiansen/laser_absres_in_alitpccaliblaser.pdf
http://www.hep.lu.se/staff/christiansen/laser_absres_in_alitpccaliblaser.odp
********************************************************************
I hope that part 1 can be committed to AliRoot. I would however hope
that if part 2 is not committed that in some other way these switches
are added that the standard AliRoot can be configured to apply only
certain corrections.

TPC/AliTPCcalibCalib.cxx
TPC/AliTPCcalibCalib.h
TPC/AliTPCcalibLaser.cxx
TPC/AliTPCcalibLaser.h

index b7ea0f8d72442dec22d4dd739aca956b186441b3..e88a285c431fab36a3e8be47a4004e7094934a9a 100644 (file)
@@ -58,6 +58,7 @@
 
 #include "AliTPCcalibDB.h"
 #include "AliTPCTransform.h"
+#include "AliTPCRecoParam.h"
 #include "AliTPCclusterMI.h"
 #include "AliTPCseed.h"
 #include "AliTPCPointCorrection.h"
 ClassImp(AliTPCcalibCalib)
 
 AliTPCcalibCalib::AliTPCcalibCalib():
-    AliTPCcalibBase()
+AliTPCcalibBase(),
+  fApplyExBCorrection(1),      // apply ExB correction
+  fApplyTOFCorrection(1),      // apply TOF correction
+  fApplyPositionCorrection(1), // apply position correction
+  fApplySectorAlignment(1),    // apply sector alignment
+  fApplyRPhiCorrection(1),     // apply R-Phi correction
+  fApplyRCorrection(1)         // apply Radial correction
+
 {
   //
   // Constructor
@@ -74,7 +82,13 @@ AliTPCcalibCalib::AliTPCcalibCalib():
 
 
 AliTPCcalibCalib::AliTPCcalibCalib(const Text_t *name, const Text_t *title) 
-  :AliTPCcalibBase()
+  :AliTPCcalibBase(),
+  fApplyExBCorrection(1),      // apply ExB correction
+  fApplyTOFCorrection(1),      // apply TOF correction
+  fApplyPositionCorrection(1), // apply position correction
+  fApplySectorAlignment(1),    // apply sector alignment
+  fApplyRPhiCorrection(1),     // apply R-Phi correction
+  fApplyRCorrection(1)         // apply Radial correction
 {  
   SetName(name);
   SetTitle(title);
@@ -82,7 +96,14 @@ AliTPCcalibCalib::AliTPCcalibCalib(const Text_t *name, const Text_t *title)
 
 
 AliTPCcalibCalib::AliTPCcalibCalib(const AliTPCcalibCalib&calib):
-  AliTPCcalibBase(calib)
+  AliTPCcalibBase(calib),
+  fApplyExBCorrection(calib.GetApplyExBCorrection()),
+  fApplyTOFCorrection(calib.GetApplyTOFCorrection()),
+  fApplyPositionCorrection(calib.GetApplyPositionCorrection()),
+  fApplySectorAlignment(calib.GetApplySectorAlignment()),
+  fApplyRPhiCorrection(calib.GetApplyRPhiCorrection()),
+  fApplyRCorrection(calib.GetApplyRCorrection())
+
 {
   //
   // copy constructor
@@ -155,6 +176,15 @@ Bool_t  AliTPCcalibCalib::RefitTrack(AliESDtrack * track, AliTPCseed *seed){
   AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
   transform->SetCurrentRun(fRun);
   transform->SetCurrentTimeStamp((UInt_t)fTime);
+  if(!fApplyExBCorrection) { // disable ExB correction in transform
+    if(transform->GetCurrentRecoParam())
+      transform->GetCurrentRecoParamNonConst()->SetUseExBCorrection(0);
+  }
+  if(!fApplyTOFCorrection) { // disable TOF correction in transform
+    if(transform->GetCurrentRecoParam())
+      transform->GetCurrentRecoParamNonConst()->SetUseTOFCorrection(kFALSE);
+  }
+
   //
   // First apply calibration
   //
@@ -165,6 +195,7 @@ Bool_t  AliTPCcalibCalib::RefitTrack(AliESDtrack * track, AliTPCseed *seed){
     AliTPCclusterMI cl0(*cluster);
     Double_t x[3]={cluster->GetRow(),cluster->GetPad(),cluster->GetTimeBin()};
     Int_t i[1]={cluster->GetDetector()};
+
     transform->Transform(x,i,0,1);
     //
     // get position correction
@@ -173,8 +204,10 @@ Bool_t  AliTPCcalibCalib::RefitTrack(AliESDtrack * track, AliTPCseed *seed){
     if (cluster->GetDetector()>35) ipad=1;
     Float_t dy =AliTPCClusterParam::SPosCorrection(0,ipad,cluster->GetPad(),cluster->GetTimeBin(),cluster->GetZ(),cluster->GetSigmaY2(),cluster->GetSigmaZ2(),cluster->GetMax());
     Float_t dz =AliTPCClusterParam::SPosCorrection(1,ipad,cluster->GetPad(),cluster->GetTimeBin(),cluster->GetZ(),cluster->GetSigmaY2(),cluster->GetSigmaZ2(),cluster->GetMax());
+    // if(fApplyPositionCorrection) {
     //x[1]-=dy;
     //x[2]-=dz;
+    // }
     //
     // Apply sector alignment
     //
@@ -184,7 +217,7 @@ Bool_t  AliTPCcalibCalib::RefitTrack(AliESDtrack * track, AliTPCseed *seed){
                                                               cluster->GetY(),cluster->GetZ());
     Double_t dzq = AliTPCPointCorrection::SGetCorrectionSector(2,cluster->GetDetector()%36,cluster->GetX(),
                                                               cluster->GetY(),cluster->GetZ());
-    if (kTRUE){
+    if (fApplySectorAlignment){
       x[0]-=dxq;
       x[1]-=dyq;
       x[2]-=dzq;
@@ -198,9 +231,11 @@ Bool_t  AliTPCcalibCalib::RefitTrack(AliESDtrack * track, AliTPCseed *seed){
     // R correction
     Double_t corrR   = corr->CorrectionOutR0(kFALSE,kFALSE,cluster->GetX(),cluster->GetY(),cluster->GetZ(),cluster->GetDetector());
 
-    if (kTRUE){
+    if (fApplyRPhiCorrection){
       if (cluster->GetY()>0) x[1]+=corrclY;  // rphi correction
       if (cluster->GetY()<0) x[1]-=corrclY;  // rphi correction
+    }
+    if (fApplyRCorrection){      
       x[0]+=corrR;                           // radial correction
     }
 
index 012e7743e344f825b52ed2ddd8d324b1c3590857..437e5d288658f500856a0c1952c61bc12dbad59d 100644 (file)
@@ -31,9 +31,31 @@ public:
   Bool_t  RejectCluster(AliTPCclusterMI* cl, AliExternalTrackParam * param);
   void     Process(AliESDtrack *track, Int_t runNo=-1){AliTPCcalibBase::Process(track,runNo);};
   void     Process(AliTPCseed *track){return AliTPCcalibBase::Process(track);}
+
+  void  SetApplyExBCorrection(Int_t flag){fApplyExBCorrection=flag;}
+  void  SetApplyTOFCorrection(Int_t flag){fApplyTOFCorrection=flag;}
+  void  SetApplyPositionCorrection(Int_t flag){fApplyPositionCorrection=flag;}
+  void  SetApplySectorAlignment(Int_t flag){fApplySectorAlignment=flag;}
+  void  SetApplyRPhiCorrection(Int_t flag){fApplyRPhiCorrection=flag;}
+  void  SetApplyRCorrection(Int_t flag){fApplyRCorrection=flag;}
+
+  //
+  Int_t GetApplyExBCorrection() const {return fApplyExBCorrection;}
+  Int_t GetApplyTOFCorrection() const {return fApplyTOFCorrection;}
+  Int_t GetApplyPositionCorrection() const {return fApplyPositionCorrection;}
+  Int_t GetApplySectorAlignment() const {return fApplySectorAlignment;}
+  Int_t GetApplyRPhiCorrection() const {return fApplyRPhiCorrection;}
+  Int_t GetApplyRCorrection() const {return fApplyRCorrection;}
+
 protected: 
+  Int_t fApplyExBCorrection;      // apply ExB correction (in AliTPCTransform)
+  Int_t fApplyTOFCorrection;      // apply TOF correction (in AliTPCTransform)
+  Int_t fApplyPositionCorrection; // apply position correction
+  Int_t fApplySectorAlignment;    // apply sector alignment
+  Int_t fApplyRPhiCorrection;     // apply R-Phi correction
+  Int_t fApplyRCorrection;        // apply Radial correction
 private:
-  ClassDef(AliTPCcalibCalib,1)
+  ClassDef(AliTPCcalibCalib,2)
 };
 
 #endif
index fcc5a6af59c6b41cdc50175153cf789de9137e7e..5613ff3e99ee8fd00485a5df87f2dff66d0d9fbd 100644 (file)
@@ -190,6 +190,8 @@ AliTPCcalibLaser::AliTPCcalibLaser():
   fDeltaZres(336),   //->2D histo fo residuals
   fDeltaYres2(336),   //->2D histo of residuals
   fDeltaZres2(336),   //->2D histo fo residuals
+  fDeltaYresAbs(336),   //->2D histo of residuals
+  fDeltaZresAbs(336),   //->2D histo of residuals
   //fDeltaYres3(336),   //->2D histo of residuals
   //fDeltaZres3(336),   //->2D histo fo residuals
   fFitAside(new TVectorD(5)),
@@ -198,7 +200,18 @@ AliTPCcalibLaser::AliTPCcalibLaser():
   fEdgeXcuts(3),    
   fEdgeYcuts(3),    
   fNClCuts(5),      
-  fNcuts(0)        
+  fNcuts(0),
+  fBeamSectorOuter(336),
+  fBeamSectorInner(336),
+  fBeamOffsetYOuter(336),
+  fBeamSlopeYOuter(336),
+  fBeamOffsetYInner(336),
+  fBeamSlopeYInner(336),
+  fBeamOffsetZOuter(336),
+  fBeamSlopeZOuter(336),
+  fBeamOffsetZInner(336),
+  fBeamSlopeZInner(336),
+  fInverseSlopeZ(kTRUE)
 {
   //
   // Constructor
@@ -274,6 +287,8 @@ AliTPCcalibLaser::AliTPCcalibLaser(const Text_t *name, const Text_t *title, Bool
   fDeltaZres(336),  
   fDeltaYres2(336),
   fDeltaZres2(336),  
+  fDeltaYresAbs(336),
+  fDeltaZresAbs(336),
   //  fDeltaYres3(336),
   //fDeltaZres3(336),  
   fFitAside(new TVectorD(5)),        // drift fit - A side
@@ -282,7 +297,18 @@ AliTPCcalibLaser::AliTPCcalibLaser(const Text_t *name, const Text_t *title, Bool
   fEdgeXcuts(3),       // cuts in local x direction; used in the refit of the laser tracks
   fEdgeYcuts(3),       // cuts in local y direction; used in the refit of the laser tracks
   fNClCuts(5),         // cuts on the number of clusters per tracklet; used in the refit of the laser tracks
-  fNcuts(0)           // number of cuts
+  fNcuts(0),           // number of cuts
+  fBeamSectorOuter(336),
+  fBeamSectorInner(336),
+  fBeamOffsetYOuter(336),
+  fBeamSlopeYOuter(336),
+  fBeamOffsetYInner(336),
+  fBeamSlopeYInner(336),
+  fBeamOffsetZOuter(336),
+  fBeamSlopeZOuter(336),
+  fBeamOffsetZInner(336),
+  fBeamSlopeZInner(336),
+  fInverseSlopeZ(kTRUE)
 {
   SetName(name);
   SetTitle(title);
@@ -359,6 +385,8 @@ AliTPCcalibLaser::AliTPCcalibLaser(const AliTPCcalibLaser& calibLaser):
   fDeltaZres(((calibLaser.fDeltaZres))),  
   fDeltaYres2(((calibLaser.fDeltaYres))),
   fDeltaZres2(((calibLaser.fDeltaZres))),  
+  fDeltaYresAbs(((calibLaser.fDeltaYresAbs))),
+  fDeltaZresAbs(((calibLaser.fDeltaZresAbs))),
   //  fDeltaYres3(((calibLaser.fDeltaYres))),
   //fDeltaZres3(((calibLaser.fDeltaZres))),  
   fFitAside(new TVectorD(5)),        // drift fit - A side
@@ -367,7 +395,19 @@ AliTPCcalibLaser::AliTPCcalibLaser(const AliTPCcalibLaser& calibLaser):
   fEdgeXcuts(3),       // cuts in local x direction; used in the refit of the laser tracks
   fEdgeYcuts(3),       // cuts in local y direction; used in the refit of the laser tracks
   fNClCuts(5),         // cuts on the number of clusters per tracklet; used in the refit of the laser tracks
-  fNcuts(0)           // number of cuts
+  fNcuts(0),           // number of cuts
+  fBeamSectorOuter(336),
+  fBeamSectorInner(336),
+  fBeamOffsetYOuter(336),
+  fBeamSlopeYOuter(336),
+  fBeamOffsetYInner(336),
+  fBeamSlopeYInner(336),
+  fBeamOffsetZOuter(336),
+  fBeamSlopeZOuter(336),
+  fBeamOffsetZInner(336),
+  fBeamSlopeZInner(336),
+  fInverseSlopeZ(calibLaser.fInverseSlopeZ)
+
 {
   //
   // copy constructor
@@ -464,7 +504,10 @@ AliTPCcalibLaser::~AliTPCcalibLaser() {
   fDeltaZres2.SetOwner();
   fDeltaZres2.Delete();
   
-
+  fDeltaYresAbs.SetOwner();
+  fDeltaYresAbs.Delete();
+  fDeltaZresAbs.SetOwner();
+  fDeltaZresAbs.Delete();
 }
 
 
@@ -1612,32 +1655,87 @@ void AliTPCcalibLaser::RefitLaserJW(Int_t id){
       TH2F *profz = (TH2F*)fDeltaZres.UncheckedAt(id);
       TH2F *profy2 = (TH2F*)fDeltaYres2.UncheckedAt(id);
       TH2F *profz2 = (TH2F*)fDeltaZres2.UncheckedAt(id);
+      TH2F *profyabs = (TH2F*)fDeltaYresAbs.UncheckedAt(id);
+      TH2F *profzabs = (TH2F*)fDeltaZresAbs.UncheckedAt(id);
       //      TH2F *profy3 = (TH2F*)fDeltaYres3.UncheckedAt(id);
       //TH2F *profz3 = (TH2F*)fDeltaZres3.UncheckedAt(id);
       //
       for (Int_t irow=158;irow>-1;--irow) {
        if (vecSec[irow]==-1)continue; // no cluster info
        if (isReject[irow]>0.5) continue;  // 
-       //Double_t x    = vecX[irow];
        Double_t ycl  = vecClY[irow];
        Double_t yfit = vecY1[irow];
        Double_t yfit2 = vecY2[irow];
+       Double_t x    = vecX[irow];
+       Double_t yabsbeam = -1000;
+       if(vecSec[irow]==outerSector && (outerSector%36)==fBeamSectorOuter[id])
+         yabsbeam = fBeamSlopeYOuter[id]*x + fBeamOffsetYOuter[id];
+       else if((innerSector%36)==fBeamSectorInner[id])
+         yabsbeam = fBeamSlopeYInner[id]*x + fBeamOffsetYInner[id];
+
        //      Double_t yfit3 = vecY2[irow];
        Double_t zcl  = vecClZ[irow];
        Double_t zfit = vecZ1[irow];
        Double_t zfit2 = vecZ2[irow];
        //Double_t zfit3 = vecZ2[irow];
 
+       // dz abs 
+       // The expressions for zcorrected has been obtained by
+       // inverting the fits in the FitDriftV() method (ignoring the
+       // global y dependence for now):
+       // A side:
+       // 250 - zmeasured = [0] + [1]*(250 - zreal) + .... (yglobal)
+       // =>
+       // zreal = (zmeasured + [0] - (1-[1])*250.0)/[1]
+       // 
+       // C side:
+       // 250 + zmeasured  = [0] + [1]*(250+zreal) + .... (yglobal)
+       // =>
+       // zreal = (zmeasured - [0] + (1 - [1])*250.0)/[1]
+
+       Double_t dzabs = -1000;
+       if (ltrp->GetSide()==0){
+         if ((*fFitAside)[1]>0.) { 
+           // ignore global y dependence for now
+           Double_t zcorrected = (zcl + (*fFitAside)[0] - 
+                                  (1.0-(*fFitAside)[1])*250.0)/(*fFitAside)[1];
+           //      zcorrected = zcl;
+           if(vecSec[irow]==outerSector && (outerSector%36)==fBeamSectorOuter[id])
+             dzabs = zcorrected -fBeamSlopeZOuter[id]*x -fBeamOffsetZOuter[id];
+           else if((innerSector%36)==fBeamSectorInner[id])
+             dzabs = zcorrected -fBeamSlopeZInner[id]*x -fBeamOffsetZInner[id];
+         }
+       } else {
+         if ((*fFitCside)[1]>0.) {
+           Double_t zcorrected = (zcl - (*fFitCside)[0] + 
+                                  (1.0-(*fFitCside)[1])*250.0)/(*fFitCside)[1];
+           //      zcorrected = zcl;
+           if(vecSec[irow]==outerSector && (outerSector%36)==fBeamSectorOuter[id])
+             dzabs = zcorrected -fBeamSlopeZOuter[id]*x -fBeamOffsetZOuter[id];
+           else if((innerSector%36)==fBeamSectorInner[id])
+             dzabs = zcorrected -fBeamSlopeZInner[id]*x -fBeamOffsetZInner[id];
+         }
+       }
+       
        if (TMath::Abs(yfit-ycl)<2&&TMath::Abs(zfit-zcl)<2){
          if (profy){       
              profy->Fill(irow,ycl-yfit);
              profy2->Fill(irow,ycl-yfit2);
+             if(yabsbeam<-100)
+               profyabs->Fill(irow,-0.99);
+             else
+               profyabs->Fill(irow,ycl-yabsbeam);
+
              //              profy3->Fill(irow,ycl-yfit3);
          }
          if (profz) {
              profz->Fill(irow,zcl-zfit);
              profz2->Fill(irow,zcl-zfit2);
              //profz3->Fill(irow,zcl-zfit3);
+             if(dzabs<-100)
+               profzabs->Fill(irow,-0.99);
+             else
+               profzabs->Fill(irow,dzabs);
          }
        }
       }
@@ -2384,6 +2482,15 @@ Long64_t AliTPCcalibLaser::Merge(TCollection *li) {
       h2m = (TH2F*)cal->fDeltaZres2.At(id);
       h2  = (TH2F*)fDeltaZres2.At(id);
       if (h2m&&h2) h2->Add(h2m);
+
+      // merge ProfileY histograms - abs
+      h2m = (TH2F*)cal->fDeltaYresAbs.At(id);
+      h2  = (TH2F*)fDeltaYresAbs.At(id);
+      if (h2m&&h2) h2->Add(h2m);
+
+      h2m = (TH2F*)cal->fDeltaZresAbs.At(id);
+      h2  = (TH2F*)fDeltaZresAbs.At(id);
+      if (h2m&&h2) h2->Add(h2m);
       // merge ProfileY histograms - 3
       //h2m = (TH2F*)cal->fDeltaYres3.At(id);
       //h2  = (TH2F*)fDeltaYres3.At(id);
@@ -2538,6 +2645,8 @@ void   AliTPCcalibLaser::MakeFitHistos(){
     TH2F *profz = (TH2F*)fDeltaZres.UncheckedAt(id);
     TH2F *profy2 = (TH2F*)fDeltaYres2.UncheckedAt(id);
     TH2F *profz2 = (TH2F*)fDeltaZres2.UncheckedAt(id);
+    TH2F *profyabs = (TH2F*)fDeltaYresAbs.UncheckedAt(id);
+    TH2F *profzabs = (TH2F*)fDeltaYresAbs.UncheckedAt(id);
     //    TH2F *profy3 = (TH2F*)fDeltaYres3.UncheckedAt(id);
     //TH2F *profz3 = (TH2F*)fDeltaZres3.UncheckedAt(id);
     if (!profy){
@@ -2547,6 +2656,9 @@ void   AliTPCcalibLaser::MakeFitHistos(){
       profy2=new TH2F(Form("pry%03d",id),Form("Y Residuals for Laser Beam %03d -Parabolic",id),160,0,160,50,-0.5,0.5);
       profy2->SetDirectory(0);
       fDeltaYres2.AddAt(profy2,id);
+      profyabs=new TH2F(Form("pryabs%03d",id),Form("Y Residuals for Laser Beam %03d -Absolute",id),160,0,160,100,-1.0,1.0); // has to be bigger based on earlier studies
+      profyabs->SetDirectory(0);
+      fDeltaYresAbs.AddAt(profyabs,id);
       //profy3=new TH2F(Form("pry%03d",id),Form("Y Residuals for Laser Beam %03d- Parabolic2",id),160,0,160,100,-0.5,0.5);
       //profy3->SetDirectory(0);
       //fDeltaYres3.AddAt(profy3,id);
@@ -2558,6 +2670,9 @@ void   AliTPCcalibLaser::MakeFitHistos(){
       profz2=new TH2F(Form("prz%03d",id),Form("Z Residuals for Laser Beam %03d - Parabolic",id),160,0,160,50,-0.5,0.5);
       profz2->SetDirectory(0);
       fDeltaZres2.AddAt(profz2,id);
+      profzabs=new TH2F(Form("przabs%03d",id),Form("Z Residuals for Laser Beam %03d -Absolute",id),160,0,160,100,-1.0,1.0); // has to be bigger based on earlier studies
+      profzabs->SetDirectory(0);
+      fDeltaZresAbs.AddAt(profzabs,id);
       //profz3=new TH2F(Form("prz%03d",id),Form("Z Residuals for Laser Beam %03d- Parabolic2",id),160,0,160,100,-0.5,0.5);
       //profz3->SetDirectory(0);
       //fDeltaZres3.AddAt(profz3,id);
@@ -2600,6 +2715,11 @@ void   AliTPCcalibLaser::MakeFitHistos(){
       fSignals.AddAt(hisSignal,id);
     }
   }
+
+  SetBeamParameters(fBeamOffsetZOuter, fBeamSlopeZOuter, fBeamSectorOuter,2);
+  SetBeamParameters(fBeamOffsetZInner, fBeamSlopeZInner, fBeamSectorInner,3);
+  SetBeamParameters(fBeamOffsetYOuter, fBeamSlopeYOuter, fBeamSectorOuter,0);
+  SetBeamParameters(fBeamOffsetYInner, fBeamSlopeYInner, fBeamSectorInner,1);
 }
 
 void AliTPCcalibLaser::MergeFitHistos(AliTPCcalibLaser * laser){
@@ -2932,6 +3052,83 @@ void AliTPCcalibLaser::DumpFitInfo(TTree * chainFit,Int_t id){
     
     delete pcstream;
   }
+}
+
+void AliTPCcalibLaser::SetBeamParameters(TVectorD& meanOffset,
+                                        TVectorD& meanSlope,
+                                        TVectorD& sectorArray,
+                                        Int_t option)
+{
+  // This method should ideally go in AliTPCLaser
+  // option == 0 (pads outer - closest to beam)
+  // option == 1 (pads inner)
+  // option == 2 (time outer)
+  // option == 3 (time inner)
+  Int_t nFailures = 0;
+
+  for(Int_t id = 0; id < 336; id++) {
+    
+    if (!AliTPCLaserTrack::GetTracks()) 
+      AliTPCLaserTrack::LoadTracks();
+    AliTPCLaserTrack *ltrp =
+      (AliTPCLaserTrack*)AliTPCLaserTrack::GetTracks()->UncheckedAt(id);
+    
+    AliExternalTrackParam trackParam(*ltrp);
+    
+    Double_t deltaangle = 10.0; // sector is 1/2 a sector away from mirror
+    if((option==1 || option==3)&& (ltrp->GetBeam()<=1 || ltrp->GetBeam()>=5))
+      deltaangle = 30.0; // inner sector is 1 sector further away from mirror
+    
+    Double_t angle = trackParam.GetAlpha();
+    if(angle<0)
+      angle += 2*TMath::Pi();
+    if(trackParam.GetSnp()>0) // track points to sector "before"
+      angle -= deltaangle*TMath::DegToRad();
+    else // track points to sector "after"
+      angle += deltaangle*TMath::DegToRad();
+    
+    Bool_t success = trackParam.Rotate(angle);
+    
+    if(!success) {
+      //      cout << "WARNING: Rotate failed for ID: " << id << endl;
+      nFailures++;
+    }
+    
+    angle *= TMath::RadToDeg();
+    
+    Int_t sector = TMath::Nint((angle-10.0)/20.0);
+    if(sector<0)
+      sector += 18;
+    else if(sector>=18)
+      sector -= 18;
+    if(ltrp->GetSide()==1) // C side
+      sector += 18;
+    sectorArray[id] = sector;
+
+    const Double_t x0 = 0;
+    
+    Double_t slopey  = TMath::Tan(TMath::ASin(trackParam.GetSnp()));
+    Double_t slopez  = trackParam.GetTgl();
+    //  One needs a factor sqrt(1+slopey**2) to take into account the
+    //  longer path length
+    slopez *= TMath::Sqrt(1.0 + slopey*slopey);    
+    if(fInverseSlopeZ) // wrong sign in database, should be fixed there
+      slopez *= -1;
+    //    Double_t offsetz = trackParam.GetZ();
+    Double_t offsety = trackParam.GetY() + slopey*(x0-trackParam.GetX());
+    Double_t offsetz = trackParam.GetZ() + slopez*(x0-trackParam.GetX());
+    if(option==2 || option==3) {
+      meanOffset[id] = offsetz; meanSlope[id] = slopez;
+    } else {
+      meanOffset[id] = offsety; meanSlope[id] = slopey;
+    }
+  }
+
+  if(nFailures>0)
+    AliWarning(Form("Rotate method failed %d times", nFailures));
+}
+
+
   /*
     TFile f("vscan.root");
    */  
@@ -3020,8 +3217,6 @@ TH2 * phisPy2In = (TH2*) gROOT->FindObject("hisPy2In0")
 
 */
 
-}  
-
 
 
 
index 1023080bb15d4f6c473bf24ca1bcc393398e7d74..728fbf3f65eea520b02337401ba6caa8f8424a1a 100644 (file)
@@ -57,6 +57,8 @@ public:
   void     Process(AliESDtrack *track, Int_t runNo=-1){AliTPCcalibBase::Process(track,runNo);};
   void     Process(AliTPCseed *track){return AliTPCcalibBase::Process(track);}
   //
+  void SetBeamParameters(TVectorD& meanOffset, TVectorD& meanSlope,
+                        TVectorD& sectorArray, Int_t option);
 
   AliESDEvent  * fESD;             //! ESD event  - not OWNER
   AliESDfriend * fESDfriend;       //! ESD event  - not OWNER
@@ -80,6 +82,7 @@ public:
   TObjArray      fDeltaPhi;        //-> array of histograms of delta Phi for each track
   TObjArray      fDeltaPhiP;       //-> array of histograms of delta Phi direction for each track
   TObjArray      fSignals;         //->Array of dedx signals
+
   //
   // Refit residuals histogram
   //
@@ -133,6 +136,8 @@ public:
   TObjArray      fDeltaZres;       //-> array of histograms of delta z residuals for each track
   TObjArray      fDeltaYres2;       //-> array of histograms of delta y residuals for each track
   TObjArray      fDeltaZres2;       //-> array of histograms of delta z residuals for each track
+  TObjArray      fDeltaYresAbs; //-> array of histograms of absolute delta y residuals for each track
+  TObjArray      fDeltaZresAbs; //-> array of histograms of absolute delta z residuals for each track
   //  TObjArray      fDeltaYres3;       //-> array of histograms of delta y residuals for each track
   //TObjArray      fDeltaZres3;       //-> array of histograms of delta z residuals for each track
 
@@ -145,6 +150,17 @@ public:
   TVectorD       fEdgeYcuts;       //! cuts in local y direction; used in the refit of the laser tracks
   TVectorD       fNClCuts;         //! cuts on the number of clusters per tracklet; used in the refit of the laser tracks
   Int_t          fNcuts;           //! number of cuts
+  TVectorD       fBeamSectorOuter;  //! sector map for beams in outer sector
+  TVectorD       fBeamSectorInner;  //! sector map for beams in inner sector
+  TVectorD       fBeamOffsetYOuter; //! absolute y beam offset in outer sector
+  TVectorD       fBeamSlopeYOuter;  //! absolute y beam slope  in outer sector
+  TVectorD       fBeamOffsetYInner; //! absolute y beam offset in inner sector
+  TVectorD       fBeamSlopeYInner;  //! absolute y beam slope  in inner sector
+  TVectorD       fBeamOffsetZOuter; //! absolute z beam offset in outer sectror 
+  TVectorD       fBeamSlopeZOuter;  //! absolute z beam slope  in outer sector 
+  TVectorD       fBeamOffsetZInner; //! absolute z beam offset in inner sectror 
+  TVectorD       fBeamSlopeZInner;  //! absolute z beam slope  in inner sector
+  Bool_t         fInverseSlopeZ;    //! invert slope in z - mismatch between database and lasers
   //
 private:
   ClassDef(AliTPCcalibLaser,3)