]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - MUON/MUONAlignment.C
- AliMUONRecoParam.cxx:
[u/mrichter/AliRoot.git] / MUON / MUONAlignment.C
index 40af52d6ca8e52f5325b5eb3a2a829071d7dd451..a418eaf1de4f4334702b8414066cae5d1c096360 100644 (file)
 
 /* $Id$ */
 
-// ---
-// Macro for MUON alignment using physics tracks. The macro uses AliMUONAlignment
-// class to calculate the alignment parameters.
-// An array for the alignment parameters is created and can be filled with
-// initial values that will be used as starting values by the alignment
-// algorithm.
-// Ny default the macro run over galice.root in the working directory. If a file list
-// of galice.root is provided as third argument the macro will run over all files.
-// The macro loop over the files, events and tracks. For each track
-// AliMUONAlignment::ProcessTrack(AliMUONTrack * track) and then
-// AliMUONAlignment::LocalFit(Int_t iTrack, Double_t *lTrackParam, Int_t
-// lSingleFit) are called. After all tracks have been procesed and fitted
-// AliMUONAlignment::GlobalFit(Double_t *parameters,Double_t *errors,Double_t *pulls)
-// is called. The array parameters contains the obatained misalignement parameters.
-// A realigned geometry is generated in a local CDB.
-//
-// Authors: B. Becker and J. Castillo
-// ---
+/// \ingroup macros
+/// \file MUONAlignment.C
+/// \brief Macro for MUON alignment using physics tracks. 
+///
+/// The macro uses the AliMUONAlignment class to calculate the alignment parameters.
+/// An array for the alignment parameters is created and can be filled with
+/// initial values that will be used as starting values by the alignment
+/// algorithm.
+///
+/// By default the macro run over galice.root in the working directory. If a file list
+/// of galice.root is provided as third argument the macro will run over all files.
+/// The macro loop over the files, events and tracks. For each track
+/// AliMUONAlignment::ProcessTrack(AliMUONTrack * track) and then
+/// AliMUONAlignment::LocalFit(Int_t iTrack, Double_t *lTrackParam, Int_t
+/// lSingleFit) are called. After all tracks have been procesed and fitted
+/// AliMUONAlignment::GlobalFit(Double_t *parameters,Double_t *errors,Double_t *pulls)
+/// is called. The array parameters contains the obatained misalignement parameters.
+/// A realigned geometry is generated in a local CDB.
+///
+/// \author: B. Becker and J. Castillo
 
 #if !defined(__CINT__) || defined(__MAKECINT__)
 
 #include "AliMUONTrackExtrap.h"
 #include "AliMUONTrackParam.h"
 #include "AliMUONGeometryTransformer.h"
-#include "AliMUONDataInterface.h"
-#include "AliMUONVTrackStore.h"
+#include "AliMUONESDInterface.h"
 
+#include "AliESDEvent.h"
+#include "AliESDMuonTrack.h"
 #include "AliMagFMaps.h"
 #include "AliTracker.h"
 #include "AliCDBManager.h"
@@ -56,6 +59,7 @@
 #include <TH1.h>
 #include <TGraphErrors.h>
 #include <TFile.h>
+#include <TTree.h>
 #include <TClonesArray.h>
 #include <Riostream.h>
 
@@ -63,7 +67,7 @@
 
 #endif
 
