1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
19 /// \file MUONAlignment.C
20 /// \brief Macro for MUON alignment using physics tracks.
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
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.
37 /// \author: B. Becker and J. Castillo
39 #if !defined(__CINT__) || defined(__MAKECINT__)
41 #include "AliMUONAlignment.h"
42 #include "AliMUONTrack.h"
43 #include "AliMUONTrackExtrap.h"
44 #include "AliMUONTrackParam.h"
45 #include "AliMUONGeometryTransformer.h"
46 #include "AliMUONESDInterface.h"
48 #include "AliESDEvent.h"
49 #include "AliESDMuonTrack.h"
51 #include "AliTracker.h"
52 #include "AliCDBManager.h"
53 #include "AliCDBMetaData.h"
55 #include "AliGeomManager.h"
60 #include <TGraphErrors.h>
63 #include <TClonesArray.h>
64 #include <Riostream.h>
70 void MUONAlignment(Int_t nEvents = 100000, char* geoFilename = "geometry.root", TString esdFileName = "AliESDs.root", TString fileList = "")
73 // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
74 if ( ! AliGeomManager::GetGeometry() ) {
75 AliGeomManager::LoadGeometry(geoFilename);
76 if (! AliGeomManager::GetGeometry() ) {
77 Error("MUONAlignment", "getting geometry from file %s failed", geoFilename);
83 // waiting for mag field in CDB
84 if (!TGeoGlobalMagField::Instance()->GetField()) {
85 printf("Loading field map...\n");
86 AliMagF* field = new AliMagF("Maps","Maps",2,1.,1., 10.,AliMagF::k5kG);
87 TGeoGlobalMagField::Instance()->SetField(field);
89 // set the magnetic field for track extrapolations
90 AliMUONTrackExtrap::SetField();
92 Double_t parameters[3*156];
93 Double_t errors[3*156];
94 Double_t pulls[3*156];
95 for(Int_t k=0;k<3*156;k++) {
101 Double_t trackParams[8] = {0.,0.,0.,0.,0.,0.,0.,0.};
103 // Set initial values here, good guess may help convergence
106 // parameters[iPar++] = 0.010300 ; parameters[iPar++] = 0.010600 ; parameters[iPar++] = 0.000396 ;
110 if (fileList.Contains(".list")) {
111 cout << "Reading file list: " << fileList.Data() << endl;
114 TString fullListName(".");
116 fullListName +=fileList;
117 sFileList.open(fileList.Data());
120 TH1F *fInvBenMom = new TH1F("fInvBenMom","fInvBenMom",200,-0.1,0.1);
121 TH1F *fBenMom = new TH1F("fBenMom","fBenMom",200,-40,40);
123 AliMUONAlignment* alig = new AliMUONAlignment();
124 alig->InitGlobalParameters(parameters);
126 AliMUONGeometryTransformer *transform = new AliMUONGeometryTransformer();
127 transform->LoadGeometryData();
128 alig->SetGeometryTransformer(transform);
130 // Do alignment with magnetic field off
131 alig->SetBFieldOn(kFALSE);
133 // Set tracking station to use
134 Bool_t bStOnOff[5] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
135 Bool_t bChOnOff[10] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
137 // Set degrees of freedom to align (see AliMUONAlignment)
138 alig->AllowVariations(bChOnOff);
140 // Fix parameters or add constraints here
141 // for (Int_t iSt=0; iSt<5; iSt++)
142 // if (!bStOnOff[iSt]) alig->FixStation(iSt+1);
143 for (Int_t iCh=0; iCh<10; iCh++)
144 if (!bChOnOff[iCh]) alig->FixChamber(iCh+1);
146 // Left and right sides of the detector are independent, one can choose to align
148 Bool_t bSpecLROnOff[2] = {kTRUE,kTRUE};
149 alig->FixHalfSpectrometer(bChOnOff,bSpecLROnOff);
151 alig->SetChOnOff(bChOnOff);
152 alig->SetSpecLROnOff(bChOnOff);
154 // Set predifined global constrains: X, Y, P, XvsZ, YvsZ, PvsZ, XvsY, YvsY, PvsY
155 Bool_t bVarXYT[9] = {kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE,kTRUE};
156 Bool_t bDetTLBR[4] = {kFALSE,kTRUE,kFALSE,kTRUE};
157 alig->AddConstraints(bChOnOff,bVarXYT,bDetTLBR,bSpecLROnOff);
162 Int_t lMaxFile = 1000;
165 bool bKeepLoop = kTRUE;
169 while(bKeepLoop && iFile<lMaxFile){
172 sFileList.getline(cFileName,100);
173 if (sFileList.eof()) bKeepLoop = kFALSE;
176 sprintf(cFileName,esdFileName.Data());
179 if (!strstr(cFileName,"AliESDs.root")) continue;
180 cout << "Using file: " << cFileName << endl;
183 TFile* esdFile = TFile::Open(cFileName); // open the file
184 if (!esdFile || !esdFile->IsOpen()) {
185 cout << "opening ESD file " << cFileName << "failed" << endl;
188 TTree* esdTree = (TTree*) esdFile->Get("esdTree"); // get the tree
190 cout << "no ESD tree found" << endl;
194 AliESDEvent* esdEvent = new AliESDEvent(); // link ESD event to the tree
195 esdEvent->ReadFromTree(esdTree);
197 Int_t nevents = esdTree->GetEntries();
198 cout << "... with " << nevents << endl;
199 for(Int_t event = 0; event < nevents; event++) {
200 if (iEvent >= nEvents){
206 if (esdTree->GetEvent(event) <= 0) {
207 cout << "fails to read ESD object for event " << event << endl;
211 Int_t nTracks = Int_t(esdEvent->GetNumberOfMuonTracks());
212 if (!event%100) cout << " there are " << nTracks << " tracks in event " << event << endl;
213 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
214 AliESDMuonTrack* esdTrack = esdEvent->GetMuonTrack(iTrack);
215 if (!esdTrack->ClustersStored()) continue;
216 Double_t invBenMom = esdTrack->GetInverseBendingMomentum();
217 fInvBenMom->Fill(invBenMom);
218 fBenMom->Fill(1./invBenMom);
219 if (TMath::Abs(invBenMom)<=1.04) {
221 AliMUONESDInterface::ESDToMUON(*esdTrack, track);
222 alig->ProcessTrack(&track);
223 alig->LocalFit(iTrackOk++,trackParams,0);
230 cout << "Processed " << iTrackTot << " Tracks so far." << endl;
232 alig->GlobalFit(parameters,errors,pulls);
234 cout << "Done with GlobalFit " << endl;
237 Double_t DEid[156] = {0};
238 Double_t MSDEx[156] = {0};
239 Double_t MSDEy[156] = {0};
240 Double_t MSDExt[156] = {0};
241 Double_t MSDEyt[156] = {0};
242 Double_t DEidErr[156] = {0};
243 Double_t MSDExErr[156] = {0};
244 Double_t MSDEyErr[156] = {0};
245 Double_t MSDExtErr[156] = {0};
246 Double_t MSDEytErr[156] = {0};
247 Int_t lNDetElem = 4*2+4*2+18*2+26*2+26*2;
248 Int_t lNDetElemCh[10] = {4,4,4,4,18,18,26,26,26,26};
249 // Int_t lSNDetElemCh[10] = {4,8,12,16,34,52,78,104,130,156};
250 Int_t idOffset = 0; // 400
251 Int_t lSDetElemCh = 0;
252 for(Int_t iDE=0; iDE<lNDetElem; iDE++){
254 DEid[iDE] = idOffset+100;
257 for(Int_t iCh=0; iCh<9; iCh++){
258 lSDetElemCh += lNDetElemCh[iCh];
259 if (iDE>=lSDetElemCh) {
261 DEid[iDE] -= lNDetElemCh[iCh];
264 MSDEx[iDE]=parameters[3*iDE+0];
265 MSDEy[iDE]=parameters[3*iDE+1];
266 MSDExt[iDE]=parameters[3*iDE+2];
267 MSDEyt[iDE]=parameters[3*iDE+2];
268 MSDExErr[iDE]=(Double_t)alig->GetParError(3*iDE+0);
269 MSDEyErr[iDE]=(Double_t)alig->GetParError(3*iDE+1);
270 MSDExtErr[iDE]=(Double_t)alig->GetParError(3*iDE+2);
271 MSDEytErr[iDE]=(Double_t)alig->GetParError(3*iDE+2);
274 cout << "Let's create graphs ... " << endl;
276 TGraphErrors *gMSDEx = new TGraphErrors(lNDetElem,DEid,MSDEx,DEidErr,MSDExErr);
277 TGraphErrors *gMSDEy = new TGraphErrors(lNDetElem,DEid,MSDEy,DEidErr,MSDEyErr);
278 TGraphErrors *gMSDExt = new TGraphErrors(lNDetElem,DEid,MSDExt,DEidErr,MSDExtErr);
279 TGraphErrors *gMSDEyt = new TGraphErrors(lNDetElem,DEid,MSDEyt,DEidErr,MSDEytErr);
281 cout << "... graphs created, open file ... " << endl;
283 TFile *hFile = new TFile("measShifts.root","RECREATE");
285 cout << "... file opened ... " << endl;
287 gMSDEx->Write("gMSDEx");
288 gMSDEy->Write("gMSDEy");
289 gMSDExt->Write("gMSDExt");
290 gMSDEyt->Write("gMSDEyt");
295 cout << "... and closed!" << endl;
297 AliMUONGeometryTransformer *newTransform = alig->ReAlign(transform,parameters,true);
298 newTransform->WriteTransformations("transform2ReAlign.dat");
300 // Generate realigned data in local cdb
301 const TClonesArray* array = newTransform->GetMisAlignmentData();
303 // 100 mum residual resolution for chamber misalignments?
304 alig->SetAlignmentResolution(array,-1,0.01,0.01,0.004,0.003);
307 AliCDBManager* cdbManager = AliCDBManager::Instance();
308 cdbManager->SetDefaultStorage("local://ReAlignCDB");
310 AliCDBMetaData* cdbData = new AliCDBMetaData();
311 cdbData->SetResponsible("Dimuon Offline project");
312 cdbData->SetComment("MUON alignment objects with residual misalignment");
313 AliCDBId id("MUON/Align/Data", 0, AliCDBRunRange::Infinity());
314 cdbManager->Put(const_cast<TClonesArray*>(array), id, cdbData);