]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Possibility to use TGeoMatrix functionality
authormarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 13 Sep 2010 22:23:54 +0000 (22:23 +0000)
committermarian <marian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 13 Sep 2010 22:23:54 +0000 (22:23 +0000)
Convinient for fitting of the laignment parameters
(Marian)

TPC/AliTPCCalibGlobalMisalignment.cxx
TPC/AliTPCCalibGlobalMisalignment.h

index f7b2049e01112ea18413f660c4d473ce9cbffc3c..a2cb18e91c2c998775eefa4a85f9b8b2fea7b703 100644 (file)
@@ -20,6 +20,8 @@
 // The class calculates the space point distortions due to simple         // 
 // misalignments like shifts in caresian coordinates or a rotation        //
 // of the TPC read out planes (A and C side)                              //
+// Optionaly possible to use it for visualization of the alignemnt form the Alignment OCDB //
+// fUseGeoManager has to be set to kTRUE to enable this option
 //                                                                        //
 // date: 06/05/2010                                                       //
 // Authors: Stefan Rossegger, Jim Thomas, Magnus Mager                    //
 
 #include "AliTPCCalibGlobalMisalignment.h"
 #include "TMath.h"
+#include "TGeoMatrix.h"
+#include "AliTPCROC.h"
+#include "AliTPCcalibDB.h"
+#include "AliTPCParam.h"
+#include <TGeoPhysicalNode.h>
+
 
 AliTPCCalibGlobalMisalignment::AliTPCCalibGlobalMisalignment()
   : AliTPCCorrection("mialign","Misalignment"),
     fXShift(0.),fYShift(0.),fZShift(0.),
     fRotPhiA(0.),fRotPhiC(0.),
