]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - ITS/AliITSReconstructor.cxx
Parameters of AliITSVertexer3D added in AliITSRecoParam (F.Prino)
[u/mrichter/AliRoot.git] / ITS / AliITSReconstructor.cxx
index 932968794a7b4bb47e08aa0b271ebeda5c9d9e6e..789d9411906c4a125c0aca1c1d6633c94aff6c66 100644 (file)
 #include "Riostream.h"
 #include "AliITSReconstructor.h"
 #include "AliRun.h"
-#include "AliRunLoader.h"
 #include "AliRawReader.h"
 #include "AliITSDetTypeRec.h"
+#include "AliITSgeom.h"
 #include "AliITSLoader.h"
 #include "AliITStrackerMI.h"
+#include "AliITStrackerV2.h"
 #include "AliITStrackerSA.h"
 #include "AliITSVertexerIons.h"
 #include "AliITSVertexerFast.h"
-#include "AliITSVertexerPPZ.h"
+#include "AliITSVertexer3D.h"
 #include "AliITSVertexerZ.h"
-#include "AliESD.h"
+#include "AliITSVertexerCosmics.h"
+#include "AliESDEvent.h"
 #include "AliITSpidESD.h"
 #include "AliITSpidESD1.h"
 #include "AliITSpidESD2.h"
 #include "AliITSInitGeometry.h"
 
+
 ClassImp(AliITSReconstructor)
 
 //___________________________________________________________________________
 AliITSReconstructor::AliITSReconstructor() : AliReconstructor(),
-fItsPID(0)
+fItsPID(0),
+fDetTypeRec(0)
 {
   // Default constructor
 }
@@ -52,10 +56,13 @@ fItsPID(0)
 AliITSReconstructor::~AliITSReconstructor(){
 // destructor
   delete fItsPID;
+  if(fDetTypeRec) delete fDetTypeRec;
 } 
 //______________________________________________________________________
 AliITSReconstructor::AliITSReconstructor(const AliITSReconstructor &ob) :AliReconstructor(ob),
-fItsPID(ob.fItsPID) 
+fItsPID(ob.fItsPID),
+fDetTypeRec(ob.fDetTypeRec)
+
 {
   // Copy constructor
 }
@@ -67,14 +74,13 @@ AliITSReconstructor& AliITSReconstructor::operator=(const AliITSReconstructor&
   new(this) AliITSReconstructor(ob);
   return *this;
 }
+
 //______________________________________________________________________