-void MUONAlignment(Int_t nEvents = 100000, char* geoFilename = "geometry.root", TString fileName = "galice.root", TString fileList = "")
+void MUONAlignment(Int_t nEvents = 100000, char* geoFilename = "geometry.root", TString esdFileName = "AliESDs.root", TString fileList = "")
 {
  
   // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
@@ -143,9 +147,6 @@ void MUONAlignment(Int_t nEvents = 100000, char* geoFilename = "geometry.root",
 
 
   char cFileName[100];  
-  AliMUONDataInterface *amdi =0x0;
-  AliMUONVTrackStore *trackStore =0x0;
-  AliMUONTrack* track =0x0;  
     
   Int_t lMaxFile = 1000;
   Int_t iFile = 0;
@@ -161,16 +162,28 @@ void MUONAlignment(Int_t nEvents = 100000, char* geoFilename = "geometry.root",
       if (sFileList.eof()) bKeepLoop = kFALSE;
     }
     else {
-      sprintf(cFileName,fileName.Data());
+      sprintf(cFileName,esdFileName.Data());
       bKeepLoop = kFALSE;
     }
-    if (!strstr(cFileName,"galice.root")) continue;      
+    if (!strstr(cFileName,"AliESDs.root")) continue;      
     cout << "Using file: " << cFileName << endl;
+    
+    // load ESD event
+    TFile* esdFile = TFile::Open(cFileName); // open the file
+    if (!esdFile || !esdFile->IsOpen()) {
+      cout << "opening ESD file " << cFileName << "failed" << endl;
+      continue;
+    }
+    TTree* esdTree = (TTree*) esdFile->Get("esdTree"); // get the tree
+    if (!esdTree) {
+      cout << "no ESD tree found" << endl;
+      esdFile->Close();
+      continue;
+    }
+    AliESDEvent* esdEvent = new AliESDEvent(); // link ESD event to the tree
+    esdEvent->ReadFromTree(esdTree);
 
-    if (!amdi) amdi = new AliMUONDataInterface(cFileName);
-    else amdi->Open(cFileName);
-
-    Int_t nevents = amdi->NumberOfEvents();
+    Int_t nevents = esdTree->GetEntries();
     cout << "... with " << nevents << endl;
     for(Int_t event = 0; event < nevents; event++) {
       if (iEvent >= nEvents){
@@ -179,25 +192,30 @@ void MUONAlignment(Int_t nEvents = 100000, char* geoFilename = "geometry.root",
       }
       iEvent++;
 
-      trackStore = amdi->TrackStore(event);
+      if (esdTree->GetEvent(event) <= 0) {
+       cout << "fails to read ESD object for event " << event << endl;
+       continue;
+      }
 
-      Int_t ntracks = trackStore->GetSize();
-      if (!event%100)
-       cout << " there are " << ntracks << " tracks in event " << event << endl;
-      TIter next(trackStore->CreateIterator());
-      while ( ( track = static_cast<AliMUONTrack*>(next()) ) ) {
-       AliMUONTrackParam trackParam(*((AliMUONTrackParam*)(track->GetTrackParamAtHit()->First())));
-       AliMUONTrackExtrap::ExtrapToVertex(&trackParam,0.,0.,0.);
-       Double_t invBenMom = trackParam.GetInverseBendingMomentum();
+      Int_t nTracks = Int_t(esdEvent->GetNumberOfMuonTracks());
+      if (!event%100) cout << " there are " << nTracks << " tracks in event " << event << endl;
+      for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
+       AliESDMuonTrack* esdTrack = esdEvent->GetMuonTrack(iTrack);
+       if (!esdTrack->ClustersStored()) continue;
+       Double_t invBenMom = esdTrack->GetInverseBendingMomentum();
        fInvBenMom->Fill(invBenMom);
        fBenMom->Fill(1./invBenMom);
        if (TMath::Abs(invBenMom)<=1.04) {
-         alig->ProcessTrack(track);
+         AliMUONTrack track;
+         AliMUONESDInterface::ESDToMUON(*esdTrack, track);
+         alig->ProcessTrack(&track);
          alig->LocalFit(iTrackOk++,trackParams,0);
        }
        iTrackTot++;
       }
     }
+    delete esdEvent;
+    esdFile->Close();
     cout << "Processed " << iTrackTot << " Tracks so far." << endl;
   }
   alig->GlobalFit(parameters,errors,pulls);
@@ -270,6 +288,9 @@ void MUONAlignment(Int_t nEvents = 100000, char* geoFilename = "geometry.root",
   
   // Generate realigned data in local cdb
   const TClonesArray* array = newTransform->GetMisAlignmentData();
+
+  // 100 mum residual resolution for chamber misalignments?
+  alig->SetAlignmentResolution(array,-1,0.01,0.01,0.004,0.003);
    
   // CDB manager
   AliCDBManager* cdbManager = AliCDBManager::Instance();
@@ -281,4 +302,5 @@ void MUONAlignment(Int_t nEvents = 100000, char* geoFilename = "geometry.root",
   AliCDBId id("MUON/Align/Data", 0, 0); 
   cdbManager->Put(const_cast<TClonesArray*>(array), id, cdbData);
 
-} 
+}
+