-    fdRPhiOffsetA(0.), fdRPhiOffsetC(0.)
+    fdRPhiOffsetA(0.), fdRPhiOffsetC(0.), 
+    fQuadrantDQ1(0), fQuadrantDQ2(0), fQuadrantQ2(0), fMatrixGlobal(0),
+    fMatrixASide(0), fMatrixCSide(0),
+    fUseGeomanager(kFALSE)
 {
   //
   // default constructor
@@ -43,8 +54,30 @@ AliTPCCalibGlobalMisalignment::~AliTPCCalibGlobalMisalignment() {
   //
   // default destructor
   //
+  delete fQuadrantDQ1;   //OROC medium pads delta ly+ - ly-
+  delete fQuadrantDQ2;   //OROC long   pads delta ly+ - ly-
+  delete fQuadrantQ2;    //OROC long   pads - OROC medium pads
+
+}
+void AliTPCCalibGlobalMisalignment::SetQuadranAlign(const TVectorD *dq1, const TVectorD *dq2, const TVectorD *q2){
+  //
+  // Set quadrant alignment
+  // 3 vectors for 36 (super) sectors
+  //
+  fQuadrantDQ1 = new TVectorD(*dq1);    //OROC medium pads delta ly+ - ly-
+  fQuadrantDQ2 = new TVectorD(*dq2);;   //OROC long   pads delta ly+ - ly-
+  fQuadrantQ2  = new TVectorD(*q2);;    //OROC long   pads - OROC medium pads
 }
 
+void AliTPCCalibGlobalMisalignment::SetGlobalAlign(const TGeoMatrix * matrixGlobal, const TGeoMatrix *matrixA, const TGeoMatrix *matrixC ){
+  //
+  // Set global misalignment as TGeoMatrix
+  // 
+  if (matrixGlobal) fMatrixGlobal = new TGeoHMatrix(*matrixGlobal);
+  if (matrixA) fMatrixASide = new TGeoHMatrix(*matrixA);
+  if (matrixC) fMatrixCSide = new TGeoHMatrix(*matrixC);
+}
 
 
 //void AliTPCCalibGlobalMisalignment::Init() {
@@ -70,25 +103,59 @@ AliTPCCalibGlobalMisalignment::~AliTPCCalibGlobalMisalignment() {
 void AliTPCCalibGlobalMisalignment::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
   //
   // Calculates the simple correction due to a shift (in x,y,z) or an rotation of the TPC (around z)
-  // 
-  Double_t r, phi;
+  //  
+  static AliTPCROC *tpcRoc =AliTPCROC::Instance();  
+  Double_t r=0, phi=0;
   r   = TMath::Sqrt( x[0]*x[0] + x[1]*x[1] );
   phi = TMath::ATan2(x[1],x[0]);
+  // Getsector number
+  Double_t sec=TMath::Nint(-0.5+(phi*9./TMath::Pi()));
+  if (sec<0) sec+=18;
+  Int_t isec = TMath::Nint(sec);    
+  if  (roc%36>=18) isec+=18;
+  //
+  // Get the point on the local coordiante frame
+  //
+  Double_t alpha=(sec+0.5)*TMath::Pi()/9;
+  Double_t pos[3]={0,0,x[2]};
+  pos[0]=  TMath::Cos(alpha)*x[0]+TMath::Sin(alpha)*x[1];
+  pos[1]= -TMath::Sin(alpha)*x[0]+TMath::Cos(alpha)*x[1];
+  if (pos[0]>tpcRoc->GetPadRowRadiiUp(0))  isec+=36;
 
-
+  //
+  // apply quadrant alignment if available - in local coordinate frame
+  //
+  Double_t posQG[3]={x[0],x[1],x[2]};
+  if (fQuadrantDQ1){
+    Double_t dly=0;
+    Bool_t isQ1 = TMath::Abs(pos[0]-161)<28;
+    Bool_t isQ2 = (pos[0]>189);
+    if (isQ1){
+      if (pos[1]>0.) dly+=(*fQuadrantDQ1)[isec];
+      if (pos[1]<0.) dly-=(*fQuadrantDQ1)[isec];      
+    }
+    if (isQ2){
+      dly+=(*fQuadrantQ2)[isec];
+      if (pos[1]>0.) dly+=(*fQuadrantDQ2)[isec];
+      if (pos[1]<0.) dly-=(*fQuadrantDQ2)[isec];
+    }
+    // Tranform the corrected point to the global frame
+    posQG[0]=  TMath::Cos(alpha)*pos[0]-TMath::Sin(alpha)*(pos[1]+dly);
+    posQG[1]=  TMath::Sin(alpha)*pos[0]+TMath::Cos(alpha)*(pos[1]+dly);
+  }
+  //
   // rotation of the read-out planes
   if  (roc%36<18) // A side
     phi += fRotPhiA;
   else         // C side
     phi += fRotPhiC;
-
+  
   // Simply adding a constant dRPHi residual. PURELY FOR CALIBRATION PURPOSES
   if  (roc%36<18) // A side
     phi += fdRPhiOffsetA/r;
   else         // C side
     phi += fdRPhiOffsetC/r;
-
+  
   dx[0] = r * TMath::Cos(phi) - x[0];
   dx[1] = r * TMath::Sin(phi) - x[1]; 
   dx[2] = 0.; 
@@ -97,7 +164,61 @@ void AliTPCCalibGlobalMisalignment::GetCorrection(const Float_t x[],const Short_
   dx[0] -= fXShift;
   dx[1] -= fYShift;
   dx[2] -= fZShift;
-
+  // quadrant shifts
+  dx[0] += (posQG[0]-x[0]);
+  dx[1] += (posQG[1]-x[1]);
+  //
+  // alignment matrix in local frame
+  //
+  if (fUseGeomanager){ //loading from the OCDB
+    Double_t posC[3] ={pos[0],pos[1],pos[2]};
+    //
+    //2. correct the point in the local frame
+    AliTPCParam *param = AliTPCcalibDB::Instance()->GetParameters();
+    if (!param){
+      //AliFatal("OCDB not initialized");
+    }
+    TGeoHMatrix  *mat = param->GetClusterMatrix(isec);
+    //
+    if (mat) mat->LocalToMaster(pos,posC);
+    Double_t posCG[3]={posC[0],posC[1],posC[2]};
+    //3. tranform the corrected point to the global frame
+    posCG[0]=  TMath::Cos(alpha)*posC[0]-TMath::Sin(alpha)*posC[1];
+    posCG[1]=  TMath::Sin(alpha)*posC[0]+TMath::Cos(alpha)*posC[1];
+    posCG[2]=  posC[2];
+    //4. Add delta
+    dx[0]+=posCG[0]-x[0];
+    dx[1]+=posCG[1]-x[1];
+    dx[2]+=posCG[2]-x[2];
+  }
+  if (fMatrixGlobal){
+    // apply global alignment matrix
+    Double_t ppos[3]={x[0],x[1],x[2]};
+    Double_t pposC[3]={x[0],x[1],x[2]};
+    fMatrixGlobal->LocalToMaster(ppos,pposC);
+    dx[0]+=pposC[0]-ppos[0];
+    dx[1]+=pposC[1]-ppos[1];
+    dx[2]+=pposC[2]-ppos[2];
+  }
+
+  if (fMatrixASide && roc%36<18){
+    // apply global alignment matrix
+    Double_t ppos[3]={x[0],x[1],x[2]};
+    Double_t pposC[3]={x[0],x[1],x[2]};
+    fMatrixASide->LocalToMaster(ppos,pposC);
+    dx[0]+=pposC[0]-ppos[0];
+    dx[1]+=pposC[1]-ppos[1];
+    dx[2]+=pposC[2]-ppos[2];
+  }
+  if (fMatrixCSide && roc%36>=18){
+    // apply global alignment matrix
+    Double_t ppos[3]={x[0],x[1],x[2]};
+    Double_t pposC[3]={x[0],x[1],x[2]};
+    fMatrixCSide->LocalToMaster(ppos,pposC);
+    dx[0]+=pposC[0]-ppos[0];
+    dx[1]+=pposC[1]-ppos[1];
+    dx[2]+=pposC[2]-ppos[2];
+  }
 }
 
 void AliTPCCalibGlobalMisalignment::Print(Option_t* /*option*/ ) const {
index 12f3056114d6a810b4561d07233efb62ffeef15c..0e2c30fe68efe932a58931320d92047471e411bc 100644 (file)
@@ -16,6 +16,9 @@
 ////////////////////////////////////////////////////////////////////////////
 
 #include "AliTPCCorrection.h"
+#include "TVectorD.h"
+class TGeoMatrix;
+
 
 class AliTPCCalibGlobalMisalignment : public AliTPCCorrection {
 public:
@@ -34,7 +37,8 @@ public:
   void SetRotPhiC(Float_t rotPhiC) {fRotPhiC=rotPhiC;}
   void SetdRPhiOffsetA(Float_t dRPhiOffsetA) {fdRPhiOffsetA=dRPhiOffsetA;}
   void SetdRPhiOffsetC(Float_t dRPhiOffsetC) {fdRPhiOffsetC=dRPhiOffsetC;}
-
+  void SetUseGeoManager(Bool_t useGeomanager) {fUseGeomanager = useGeomanager;}
+  
   Float_t GetXShift() const {return fXShift;}
   Float_t GetYShift() const {return fYShift;}
   Float_t GetZShift() const {return fZShift;}
@@ -42,9 +46,10 @@ public:
   Float_t GetRotPhiC() const {return fRotPhiC;}
   Float_t GetdRPhiOffsetA() const {return fdRPhiOffsetA;}
   Float_t GetdRPhiOffsetC() const {return fdRPhiOffsetC;}
-
+  Bool_t  GetUseGeoManager() const { return fUseGeomanager;}
   virtual void Print(Option_t* option="") const;
-
+  void SetQuadranAlign(const TVectorD *dq1, const TVectorD *dq2, const TVectorD *q2); 
+  void SetGlobalAlign(const TGeoMatrix * matrixGlobal, const TGeoMatrix *matrixA, const TGeoMatrix *matrixC );
 protected:
   virtual void GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]);
 
@@ -57,8 +62,23 @@ private:
   Float_t fRotPhiC;      // simple rotation of C side read-out plane around the Z axis [rad]
   Float_t fdRPhiOffsetA;  // add a constant offset of dRPhi (or local Y) in [cm]: purely for calibration purposes!
   Float_t fdRPhiOffsetC;  // add a constant offset of dRPhi (or local Y) in [cm]: purely for calibration purposes!
-
-  ClassDef(AliTPCCalibGlobalMisalignment,1);
+  //
+  // Quadrant alignment
+  //
+  TVectorD *fQuadrantDQ1;   //OROC medium pads delta ly+ - ly-
+  TVectorD *fQuadrantDQ2;   //OROC long   pads delta ly+ - ly-
+  TVectorD *fQuadrantQ2;   //OROC long   pads - OROC medium pads
+  //
+  // Global alignment - use native ROOT representation
+  //
+  TGeoMatrix * fMatrixGlobal; // global Alignment common
+  TGeoMatrix * fMatrixASide; // global Alignment A side
+  TGeoMatrix * fMatrixCSide; // global Alignment C side
+  //
+  Bool_t  fUseGeomanager;  // switch to use GeoManager - for visualization purposes
+  AliTPCCalibGlobalMisalignment& operator=(const AliTPCCalibGlobalMisalignment&);
+  AliTPCCalibGlobalMisalignment(const AliTPCCalibGlobalMisalignment&);
+  ClassDef(AliTPCCalibGlobalMisalignment,2);
 };
 
 #endif