]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONAlignment.C
- Set alignment resolution in alignment objects
[u/mrichter/AliRoot.git] / MUON / MUONAlignment.C
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 /// \ingroup macros
19 /// \file MUONAlignment.C
20 /// \brief Macro for MUON alignment using physics tracks. 
21 ///
22 /// The macro uses the AliMUONAlignment class to calculate the alignment parameters.
23 /// An array for the alignment parameters is created and can be filled with
24 /// initial values that will be used as starting values by the alignment
25 /// algorithm.
26 ///
27 /// By default the macro run over galice.root in the working directory. If a file list
28 /// of galice.root is provided as third argument the macro will run over all files.
29 /// The macro loop over the files, events and tracks. For each track
30 /// AliMUONAlignment::ProcessTrack(AliMUONTrack * track) and then
31 /// AliMUONAlignment::LocalFit(Int_t iTrack, Double_t *lTrackParam, Int_t
32 /// lSingleFit) are called. After all tracks have been procesed and fitted
33 /// AliMUONAlignment::GlobalFit(Double_t *parameters,Double_t *errors,Double_t *pulls)
34 /// is called. The array parameters contains the obatained misalignement parameters.
35 /// A realigned geometry is generated in a local CDB.
36 ///
37 /// \author: B. Becker and J. Castillo
38
39 #if !defined(__CINT__) || defined(__MAKECINT__)
40
41 #include "AliMUONAlignment.h"
42 #include "AliMUONTrack.h"
43 #include "AliMUONTrackExtrap.h"
44 #include "AliMUONTrackParam.h"
45 #include "AliMUONGeometryTransformer.h"
46 #include "AliMUONESDInterface.h"
47 #include "AliMUONRecoParam.h"
48
49 #include "AliESDEvent.h"
50 #include "AliESDMuonTrack.h"
51 #include "AliMagFMaps.h"
52 #include "AliTracker.h"
53 #include "AliCDBManager.h"
54 #include "AliCDBMetaData.h"
55 #include "AliCDBId.h"
56 #include "AliGeomManager.h"
57
58 #include <TString.h>
59 #include <TError.h>
60 #include <TH1.h>
61 #include <TGraphErrors.h>
62 #include <TFile.h>
63 #include <TTree.h>
64 #include <TClonesArray.h>
65 #include <Riostream.h>
66
67 #include <fstream>
68
69 #endif
70
71 void MUONAlignment(Int_t nEvents = 100000, char* geoFilename = "geometry.root", TString esdFileName = "AliESDs.root", TString fileList = "")
72 {
73  
74   // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
75   if ( ! AliGeomManager::GetGeometry() ) {
76     AliGeomManager::LoadGeometry(geoFilename);
77     if (! AliGeomManager::GetGeometry() ) {
78       Error("MUONAlignment", "getting geometry from file %s failed", geoFilename);
79       return;
80     }
81   }
82   
83   // set  mag field 
84   // waiting for mag field in CDB 
85   printf("Loading field map...\n");
86   AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG);
87   AliTracker::SetFieldMap(field, kFALSE);
88   // set the magnetic field for track extrapolations
89   AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap());
90
91   Double_t parameters[3*156];
92   Double_t errors[3*156];
93   Double_t pulls[3*156];
94   for(Int_t k=0;k<3*156;k++) {
95     parameters[k]=0.;
96     errors[k]=0.;
97     pulls[k]=0.;
98   }
99
100   Double_t trackParams[8] = {0.,0.,0.,0.,0.,0.,0.,0.};
101
102   // Set initial values here, good guess may help convergence
103   // St 1 
104   //  Int_t iPar = 0;
105   //  parameters[iPar++] =  0.010300 ;  parameters[iPar++] =  0.010600 ;  parameters[iPar++] =  0.000396 ;  
106
107   bool bLoop = kFALSE;
108   ifstream sFileList;
109   if (fileList.Contains(".list")) {
110     cout << "Reading file list: " << fileList.Data() << endl;
111     bLoop = kTRUE;
112
113     TString fullListName(".");
114     fullListName +="/";
115     fullListName +=fileList;
116     sFileList.open(fileList.Data());
117   }
118   
119   TH1F *fInvBenMom = new TH1F("fInvBenMom","fInvBenMom",200,-0.1,0.1); 
120   TH1F *fBenMom = new TH1F("fBenMom","fBenMom",200,-40,40); 
121
122   AliMUONAlignment* alig = new AliMUONAlignment();
123   alig->InitGlobalParameters(parameters);
124
125   AliMUONGeometryTransformer *transform = new AliMUONGeometryTransformer();
126   transform->LoadGeometryData();
127   alig->SetGeometryTransformer(transform);
128
129   // Do alignment with magnetic field off
130   alig->SetBFieldOn(kFALSE);
131   
132   // Set tracking station to use
133   Bool_t bStOnOff[5] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
134
135   // Fix parameters or add constraints here
136   for (Int_t iSt=0; iSt<5; iSt++)
137     if (!bStOnOff[iSt]) alig->FixStation(iSt+1);
138
139   // Left and right sides of the detector are independent, one can choose to align 
140   // only one side
141   Bool_t bSpecLROnOff[2] = {kTRUE,kTRUE};
142   alig->FixHalfSpectrometer(bStOnOff,bSpecLROnOff);
143
144   // Set predifined global constrains: X, Y, P, XvsZ, YvsZ, PvsZ, XvsY, YvsY, PvsY
145   Bool_t bVarXYT[9] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
146   Bool_t bDetTLBR[4] = {kFALSE,kTRUE,kFALSE,kTRUE};
147   alig->AddConstraints(bStOnOff,bVarXYT,bDetTLBR,bSpecLROnOff);
148
149
150   char cFileName[100];  
151     
152   Int_t lMaxFile = 1000;
153   Int_t iFile = 0;
154   Int_t iEvent = 0;
155   bool bKeepLoop = kTRUE;
156   Int_t iTrackTot=0;
157   Int_t iTrackOk=0;
158
159   while(bKeepLoop && iFile<lMaxFile){
160     iFile++;
161     if (bLoop) {
162       sFileList.getline(cFileName,100);
163       if (sFileList.eof()) bKeepLoop = kFALSE;
164     }
165     else {
166       sprintf(cFileName,esdFileName.Data());
167       bKeepLoop = kFALSE;
168     }
169     if (!strstr(cFileName,"AliESDs.root")) continue;      
170     cout << "Using file: " << cFileName << endl;
171     
172     // load ESD event
173     TFile* esdFile = TFile::Open(cFileName); // open the file
174     if (!esdFile || !esdFile->IsOpen()) {
175       cout << "opening ESD file " << cFileName << "failed" << endl;
176       continue;
177     }
178     TTree* esdTree = (TTree*) esdFile->Get("esdTree"); // get the tree
179     if (!esdTree) {
180       cout << "no ESD tree found" << endl;
181       esdFile->Close();
182       continue;
183     }
184     AliESDEvent* esdEvent = new AliESDEvent(); // link ESD event to the tree
185     esdEvent->ReadFromTree(esdTree);
186
187     Int_t nevents = esdTree->GetEntries();
188     cout << "... with " << nevents << endl;
189     for(Int_t event = 0; event < nevents; event++) {
190       if (iEvent >= nEvents){
191         bKeepLoop = kFALSE;
192         break;
193       }
194       iEvent++;
195
196       if (esdTree->GetEvent(event) <= 0) {
197         cout << "fails to read ESD object for event " << event << endl;
198         continue;
199       }
200
201       Int_t nTracks = Int_t(esdEvent->GetNumberOfMuonTracks());
202       if (!event%100) cout << " there are " << nTracks << " tracks in event " << event << endl;
203       for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
204         AliESDMuonTrack* esdTrack = esdEvent->GetMuonTrack(iTrack);
205         if (!esdTrack->ClustersStored()) continue;
206         Double_t invBenMom = esdTrack->GetInverseBendingMomentum();
207         fInvBenMom->Fill(invBenMom);
208         fBenMom->Fill(1./invBenMom);
209         if (TMath::Abs(invBenMom)<=1.04) {
210           AliMUONTrack track;
211           AliMUONESDInterface::ESDToMUON(AliMUONRecoParam::GetLowFluxParam(),*esdTrack, track);
212           alig->ProcessTrack(&track);
213           alig->LocalFit(iTrackOk++,trackParams,0);
214         }
215         iTrackTot++;
216       }
217     }
218     delete esdEvent;
219     esdFile->Close();
220     cout << "Processed " << iTrackTot << " Tracks so far." << endl;
221   }
222   alig->GlobalFit(parameters,errors,pulls);
223
224   cout << "Done with GlobalFit " << endl;
225
226   // Store results
227   Double_t DEid[156] = {0};
228   Double_t MSDEx[156] = {0};
229   Double_t MSDEy[156] = {0};
230   Double_t MSDExt[156] = {0};
231   Double_t MSDEyt[156] = {0};
232   Double_t DEidErr[156] = {0};
233   Double_t MSDExErr[156] = {0};
234   Double_t MSDEyErr[156] = {0};
235   Double_t MSDExtErr[156] = {0};
236   Double_t MSDEytErr[156] = {0};
237   Int_t lNDetElem = 4*2+4*2+18*2+26*2+26*2;
238   Int_t lNDetElemCh[10] = {4,4,4,4,18,18,26,26,26,26};
239   // Int_t lSNDetElemCh[10] = {4,8,12,16,34,52,78,104,130,156};
240   Int_t idOffset = 0; // 400
241   Int_t lSDetElemCh = 0;
242   for(Int_t iDE=0; iDE<lNDetElem; iDE++){
243     DEidErr[iDE] = 0.;
244     DEid[iDE] = idOffset+100;
245     DEid[iDE] += iDE; 
246     lSDetElemCh = 0;
247     for(Int_t iCh=0; iCh<9; iCh++){
248       lSDetElemCh += lNDetElemCh[iCh];
249       if (iDE>=lSDetElemCh) {
250         DEid[iDE] += 100;
251         DEid[iDE] -= lNDetElemCh[iCh];
252       }
253     }
254     MSDEx[iDE]=parameters[3*iDE+0];
255     MSDEy[iDE]=parameters[3*iDE+1];
256     MSDExt[iDE]=parameters[3*iDE+2];
257     MSDEyt[iDE]=parameters[3*iDE+2];
258     MSDExErr[iDE]=(Double_t)alig->GetParError(3*iDE+0);
259     MSDEyErr[iDE]=(Double_t)alig->GetParError(3*iDE+1);
260     MSDExtErr[iDE]=(Double_t)alig->GetParError(3*iDE+2);
261     MSDEytErr[iDE]=(Double_t)alig->GetParError(3*iDE+2);
262   }
263
264   cout << "Let's create graphs ...  " << endl;
265
266   TGraphErrors *gMSDEx = new TGraphErrors(lNDetElem,DEid,MSDEx,DEidErr,MSDExErr); 
267   TGraphErrors *gMSDEy = new TGraphErrors(lNDetElem,DEid,MSDEy,DEidErr,MSDEyErr); 
268   TGraphErrors *gMSDExt = new TGraphErrors(lNDetElem,DEid,MSDExt,DEidErr,MSDExtErr); 
269   TGraphErrors *gMSDEyt = new TGraphErrors(lNDetElem,DEid,MSDEyt,DEidErr,MSDEytErr); 
270
271   cout << "... graphs created, open file ...  " << endl;
272
273   TFile *hFile = new TFile("measShifts.root","RECREATE");
274
275   cout << "... file opened ...  " << endl;
276
277   gMSDEx->Write("gMSDEx");
278   gMSDEy->Write("gMSDEy");
279   gMSDExt->Write("gMSDExt");
280   gMSDEyt->Write("gMSDEyt");
281   fInvBenMom->Write();
282   fBenMom->Write();
283   hFile->Close();
284   
285   cout << "... and closed!" << endl;
286   // Re Align
287   AliMUONGeometryTransformer *newTransform = alig->ReAlign(transform,parameters,true); 
288   newTransform->WriteTransformations("transform2ReAlign.dat");
289   
290   // Generate realigned data in local cdb
291   const TClonesArray* array = newTransform->GetMisAlignmentData();
292
293   // 100 mum residual resolution for chamber misalignments?
294   alig->SetAlignmentResolution(array,-1,0.01,0.01,0.004,0.003);
295    
296   // CDB manager
297   AliCDBManager* cdbManager = AliCDBManager::Instance();
298   cdbManager->SetDefaultStorage("local://ReAlignCDB");
299   
300   AliCDBMetaData* cdbData = new AliCDBMetaData();
301   cdbData->SetResponsible("Dimuon Offline project");
302   cdbData->SetComment("MUON alignment objects with residual misalignment");
303   AliCDBId id("MUON/Align/Data", 0, 0); 
304   cdbManager->Put(const_cast<TClonesArray*>(array), id, cdbData);
305
306 }
307