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