-void AliITSReconstructor::Init(AliRunLoader *runLoader) const{
+void AliITSReconstructor::Init({
     // Initalize this constructor bet getting/creating the objects
     // nesseary for a proper ITS reconstruction.
     // Inputs:
-    //    AliRunLoader *runLoader   Pointer to the run loader to allow
-    //                              the getting of files/folders data
-    //                              needed to do reconstruction
+    //   none.
     // Output:
     //   none.
     // Return:
@@ -83,177 +89,133 @@ void AliITSReconstructor::Init(AliRunLoader *runLoader) const{
     AliITSInitGeometry initgeom;
     AliITSgeom *geom = initgeom.CreateAliITSgeom();
     AliInfo(Form("Geometry name: %s",(initgeom.GetGeometryName()).Data()));
-    AliITSLoader* loader = static_cast<AliITSLoader*>
-       (runLoader->GetLoader("ITSLoader"));
-    if (!loader) {
-       Error("Init", "ITS loader not found");
-       return;
-    }
-    loader->SetITSgeom(geom);
+
+    fDetTypeRec = new AliITSDetTypeRec();
+    fDetTypeRec->SetITSgeom(geom);
+    fDetTypeRec->SetDefaults();
+
     return;
 }
+
 //_____________________________________________________________________________
-void AliITSReconstructor::Reconstruct(AliRunLoader* runLoader) const
+void AliITSReconstructor::Reconstruct(TTree *digitsTree, TTree *clustersTree) const
 {
 // reconstruct clusters
 
+  Int_t cluFindOpt = GetRecoParam()->GetClusterFinder();
+  Bool_t useV2=kTRUE;   // Default: V2 cluster finder
+  if(cluFindOpt==1) useV2=kFALSE;
 
-  AliITSLoader* loader = static_cast<AliITSLoader*>(runLoader->GetLoader("ITSLoader"));
-  if (!loader) {
-    Error("Reconstruct", "ITS loader not found");
-    return;
-  }
-  AliITSDetTypeRec* rec = new AliITSDetTypeRec(loader);
-  rec->SetDefaults();
-
-  loader->LoadRecPoints("recreate");
-  loader->LoadDigits("read");
-  runLoader->LoadKinematics();
-  TString option = GetOption();
-  Bool_t clusfinder=kTRUE;   // Default: V2 cluster finder
-  if(option.Contains("OrigCF"))clusfinder=kFALSE;
-
-  Int_t nEvents = runLoader->GetNumberOfEvents();
-
-  for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
-    runLoader->GetEvent(iEvent);
-    if(loader->TreeR()==0x0) loader->MakeTree("R");
-    rec->MakeBranch("R");
-    rec->SetTreeAddress();
-    rec->DigitsToRecPoints(iEvent,0,"All",clusfinder);    
-  }
-
-  loader->UnloadRecPoints();
-  loader->UnloadDigits();
-  runLoader->UnloadKinematics();
+  fDetTypeRec->SetTreeAddressD(digitsTree);
+  fDetTypeRec->MakeBranch(clustersTree,"R");
+  fDetTypeRec->SetTreeAddressR(clustersTree);
+  fDetTypeRec->DigitsToRecPoints(digitsTree,clustersTree,0,"All",useV2);    
 }
 
 //_________________________________________________________________
-void AliITSReconstructor::Reconstruct(AliRunLoader* runLoader, 
-                                      AliRawReader* rawReader) const
+void AliITSReconstructor::Reconstruct(AliRawReader* rawReader, TTree *clustersTree) const
 {
-// reconstruct clusters
-
+  // reconstruct clusters from raw data
  
-  AliITSLoader* loader = static_cast<AliITSLoader*>(runLoader->GetLoader("ITSLoader"));
-  if (!loader) {
-    Error("Reconstruct", "ITS loader not found");
-    return;
-  }
-
-  AliITSDetTypeRec* rec = new AliITSDetTypeRec(loader);
-  rec->SetDefaults();
-  rec->SetDefaultClusterFindersV2(kTRUE);
-
-  loader->LoadRecPoints("recreate");
-
-  Int_t iEvent = 0;
-
-  while(rawReader->NextEvent()) {
-    runLoader->GetEvent(iEvent++);
-    if(loader->TreeR()==0x0) loader->MakeTree("R");
-    rec->DigitsToRecPoints(rawReader);
-  }
-
-  loader->UnloadRecPoints();
+  fDetTypeRec->SetDefaultClusterFindersV2(kTRUE);
+  fDetTypeRec->DigitsToRecPoints(rawReader,clustersTree);
 }
 
 //_____________________________________________________________________________
-AliTracker* AliITSReconstructor::CreateTracker(AliRunLoader* runLoader)const
+AliTracker* AliITSReconstructor::CreateTracker(const
 {
 // create a ITS tracker
 
-  
-  TString selectedTracker = GetOption();
+  Int_t trackerOpt = GetRecoParam()->GetTracker();
   AliTracker* tracker;    
-  if (selectedTracker.Contains("MI")) {
+  if (trackerOpt==1) {
     tracker = new AliITStrackerMI(0);
+    AliITStrackerMI *mit=(AliITStrackerMI*)tracker;
+    mit->SetDetTypeRec(fDetTypeRec);
+  }  
+  else if (trackerOpt==2) {
+    tracker = new AliITStrackerV2(0);
   }
   else {
     tracker =  new AliITStrackerSA(0);  // inherits from AliITStrackerMI
+    AliITStrackerSA *sat=(AliITStrackerSA*)tracker;
+    sat->SetDetTypeRec(fDetTypeRec);
+    if(GetRecoParam()->GetTrackerSAOnly()) sat->SetSAFlag(kTRUE);
+    if(sat->GetSAFlag())AliDebug(1,"Tracking Performed in ITS only\n");
+    sat->SetOuterStartLayer(GetRecoParam()->GetOuterStartLayerSA());
   }
 
-  TString selectedPIDmethod = GetOption();
-  AliITSLoader *loader = (AliITSLoader*)runLoader->GetLoader("ITSLoader");
-  if (!loader) {
-    Error("CreateTracker", "ITS loader not found");
-  }
-  if(selectedPIDmethod.Contains("LandauFitPID")){
-    loader->AdoptITSpid(new AliITSpidESD2((AliITStrackerMI*)tracker,loader));
+  Int_t pidOpt = GetRecoParam()->GetPID();
+
+  AliITSReconstructor* nc = const_cast<AliITSReconstructor*>(this);
+  if(pidOpt==1){
+    Info("FillESD","ITS LandauFitPID option has been selected\n");
+    nc->fItsPID = new AliITSpidESD2((AliITStrackerMI*)tracker);
   }
   else{
-    Double_t parITS[] = {34., 0.15, 10.};
-    loader->AdoptITSpid(new AliITSpidESD1(parITS));
+    Info("FillESD","ITS default PID\n");
+    Double_t parITS[] = {79.,0.13, 5.}; //IB: this is  "pp tuning"
+    nc->fItsPID = new AliITSpidESD1(parITS);
   }
   return tracker;
   
 }
 
 //_____________________________________________________________________________
-AliVertexer* AliITSReconstructor::CreateVertexer(AliRunLoader* /*runLoader*/) const
+AliVertexer* AliITSReconstructor::CreateVertexer() const
 {
 // create a ITS vertexer
 
-  TString selectedVertexer = GetOption();
-  if(selectedVertexer.Contains("ions") || selectedVertexer.Contains("IONS")){
+  Int_t vtxOpt = GetRecoParam()->GetVertexer();
+  if(vtxOpt==3){
     Info("CreateVertexer","a AliITSVertexerIons object has been selected\n");
-    return new AliITSVertexerIons("null");
+    return new AliITSVertexerIons();
   }
-  if(selectedVertexer.Contains("smear") || selectedVertexer.Contains("SMEAR")){
+  if(vtxOpt==4){
     Double_t smear[3]={0.005,0.005,0.01};
     Info("CreateVertexer","a AliITSVertexerFast object has been selected\n"); 
     return new AliITSVertexerFast(smear);
   }
-  if(selectedVertexer.Contains("ppz") || selectedVertexer.Contains("PPZ")){
-    Info("CreateVertexer","a AliITSVertexerPPZ object has been selected\n");
-    return new AliITSVertexerPPZ("null");
-  }
-  // by default an AliITSVertexerZ object is instatiated
-  Info("CreateVertexer","a AliITSVertexerZ object has been selected\n");
-  return new AliITSVertexerZ("null");
+  if(vtxOpt==1){
+    Info("CreateVertexer","a AliITSVertexerZ object has been selected\n");
+    return new AliITSVertexerZ();
+  }
+  if(vtxOpt==2){
+    Info("CreateVertexer","a AliITSVertexerCosmics object has been selected\n");
+    return new AliITSVertexerCosmics();
+  }
+  // by default an AliITSVertexer3D object is instatiated
+  Info("CreateVertexer","a AliITSVertexer3D object has been selected\n");
+  AliITSVertexer3D*  vtxr = new AliITSVertexer3D();
+  Float_t dzw=GetRecoParam()->GetVertexer3DWideFiducialRegionZ();
+  Float_t drw=GetRecoParam()->GetVertexer3DWideFiducialRegionR();
+  vtxr->SetWideFiducialRegion(dzw,drw);
+  Float_t dzn=GetRecoParam()->GetVertexer3DNarrowFiducialRegionZ();
+  Float_t drn=GetRecoParam()->GetVertexer3DNarrowFiducialRegionR();
+  vtxr->SetNarrowFiducialRegion(dzn,drn);
+  Float_t dphil=GetRecoParam()->GetVertexer3DLooseDeltaPhiCut();
+  Float_t dphit=GetRecoParam()->GetVertexer3DTightDeltaPhiCut();
+  vtxr->SetDeltaPhiCuts(dphil,dphit);
+  Float_t dcacut=GetRecoParam()->GetVertexer3DDCACut();
+  vtxr->SetDCACut(dcacut);
+  return vtxr;
 }
 
 //_____________________________________________________________________________
-void AliITSReconstructor::FillESD(AliRunLoader* runLoader
-                                 AliESD* esd) const
+void AliITSReconstructor::FillESD(TTree * /*digitsTree*/, TTree *clustersTree
+                                 AliESDEvent* esd) const
 {
 // make PID, find V0s and cascade
-  AliITSLoader *loader = (AliITSLoader*)runLoader->GetLoader("ITSLoader");
-  AliITSpidESD *pidESD = 0;
-  TString selectedPIDmethod = GetOption();
-  if(selectedPIDmethod.Contains("LandauFitPID")){
-    Info("FillESD","ITS LandauFitPID option has been selected\n");
-    pidESD=loader->GetITSpid();
-  }
-  else{
-    Info("FillESD","ITS default PID\n");
-    pidESD=loader->GetITSpid();
-  }
-  if(pidESD!=0){
-    pidESD->MakePID(esd);
+  if(fItsPID!=0) {
+    Int_t pidOpt = GetRecoParam()->GetPID();
+    if(pidOpt==1){
+      fItsPID->MakePID(clustersTree,esd);
+    }else{
+      fItsPID->MakePID(esd);
+    }
   }
   else {
     Error("FillESD","!! cannot do the PID !!\n");
   }
 }
-
-
-//_____________________________________________________________________________
-AliITSgeom* AliITSReconstructor::GetITSgeom(AliRunLoader* runLoader) const
-{
-// get the ITS geometry
-
-  if (!runLoader->GetAliRun()) runLoader->LoadgAlice();
-  if (!runLoader->GetAliRun()) {
-    Error("GetITSgeom", "couldn't get AliRun object");
-    return NULL;
-  }
-  AliITSLoader *loader = (AliITSLoader*)runLoader->GetLoader("ITSLoader");
-  AliITSgeom* geom = (AliITSgeom*)loader->GetITSgeom();
-  if(!geom){
-    Error("GetITSgeom","no ITS geometry available");
-    return NULL;
-  }
-  
-  return geom;
-}