]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Added Lorentz angle correction
authorshahoian <shahoian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 15 Nov 2012 02:56:58 +0000 (02:56 +0000)
committershahoian <shahoian@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 15 Nov 2012 02:56:58 +0000 (02:56 +0000)
ITS/UPGRADE/AliITSUClusterizer.cxx
ITS/UPGRADE/AliITSUClusterizer.h
ITS/UPGRADE/AliITSUReconstructor.cxx
ITS/UPGRADE/AliITSUReconstructor.h

index a02fbf70344504bd59634381e141ee9b9f202077..d5e9da94569b4c48acb9168c0b5438db359cdd81 100644 (file)
@@ -8,7 +8,9 @@
 #include "AliITSUSegmentationPix.h"
 #include "AliITSdigit.h"
 #include "AliITSURecoParam.h"
+#include "AliITSUAux.h"
 using namespace TMath;
+using namespace AliITSUAux;
 
 ClassImp(AliITSUClusterizer)
 
@@ -19,6 +21,8 @@ AliITSUClusterizer::AliITSUClusterizer(Int_t initNRow)
   ,fRecoParam(0)
   ,fInputDigits(0)
   ,fInputDigitsReadIndex(0)
+  ,fLayerID(0)
+  ,fLorAngCorrection(0)
   ,fOutputClusters(0)
   ,fDigitFreelist(0)
   ,fPartFreelist(0)
@@ -172,6 +176,7 @@ void AliITSUClusterizer::Transform(AliITSUClusterPix *cluster,AliITSUClusterizer
     nx = 1+Nint(dx/px);
     nz = 1+Nint(dz/pz);
   }
+  x -= fLorAngCorrection;  // LorentzAngle correction
   cluster->SetX(x);
   cluster->SetZ(z);
   cluster->SetY(0);
@@ -295,3 +300,11 @@ void AliITSUClusterizer::Clusterize()
   CloseRemainingParts(iNextRowBegin);
   return;
 }
+
+//______________________________________________________________________________
+void AliITSUClusterizer::PrepareLorentzAngleCorrection(Double_t bz)
+{
+  // calculate parameters for Lorentz Angle correction. Must be called 
+  // after setting segmentation and recoparams
+  fLorAngCorrection = 0.5*fRecoParam->GetTanLorentzAngle(fLayerID)*bz/kNominalBz*fSegm->Dy();
+}
index 6a07e24cff4eb386cda8123b5955256581de154c..0ab05966059c1e08e7ad76aedeac545629e2acae 100644 (file)
@@ -21,8 +21,11 @@ class AliITSUClusterizer : public TObject
   void Clusterize();\r
   void SetSegmentation(const AliITSUSegmentationPix *segm);\r
   void SetRecoParam(const AliITSURecoParam* param)     {fRecoParam = param;}\r
+  void SetLayerID(Int_t id)                            {fLayerID = id;}\r
   void SetVolID(Int_t id)                              {fVolID = id;}\r
   void SetNRow(Int_t nrow);\r
+  void PrepareLorentzAngleCorrection(Double_t bz);\r
+  //\r
   // interface methods\r
   void MakeRecPointBranch(TTree */*treeR*/)            {};\r
   void SetRecPointTreeAddress(TTree */*treeR*/)        {};\r
@@ -94,6 +97,8 @@ class AliITSUClusterizer : public TObject
   // Digit Input\r
   const TClonesArray *fInputDigits;         // supplied digits\r
   Int_t         fInputDigitsReadIndex;      // digits counter\r
+  Int_t         fLayerID;                   // current layer id\r
+  Double_t      fLorAngCorrection;          // Lorentz Angle correction for current layer\r
   // Cluster Output\r
   TClonesArray *fOutputClusters;            // external container to store clusters\r
   //\r
index aee4b10255ff49b328c21deb50deb90f91ed8f32..d4ce30e6efce593703e91cc141809ef305221f52 100644 (file)
 #include "AliITSUDigitPix.h"
 #include "AliITSUClusterizer.h"
 #include "AliITSUClusterPix.h"
+#include "AliMagF.h"
 
 ClassImp(AliITSUReconstructor)
 
 //___________________________________________________________________________
 AliITSUReconstructor::AliITSUReconstructor() 
 :  AliReconstructor()
-  ,fGM(0)
+  ,fGeom(0)
   ,fClusterFinders(0)
