1 /**************************************************************************
2 * Copyright(c) 2004, 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 **************************************************************************/
16 // Thil class computes background corrections for the FMD. The correction is computed
17 // in eta,phi cells and the objects stored can be put into alien to use with analysis.
18 // It is based on the AliFMDInput class that is used to loop over hits and primaries.
20 // Author: Hans Hjersing Dalsgaard, NBI, hans.dalsgaard@cern.ch
24 #include "AliSimulation.h"
25 #include "TStopwatch.h"
28 #include "AliRunLoader.h"
29 #include "AliGeomManager.h"
30 #include "AliFMDGeometry.h"
32 #include "TParticle.h"
33 #include "TDatabasePDG.h"
34 #include "TParticlePDG.h"
36 #include "AliFMDBackgroundCorrection.h"
38 #include "AliCDBManager.h"
40 #include "TClonesArray.h"
42 #include "AliFMDHit.h"
43 #include "AliLoader.h"
46 #include "AliGenEventHeader.h"
47 #include "AliHeader.h"
51 #include "AliCDBMetaData.h"
53 #include "AliFMDParameters.h"
58 ClassImp(AliFMDBackgroundCorrection)
59 //_____________________________________________________________________
60 AliFMDBackgroundCorrection::AliFMDBackgroundCorrection() :
65 //_____________________________________________________________________
66 AliFMDBackgroundCorrection::AliFMDInputBG::AliFMDInputBG() :
83 //_____________________________________________________________________
85 void AliFMDBackgroundCorrection::GenerateBackgroundCorrection(Int_t nvtxbins,
91 const Char_t* filename,
95 TGrid::Connect("alien:",0,0,"t");
99 AliCDBManager::Instance()->SetDefaultStorage("alien://Folder=/alice/data/2008/LHC08d/OCDB/");
100 AliCDBManager::Instance()->SetRun(runNo);
102 #if defined(__CINT__)
103 gSystem->Load("liblhapdf");
104 gSystem->Load("libEGPythia6");
105 gSystem->Load("libpythia6");
106 gSystem->Load("libAliPythia6");
107 gSystem->Load("libgeant321");
111 //Setting up the geometry
112 //-----------------------------------------------
113 if (AliGeomManager::GetGeometry() == NULL)
114 AliGeomManager::LoadGeometry();
116 AliFMDGeometry* geo = AliFMDGeometry::Instance();
118 geo->InitTransformations();
122 AliInfo("Processing hits and primaries ");
125 input.SetVtxCutZ(zvtxcut);
126 input.SetNvtxBins(nvtxbins);
127 input.SetNbinsEta(nBinsEta);
130 AliInfo(Form("Found %d primaries and %d hits.", input.GetNprim(),input.GetNhits()));
132 TObjArray* hitArray = input.GetHits();
133 TObjArray* primaryArray = input.GetPrimaries();
134 fCorrectionArray.SetName("FMD_bg_correction");
135 fCorrectionArray.SetOwner();
136 for(Int_t det= 1; det <=3; det++) {
137 Int_t nRings = (det==1 ? 1 : 2);
139 TObjArray* detArrayCorrection = new TObjArray();
140 detArrayCorrection->SetName(Form("FMD%d",det));
141 fCorrectionArray.AddAtAndExpand(detArrayCorrection,det);
144 for(Int_t iring = 0; iring<nRings; iring++) {
145 TObjArray* primRingArray = (TObjArray*)primaryArray->At(iring);
146 Char_t ringChar = (iring == 0 ? 'I' : 'O');
147 TObjArray* vtxArrayCorrection = new TObjArray();
148 vtxArrayCorrection->SetName(Form("FMD%d%c",det,ringChar));
149 detArrayCorrection->AddAtAndExpand(vtxArrayCorrection,iring);
151 for(Int_t vertexBin=0;vertexBin<nvtxbins;vertexBin++) {
152 TObjArray* detArray = (TObjArray*)hitArray->At(det);
153 TObjArray* vtxArray = (TObjArray*)detArray->At(iring);
154 TH2F* hHits = (TH2F*)vtxArray->At(vertexBin);
155 TH2F* hPrimary = (TH2F*)primRingArray->At(vertexBin);
156 TH2F* hCorrection = (TH2F*)hHits->Clone(Form("FMD%d%c_vtxbin_%d_correction",det,ringChar,vertexBin));
157 hCorrection->Divide(hPrimary);
158 vtxArrayCorrection->AddAtAndExpand(hCorrection,vertexBin);
164 TAxis refAxis(nvtxbins,-1*zvtxcut,zvtxcut);
169 TFile fout(filename,"RECREATE");
170 refAxis.Write("vertexbins");
174 fout.mkdir("Primaries");
175 fout.cd("Primaries");
176 primaryArray->Write();
177 fout.mkdir("Correction");
178 fout.cd("Correction");
179 fCorrectionArray.Write();
182 AliCDBManager* cdb = AliCDBManager::Instance();
183 cdb->SetDefaultStorage("local://$ALICE_ROOT");
184 AliCDBId id("FMD/AnalysisCalib/Background",runNo,endRunNo);
186 AliCDBMetaData* meta = new AliCDBMetaData;
187 meta->SetResponsible(gSystem->GetUserInfo()->fRealName.Data());
188 meta->SetAliRootVersion(gROOT->GetVersion());
189 meta->SetBeamPeriod(1);
190 meta->SetComment("Background Correction for FMD");
191 meta->SetProperty("key1", &fout );
192 cdb->Put(&fout, id, meta);
200 //_____________________________________________________________________
201 void AliFMDBackgroundCorrection::Simulate(Int_t nEvents) {
205 TGrid::Connect("alien:",0,0,"t");
206 sim.SetDefaultStorage("alien://Folder=/alice/data/2008/LHC08d/OCDB/");
207 sim.SetConfigFile("Config.C");
208 sim.SetRunQA("FMD:");
211 sim.RunSimulation(nEvents);
217 //_____________________________________________________________________
218 Bool_t AliFMDBackgroundCorrection::AliFMDInputBG::ProcessTrack(Int_t i , TParticle* p, AliFMDHit* h) {
227 // if(h->Detector() == 3)
228 // std::cout<<"Track "<<i<<" hit "<<"FMD3"<<h->Ring()<<std::endl;
232 AliGenEventHeader* genHeader = fLoader->GetHeader()->GenEventHeader();
234 genHeader->PrimaryVertex(vertex);
236 if(TMath::Abs(vertex.At(2)) > fZvtxCut)
239 Double_t delta = 2*fZvtxCut/fNvtxBins;
240 Double_t vertexBinDouble = (vertex.At(2) + fZvtxCut) / delta;
241 Int_t vertexBin = (Int_t)vertexBinDouble;
244 if(h->Q() != 0 && (i != fPrevTrack || h->Detector() != fPrevDetector || h->Ring() != fPrevRing)) {
246 Int_t det = h->Detector();
247 Char_t ring = h->Ring();
248 // if(h->Detector() == 3)
249 // std::cout<<"Detected a hit in " <<"FMD"<<det<<ring<<" ! "<<std::endl;
251 Int_t iring = (ring == 'I' ? 0 : 1);
253 TObjArray* detArray = (TObjArray*)fHitArray.At(det);
254 TObjArray* vtxArray = (TObjArray*)detArray->At(iring);
255 TH2F* hHits = (TH2F*)vtxArray->At(vertexBin);
257 Float_t phi = TMath::ATan2(h->Py(),h->Px());
259 phi = phi+2*TMath::Pi();
260 Float_t theta = TMath::ATan2(TMath::Sqrt(TMath::Power(h->X(),2)+TMath::Power(h->Y(),2)),h->Z()+vertex.At(2));
261 Float_t eta = -1*TMath::Log(TMath::Tan(0.5*theta));
262 hHits->Fill(eta,phi);
268 if(p && i != fPrevTrack) {
270 TDatabasePDG* pdgDB = TDatabasePDG::Instance();
271 TParticlePDG* pdgPart = pdgDB->GetParticle(p->GetPdgCode());
272 Float_t charge = (pdgPart ? pdgPart->Charge() : 0);
273 Float_t phi = TMath::ATan2(p->Py(),p->Px());
276 phi = phi+2*TMath::Pi();
279 Bool_t primary = (charge!=0)&&(TMath::Abs(p->Vx() - vertex.At(0))<0.1)&&(TMath::Abs(p->Vy() - vertex.At(1))<0.1)&&(TMath::Abs(p->Vz() - vertex.At(2))<0.1);
282 TObjArray* innerArray = (TObjArray*)fPrimaryArray.At(0);
283 TObjArray* outerArray = (TObjArray*)fPrimaryArray.At(1);
285 TH2F* hPrimaryInner = (TH2F*)innerArray->At(vertexBin);
286 TH2F* hPrimaryOuter = (TH2F*)outerArray->At(vertexBin);
287 hPrimaryInner->Fill(p->Eta(),phi);
288 hPrimaryOuter->Fill(p->Eta(),phi);
297 //_____________________________________________________________________
298 Bool_t AliFMDBackgroundCorrection::AliFMDInputBG::Init() {
300 fPrimaryArray.SetOwner();
301 fPrimaryArray.SetName("FMD_primary");
303 for(Int_t iring = 0; iring<2;iring++) {
304 Char_t ringChar = (iring == 0 ? 'I' : 'O');
305 TObjArray* ringArray = new TObjArray();
306 ringArray->SetName(Form("FMD_%c",ringChar));
307 fPrimaryArray.AddAtAndExpand(ringArray,iring);
308 Int_t nSec = (iring == 1 ? 40 : 20);
309 for(Int_t v=0; v<fNvtxBins;v++) {
311 TH2F* hPrimary = new TH2F(Form("hPrimary_FMD_%c_vtx%d",ringChar,v),
312 Form("hPrimary_FMD_%c_vtx%d",ringChar,v),
313 fNbinsEta, -6,6, nSec, 0,2*TMath::Pi());
314 ringArray->AddAtAndExpand(hPrimary,v);
319 fHitArray.SetOwner();
320 fHitArray.SetName("FMD_hits");
322 for(Int_t det =1; det<=3;det++)
324 TObjArray* detArrayHits = new TObjArray();
325 detArrayHits->SetName(Form("FMD%d",det));
326 fHitArray.AddAtAndExpand(detArrayHits,det);
327 Int_t nRings = (det==1 ? 1 : 2);
328 for(Int_t ring = 0;ring<nRings;ring++)
330 Int_t nSec = (ring == 1 ? 40 : 20);
331 Char_t ringChar = (ring == 0 ? 'I' : 'O');
332 TObjArray* vtxArrayHits = new TObjArray();
333 vtxArrayHits->SetName(Form("FMD%d%c",det,ringChar));
334 detArrayHits->AddAtAndExpand(vtxArrayHits,ring);
335 for(Int_t v=0; v<fNvtxBins;v++)
338 TH2F* hHits = new TH2F(Form("hHits_FMD%d%c_vtx%d",det,ringChar,v),
339 Form("hHits_FMD%d%c_vtx%d",det,ringChar,v),
340 fNbinsEta, -6,6, nSec, 0, 2*TMath::Pi());
341 vtxArrayHits->AddAtAndExpand(hHits,v);