]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONAlignment.C
effc++ warnings
[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 "AliESDEvent.h"
47 #include "AliESDMuonTrack.h"
48
49 #include "AliMagFMaps.h"
50 #include "AliTracker.h"
51 #include "AliCDBManager.h"
52 #include "AliCDBMetaData.h"
53 #include "AliCDBId.h"
54 #include "AliGeomManager.h"
55
56 #include <TString.h>
57 #include <TError.h>
58 #include <TH1.h>
59 #include <TGraphErrors.h>
60 #include <TFile.h>
61 #include <TTree.h>
62 #include <TClonesArray.h>
63 #include <Riostream.h>
64
65 #include <fstream>
66
67 #endif
68
69 void MUONAlignment(Int_t nEvents = 100000, char* geoFilename = "geometry.root", TString esdFileName = "AliESDs.root", TString fileList = "")
70 {
71  
72   // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
73   if ( ! AliGeomManager::GetGeometry() ) {
74     AliGeomManager::LoadGeometry(geoFilename);
75     if (! AliGeomManager::GetGeometry() ) {
76       Error("MUONAlignment", "getting geometry from file %s failed", geoFilename);
77       return;
78     }
79   }
80   
81   // set  mag field 
82   // waiting for mag field in CDB 
83   printf("Loading field map...\n");
84   AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG);
85   AliTracker::SetFieldMap(field, kFALSE);
86   // set the magnetic field for track extrapolations
87   AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap());
88
89   Double_t parameters[3*156];
90   Double_t errors[3*156];
91   Double_t pulls[3*156];
92   for(Int_t k=0;k<3*156;k++) {
93     parameters[k]=0.;
94     errors[k]=0.;
95     pulls[k]=0.;
96   }
97
98   Double_t trackParams[8] = {0.,0.,0.,0.,0.,0.,0.,0.};
99
100   // Set initial values here, good guess may help convergence
101   // St 1 
102   //  Int_t iPar = 0;
103   //  parameters[iPar++] =  0.010300 ;  parameters[iPar++] =  0.010600 ;  parameters[iPar++] =  0.000396 ;  
104
105   bool bLoop = kFALSE;
106   ifstream sFileList;
107   if (fileList.Contains(".list")) {
108     cout << "Reading file list: " << fileList.Data() << endl;
109     bLoop = kTRUE;
110
111     TString fullListName(".");
112     fullListName +="/";
113     fullListName +=fileList;
114     sFileList.open(fileList.Data());
115   }
116   
117   TH1F *fInvBenMom = new TH1F("fInvBenMom","fInvBenMom",200,-0.1,0.1); 
118   TH1F *fBenMom = new TH1F("fBenMom","fBenMom",200,-40,40); 
119
120   AliMUONAlignment* alig = new AliMUONAlignment();
121   alig->InitGlobalParameters(parameters);
122
123   AliMUONGeometryTransformer *transform = new AliMUONGeometryTransformer();
124   transform->LoadGeometryData();
125   alig->SetGeometryTransformer(transform);
126
127   // Do alignment with magnetic field off
128   alig->SetBFieldOn(kFALSE);
129   
130   // Set tracking station to use
131   Bool_t bStOnOff[5] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
132
133   // Fix parameters or add constraints here
134   for (Int_t iSt=0; iSt<5; iSt++)
135     if (!bStOnOff[iSt]) alig->FixStation(iSt+1);
136
137   // Left and right sides of the detector are independent, one can choose to align 
138   // only one side
139   Bool_t bSpecLROnOff[2] = {kTRUE,kTRUE};
140   alig->FixHalfSpectrometer(bStOnOff,bSpecLROnOff);
141
142   // Set predifined global constrains: X, Y, P, XvsZ, YvsZ, PvsZ, XvsY, YvsY, PvsY
143   Bool_t bVarXYT[9] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
144   Bool_t bDetTLBR[4] = {kFALSE,kTRUE,kFALSE,kTRUE};
145   alig->AddConstraints(bStOnOff,bVarXYT,bDetTLBR,bSpecLROnOff);
146
147
148   char cFileName[100];  
149     
150   Int_t lMaxFile = 1000;
151   Int_t iFile = 0;
152   Int_t iEvent = 0;
153   bool bKeepLoop = kTRUE;
154   Int_t iTrackTot=0;
155   Int_t iTrackOk=0;
156
157   while(bKeepLoop && iFile<lMaxFile){
158     iFile++;
159     if (bLoop) {
160       sFileList.getline(cFileName,100);
161       if (sFileList.eof()) bKeepLoop = kFALSE;
162     }
163     else {
164       sprintf(cFileName,esdFileName.Data());
165       bKeepLoop = kFALSE;
166     }
167     if (!strstr(cFileName,"AliESDs.root")) continue;      
168     cout << "Using file: " << cFileName << endl;
169     
170     // load ESD event
171     TFile* esdFile = TFile::Open(cFileName); // open the file
172     if (!esdFile || !esdFile->IsOpen()) {
173       cout << "opening ESD file " << cFileName << "failed" << endl;
174       continue;
175     }
176     TTree* esdTree = (TTree*) esdFile->Get("esdTree"); // get the tree
177     if (!esdTree) {
178       cout << "no ESD tree found" << endl;
179       esdFile->Close();
180       continue;
181     }
182     AliESDEvent* esdEvent = new AliESDEvent(); // link ESD event to the tree
183     esdEvent->ReadFromTree(esdTree);
184
185     Int_t nevents = esdTree->GetEntries();
186     cout << "... with " << nevents << endl;
187     for(Int_t event = 0; event < nevents; event++) {
188       if (iEvent >= nEvents){
189         bKeepLoop = kFALSE;
190         break;
191       }
192       iEvent++;
193
194       if (esdTree->GetEvent(event) <= 0) {
195         cout << "fails to read ESD object for event " << event << endl;
196         continue;
197       }
198
199       Int_t nTracks = Int_t(esdEvent->GetNumberOfMuonTracks());
200       if (!event%100) cout << " there are " << nTracks << " tracks in event " << event << endl;
201       for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
202         AliESDMuonTrack* esdTrack = esdEvent->GetMuonTrack(iTrack);
203         if (!esdTrack->ClustersStored()) continue;
204         Double_t invBenMom = esdTrack->GetInverseBendingMomentum();
205         fInvBenMom->Fill(invBenMom);
206         fBenMom->Fill(1./invBenMom);
207         if (TMath::Abs(invBenMom)<=1.04) {
208           AliMUONTrack track(*esdTrack);
209           alig->ProcessTrack(&track);
210           alig->LocalFit(iTrackOk++,trackParams,0);
211         }
212         iTrackTot++;
213       }
214     }
215     delete esdEvent;
216     esdFile->Close();
217     cout << "Processed " << iTrackTot << " Tracks so far." << endl;
218   }
219   alig->GlobalFit(parameters,errors,pulls);
220
221   cout << "Done with GlobalFit " << endl;
222
223   // Store results
224   Double_t DEid[156] = {0};
225   Double_t MSDEx[156] = {0};
226   Double_t MSDEy[156] = {0};
227   Double_t MSDExt[156] = {0};
228   Double_t MSDEyt[156] = {0};
229   Double_t DEidErr[156] = {0};
230   Double_t MSDExErr[156] = {0};
231   Double_t MSDEyErr[156] = {0};
232   Double_t MSDExtErr[156] = {0};
233   Double_t MSDEytErr[156] = {0};
234   Int_t lNDetElem = 4*2+4*2+18*2+26*2+26*2;
235   Int_t lNDetElemCh[10] = {4,4,4,4,18,18,26,26,26,26};
236   // Int_t lSNDetElemCh[10] = {4,8,12,16,34,52,78,104,130,156};
237   Int_t idOffset = 0; // 400
238   Int_t lSDetElemCh = 0;
239   for(Int_t iDE=0; iDE<lNDetElem; iDE++){
240     DEidErr[iDE] = 0.;
241     DEid[iDE] = idOffset+100;
242     DEid[iDE] += iDE; 
243     lSDetElemCh = 0;
244     for(Int_t iCh=0; iCh<9; iCh++){
245       lSDetElemCh += lNDetElemCh[iCh];
246       if (iDE>=lSDetElemCh) {
247         DEid[iDE] += 100;
248         DEid[iDE] -= lNDetElemCh[iCh];
249       }
250     }
251     MSDEx[iDE]=parameters[3*iDE+0];
252     MSDEy[iDE]=parameters[3*iDE+1];
253     MSDExt[iDE]=parameters[3*iDE+2];
254     MSDEyt[iDE]=parameters[3*iDE+2];
255     MSDExErr[iDE]=(Double_t)alig->GetParError(3*iDE+0);
256     MSDEyErr[iDE]=(Double_t)alig->GetParError(3*iDE+1);
257     MSDExtErr[iDE]=(Double_t)alig->GetParError(3*iDE+2);
258     MSDEytErr[iDE]=(Double_t)alig->GetParError(3*iDE+2);
259   }
260
261   cout << "Let's create graphs ...  " << endl;
262
263   TGraphErrors *gMSDEx = new TGraphErrors(lNDetElem,DEid,MSDEx,DEidErr,MSDExErr); 
264   TGraphErrors *gMSDEy = new TGraphErrors(lNDetElem,DEid,MSDEy,DEidErr,MSDEyErr); 
265   TGraphErrors *gMSDExt = new TGraphErrors(lNDetElem,DEid,MSDExt,DEidErr,MSDExtErr); 
266   TGraphErrors *gMSDEyt = new TGraphErrors(lNDetElem,DEid,MSDEyt,DEidErr,MSDEytErr); 
267
268   cout << "... graphs created, open file ...  " << endl;
269
270   TFile *hFile = new TFile("measShifts.root","RECREATE");
271
272   cout << "... file opened ...  " << endl;
273
274   gMSDEx->Write("gMSDEx");
275   gMSDEy->Write("gMSDEy");
276   gMSDExt->Write("gMSDExt");
277   gMSDEyt->Write("gMSDEyt");
278   fInvBenMom->Write();
279   fBenMom->Write();
280   hFile->Close();
281   
282   cout << "... and closed!" << endl;
283   // Re Align
284   AliMUONGeometryTransformer *newTransform = alig->ReAlign(transform,parameters,true); 
285   newTransform->WriteTransformations("transform2ReAlign.dat");
286   
287   // Generate realigned data in local cdb
288   const TClonesArray* array = newTransform->GetMisAlignmentData();
289    
290   // CDB manager
291   AliCDBManager* cdbManager = AliCDBManager::Instance();
292   cdbManager->SetDefaultStorage("local://ReAlignCDB");
293   
294   AliCDBMetaData* cdbData = new AliCDBMetaData();
295   cdbData->SetResponsible("Dimuon Offline project");
296   cdbData->SetComment("MUON alignment objects with residual misalignment");
297   AliCDBId id("MUON/Align/Data", 0, 0); 
298   cdbManager->Put(const_cast<TClonesArray*>(array), id, cdbData);
299
300 }
301