-  ,fRecPoints(0)
+  ,fClusters(0)
 {
   // Default constructor
 
@@ -54,49 +55,45 @@ AliITSUReconstructor::~AliITSUReconstructor()
 {
   // destructor
   //
-  if (!fGM) return; // was not initialized
+  if (!fGeom) return; // was not initialized
   //
   // same cluster finders and recpoint arrays might be attached to different layers
-  for (int i=fGM->GetNLayers();i--;) {
+  for (int i=fGeom->GetNLayers();i--;) {
     TObject* clFinder = fClusterFinders.At(i);
     if (clFinder) {
       while (fClusterFinders.Remove(clFinder)) {}
       delete clFinder;
     }
     //
-    TObject* arrRP = fRecPoints.At(i);
-    if (arrRP) {
-      while (fRecPoints.Remove(arrRP)) {}
-      delete arrRP;
-    }
+    delete[] fClusters;
   }
   //
-  delete fGM;
+  delete fGeom;
 } 
 
 //______________________________________________________________________
 void AliITSUReconstructor::Init() 
 {
   // Initalize this constructor 
-  if (fGM) AliFatal("was already done, something is wrong...");
+  AliInfo("Initializing");
+  if (fGeom) AliFatal("was already done, something is wrong...");
   //
-  fGM = new AliITSUGeomTGeo(kTRUE,kTRUE);
-  AliITSUClusterPix::SetGeom(fGM);
+  fGeom = new AliITSUGeomTGeo(kTRUE,kTRUE);
+  AliITSUClusterPix::SetGeom(fGeom);
   //  
   AliITSUClusterizer* clusPIX = 0;
-  TClonesArray* rpArrayPix = 0;
+  fClusters = new TClonesArray*[fGeom->GetNLayers()];
   //
-  for (int ilr=fGM->GetNLayers();ilr--;) {
-    int tpDet = fGM->GetLayerDetTypeID(ilr)/AliITSUGeomTGeo::kMaxSegmPerDetType;
+  for (int ilr=fGeom->GetNLayers();ilr--;) {
+    fClusters[ilr] = 0;
+    int tpDet = fGeom->GetLayerDetTypeID(ilr)/AliITSUGeomTGeo::kMaxSegmPerDetType;
     if (tpDet == AliITSUGeomTGeo::kDetTypePix) {
       if (!clusPIX)    clusPIX    = new AliITSUClusterizer();
-      if (!rpArrayPix) rpArrayPix = new TClonesArray(AliITSUClusterPix::Class());
-      //
       fClusterFinders.AddAtAndExpand(clusPIX, ilr);
-      fRecPoints.AddAtAndExpand(rpArrayPix, ilr);
+      fClusters[ilr] = new TClonesArray(AliITSUClusterPix::Class());
       //
       // to expand the buffers to max.size
-      clusPIX->SetSegmentation((AliITSUSegmentationPix*)fGM->GetSegmentation(ilr)); 
+      clusPIX->SetSegmentation((AliITSUSegmentationPix*)fGeom->GetSegmentation(ilr)); 
       continue;
     }
     else {
@@ -120,15 +117,13 @@ void AliITSUReconstructor::Reconstruct(TTree *digitsTree, TTree *clustersTree) c
   digitsTree->SetBranchAddress("ITSDigitsPix",&digArrPix);
   //
   // a new tree is created for each event: add each layer as separate branch
-  TBranch *lrBranch[fGM->GetNLayers()];
-  TClonesArray *rpClones[fGM->GetNLayers()];
+  TBranch *lrBranch[fGeom->GetNLayers()];
   //
-  for (int ilr=0;ilr<fGM->GetNLayers();ilr++) {
-    rpClones[ilr] = (TClonesArray*)fRecPoints.At(ilr);
+  for (int ilr=0;ilr<fGeom->GetNLayers();ilr++) {
     if (clustersTree) { // do we write clusters tree?
-      int tp = fGM->GetLayerDetTypeID(ilr)/AliITSUGeomTGeo::kMaxSegmPerDetType;
+      int tp = fGeom->GetLayerDetTypeID(ilr)/AliITSUGeomTGeo::kMaxSegmPerDetType;
       if (tp==AliITSUGeomTGeo::kDetTypePix) {
-       lrBranch[ilr] = clustersTree->Bronch(Form("ITSRecPoints%d",ilr),"TClonesArray",&rpClones[ilr]);
+       lrBranch[ilr] = clustersTree->Bronch(Form("ITSRecPoints%d",ilr),"TClonesArray",&fClusters[ilr]);
       }
       else {
        AliFatal(Form("Detector type %d is not defined",tp));
@@ -137,17 +132,23 @@ void AliITSUReconstructor::Reconstruct(TTree *digitsTree, TTree *clustersTree) c
   }
   //
   AliITSUClusterizer* clFinder = 0;
+  AliMagF* field = dynamic_cast<AliMagF*>(TGeoGlobalMagField::Instance()->GetField());
+  double bz = 0;
+  if (field == 0) AliError("Cannot get magnetic field from TGeoGlobalMagField");
+  else bz = field->SolenoidField();
   //
-  for (int ilr=0;ilr<fGM->GetNLayers();ilr++) {
+  for (int ilr=0;ilr<fGeom->GetNLayers();ilr++) {
     //
-    rpClones[ilr]->Clear();
+    fClusters[ilr]->Clear();
     clFinder = (AliITSUClusterizer*)fClusterFinders[ilr];
-    clFinder->SetSegmentation((AliITSUSegmentationPix*)fGM->GetSegmentation(ilr));
-    clFinder->SetClusters(rpClones[ilr]);
+    clFinder->SetSegmentation((AliITSUSegmentationPix*)fGeom->GetSegmentation(ilr));
+    clFinder->SetLayerID(ilr);
+    clFinder->SetClusters(fClusters[ilr]);
     clFinder->SetRecoParam(GetRecoParam()); // RS: Do we need to set it for every event?
+    clFinder->PrepareLorentzAngleCorrection(bz);
     //
-    int modF=fGM->GetFirstModIndex(ilr);
-    int modL=fGM->GetLastModIndex(ilr)+1;
+    int modF=fGeom->GetFirstModIndex(ilr);
+    int modL=fGeom->GetLastModIndex(ilr)+1;
     for (int imod=modF;imod<modL;imod++) {
       digitsTree->GetEntry(imod);   
       int ndig  = digArrPix->GetEntries();
@@ -158,8 +159,8 @@ void AliITSUReconstructor::Reconstruct(TTree *digitsTree, TTree *clustersTree) c
     }
     //
     AliITSUClusterPix::SetSortMode( AliITSUClusterPix::SortModeIdTrkYZ());
-    rpClones[ilr]->Sort();
-    AliDebug(1,Form(" -> Lr%d : %d Cluster",ilr,rpClones[ilr]->GetEntries()));
+    fClusters[ilr]->Sort();
+    AliDebug(1,Form(" -> Lr%d : %d Cluster",ilr,fClusters[ilr]->GetEntries()));
     if (clustersTree) lrBranch[ilr]->Fill();
   }
   if (clustersTree) clustersTree->SetEntries();
@@ -240,3 +241,17 @@ AliTracker* AliITSUReconstructor::CreateTrackleter() const
 
 }
 
+//_____________________________________________________________________________
+Int_t AliITSUReconstructor::LoadClusters(TTree* treeRP) 
+{
+  // read clusters from the tree, if it is provided
+  if (!treeRP) return 0;
+  for (int ilr=fGeom->GetNLayers();ilr--;) {
+    if (!fClusters[ilr]) AliFatal(Form("Clusters array for layer %d is not defined",ilr)); 
+    TBranch* br = treeRP->GetBranch(Form("ITSRecPoints%d",ilr));
+    if (!br) AliFatal(Form("Provided cluster tree does not contain branch for layer %d",ilr));
+    br->SetAddress(&fClusters[ilr]);
+  }
+  treeRP->GetEntry(0); // we are still in 1 ev/tree mode...
+  return 1;
+}
index 0eb3ad589a27b97c8a4547665c33c13b12eb271b..1dafcef2327426ac9247f07ccaa38821b53b9f6b 100644 (file)
@@ -21,16 +21,22 @@ class AliITSUReconstructor: public AliReconstructor {
 public:
   AliITSUReconstructor();
   virtual ~AliITSUReconstructor();
-  virtual void         Init();
-  virtual void         Reconstruct(TTree *digitsTree, TTree *clustersTree) const;
-  virtual void         Reconstruct(AliRawReader*, TTree*) const {};
-
-  virtual AliTracker*    CreateTracker() const;
-  virtual AliVertexer*   CreateVertexer() const;
+  virtual void          Init();
+  virtual void          Reconstruct(TTree *digitsTree, TTree *clustersTree) const;
+  virtual void          Reconstruct(AliRawReader*, TTree*) const {};
+  //
+  virtual AliTracker*    CreateTracker()    const;
+  virtual AliVertexer*   CreateVertexer()   const;
   virtual AliTrackleter* CreateMultFinder() const;
   virtual AliTracker*    CreateTrackleter() const;
-
+  //
   virtual const char*    GetDetectorName() const {return "ITS";}
+  //
+  TClonesArray*          GetClusters(Int_t lrID)            const {return fClusters ? fClusters[lrID] : 0;}
+  AliITSUGeomTGeo*       GetGeom()                          const {return (AliITSUGeomTGeo*)fGeom;}
+  //
+  Int_t                  LoadClusters(TTree* treeRP);
+  //
 
   static const AliITSURecoParam* GetRecoParam() { 
     return dynamic_cast<const AliITSURecoParam*>(AliReconstructor::GetRecoParam(0)); }
@@ -38,10 +44,9 @@ public:
 private:
   AliITSUReconstructor(const AliITSUReconstructor &); //Not implemented
   AliITSUReconstructor& operator=(const AliITSUReconstructor &); //Not implemented
-
-  AliITSUGeomTGeo* fGM;          // geometry wrapper
+  AliITSUGeomTGeo* fGeom;          // geometry wrapper
   TObjArray        fClusterFinders; // array of clusterfinders per layer
-  TObjArray        fRecPoints;      // container for recpoints TClonesArrays
+  TClonesArray**   fClusters;      // container for recpoints TClonesArrays
   //
   ClassDef(AliITSUReconstructor, 0)   // class for the ITSU reconstruction
 };