]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/MUONAlignment.C
New classes for shuttle (Laurent)
[u/mrichter/AliRoot.git] / MUON / MUONAlignment.C
CommitLineData
af9df62e 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
37void MUONAlignment(Int_t nEvents = 100000, TString fileName = "galice.root", TString fileList = "")
38{
39
40 Double_t parameters[3*156];
41 Double_t errors[3*156];
42 Double_t pulls[3*156];
43 for(Int_t k=0;k<3*156;k++) {
44 parameters[k]=0.;
45 errors[k]=0.;
46 pulls[k]=0.;
47 }
48
49 Double_t trackParams[8] = {0.,0.,0.,0.,0.,0.,0.,0.};
50
51 Int_t iPar = 0;
52 // Set initial values here, good guess may help convergence
53 // St 1
54 // parameters[iPar++] = 0.010300 ; parameters[iPar++] = 0.010600 ; parameters[iPar++] = 0.000396 ;
55
56 bool bLoop = kFALSE;
57 ifstream sFileList;
58 if (fileList.Contains(".list")) {
59 cout << "Reading file list: " << fileList.Data() << endl;
60 bLoop = kTRUE;
61
62 TString fullListName(".");
63 fullListName +="/";
64 fullListName +=fileList;
65 sFileList.open(fileList.Data());
66 }
67
68 TH1F *fInvBenMom = new TH1F("fInvBenMom","fInvBenMom",200,-0.1,0.1);
69 TH1F *fBenMom = new TH1F("fBenMom","fBenMom",200,-40,40);
70
71 AliMUONAlignment* alig = new AliMUONAlignment();
72 alig->InitGlobalParameters(parameters);
73
74 AliMUONGeometryTransformer *transform = new AliMUONGeometryTransformer(true);
75 transform->ReadGeometryData("volpath.dat", "transform.dat");
76 alig->SetGeometryTransformer(transform);
77
78 char cFileName[100];
79 AliMUONDataInterface amdi;
80
81 Int_t lMaxFile = 1000;
82 Int_t iFile = 0;
83 bool bKeepLoop = kTRUE;
84 while(bKeepLoop && iFile<lMaxFile){
85 iFile++;
86 if (bLoop) {
87 sFileList.getline(cFileName,100);
88 if (sFileList.eof()) bKeepLoop = kFALSE;
89 }
90 else {
91 sprintf(cFileName,fileName.Data());
92 bKeepLoop = kFALSE;
93 }
94 if (!strstr(cFileName,"galice.root")) continue;
95 cout << "Using file: " << cFileName << endl;
96 amdi.SetFile(cFileName);
97 Int_t nevents = amdi.NumberOfEvents();
98 for(Int_t event = 0; event < nevents; event++) {
99 amdi.GetEvent(event);
100 Int_t ntracks = amdi.NumberOfRecTracks();
101 cout << " there are " << ntracks << " tracks in event " << event << endl;
102 Int_t iTrack=0;
103 Int_t iTrackOk=0;
104 AliMUONTrack* track = amdi.RecTrack(iTrack);
105 while(track) {
106 Double_t invBenMom = track->GetTrackParamAtVertex()->GetInverseBendingMomentum();
107 fInvBenMom->Fill(invBenMom);
108 fBenMom->Fill(1./invBenMom);
109 if (TMath::Abs(invBenMom)<=1.04) {
110 cout << "Track " << iTrack << endl;
111 alig->ProcessTrack(track);
112 cout << "Calling LocalFit" << endl;
113 alig->LocalFit(iTrackOk++,trackParams,0);
114 }
115 iTrack++;
116 track = amdi.RecTrack(iTrack);
117 }
118 }
119 }
120 alig->GlobalFit(parameters,errors,pulls);
121
122 cout << "Done with GlobalFit " << endl;
123
124 // Store results
125 Double_t DEid[156] = {0};
126 Double_t MSDEx[156] = {0};
127 Double_t MSDEy[156] = {0};
128 Double_t MSDExt[156] = {0};
129 Double_t MSDEyt[156] = {0};
130 Double_t DEidErr[156] = {0};
131 Double_t MSDExErr[156] = {0};
132 Double_t MSDEyErr[156] = {0};
133 Double_t MSDExtErr[156] = {0};
134 Double_t MSDEytErr[156] = {0};
135 Int_t lNDetElem = 4*2+4*2+18*2+26*2+26*2;
136 Int_t lNDetElemCh[10] = {4,4,4,4,18,18,26,26,26,26};
137 Int_t lSNDetElemCh[10] = {4,8,12,16,34,52,78,104,130,156};
138 Int_t idOffset = 0; // 400
139 Int_t lSDetElemCh = 0;
140 for(Int_t iDE=0; iDE<lNDetElem; iDE++){
141 DEidErr[iDE] = 0.;
142 DEid[iDE] = idOffset+100;
143 DEid[iDE] += iDE;
144 lSDetElemCh = 0;
145 for(Int_t iCh=0; iCh<9; iCh++){
146 lSDetElemCh += lNDetElemCh[iCh];
147 if (iDE>=lSDetElemCh) {
148 DEid[iDE] += 100;
149 DEid[iDE] -= lNDetElemCh[iCh];
150 }
151 }
152 MSDEx[iDE]=parameters[3*iDE+0];
153 MSDEy[iDE]=parameters[3*iDE+1];
154 MSDExt[iDE]=parameters[3*iDE+2];
155 MSDEyt[iDE]=parameters[3*iDE+2];
156 MSDExErr[iDE]=(Double_t)alig->GetParError(3*iDE+0);
157 MSDEyErr[iDE]=(Double_t)alig->GetParError(3*iDE+1);
158 MSDExtErr[iDE]=(Double_t)alig->GetParError(3*iDE+2);
159 MSDEytErr[iDE]=(Double_t)alig->GetParError(3*iDE+2);
160 }
161
162 cout << "Let's create graphs ... " << endl;
163
164 TGraphErrors *gMSDEx = new TGraphErrors(lNDetElem,DEid,MSDEx,DEidErr,MSDExErr);
165 TGraphErrors *gMSDEy = new TGraphErrors(lNDetElem,DEid,MSDEy,DEidErr,MSDEyErr);
166 TGraphErrors *gMSDExt = new TGraphErrors(lNDetElem,DEid,MSDExt,DEidErr,MSDExtErr);
167 TGraphErrors *gMSDEyt = new TGraphErrors(lNDetElem,DEid,MSDEyt,DEidErr,MSDEytErr);
168
169 cout << "... graphs created, open file ... " << endl;
170
171 TFile *hFile = new TFile("measShifts.root","RECREATE");
172
173 cout << "... file opened ... " << endl;
174
175 gMSDEx->Write("gMSDEx");
176 gMSDEy->Write("gMSDEy");
177 gMSDExt->Write("gMSDExt");
178 gMSDEyt->Write("gMSDEyt");
179 fInvBenMom->Write();
180 fBenMom->Write();
181 hFile->Close();
182
183 cout << "... and closed!" << endl;
184 // Re Align
185 AliMUONGeometryTransformer *newTransform = alig->ReAlign(transform,parameters,true);
186 newTransform->WriteTransformations("transform2ReAlign.dat");
187
188 // Generate realigned data in local cdb
189 TClonesArray* array = newTransform->GetMisAlignmentData();
190
191 // CDB manager
192 AliCDBManager* cdbManager = AliCDBManager::Instance();
193 cdbManager->SetDefaultStorage("local://ReAlignCDB");
194
195 AliCDBMetaData* cdbData = new AliCDBMetaData();
196 cdbData->SetResponsible("Dimuon Offline project");
197 cdbData->SetComment("MUON alignment objects with residual misalignment");
198 AliCDBId id("MUON/Align/Data", 0, 0);
199 cdbManager->Put(array, id, cdbData);
200
201}