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