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