ATLAS Offline Software
Loading...
Searching...
No Matches
CheckFlow_New.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// File: Generators/FlowAfterburnber/CheckFlow_New.h
6// Description:
7// This is a simple algorithm to histogram particle properties
8// for diagnosing of flow generation
9//
10//************************** THIS PROGRAM ANALYZES FILES WITH FIXED b_imp*********************
11//
12// AuthorList:
13// Andrzej Olszewski: Initial Code February 2006
14// Andrzej Olszewski: Converted to ROOT histograms July 2007
15// Soumya Mohapatra : Re-written to check the new Flow implementations (JUNE 2011)
16
19//
21
22#include "GaudiKernel/SmartDataPtr.h"
23#include "GaudiKernel/DataSvc.h"
24
25#include "GaudiKernel/ITHistSvc.h"
26
27#include "TH1.h"
28#include "TProfile.h"
29
30#include "AtlasHepMC/GenEvent.h"
33
35#include <cmath>
36
38CheckFlow_New::CheckFlow_New(const std::string& name, ISvcLocator* pSvcLocator) :
39 AthAlgorithm(name, pSvcLocator)
40{
41}
42
43
45 ATH_CHECK(m_hijingKey.initialize());
46
47 ATH_MSG_INFO(">>> CheckFlow_New from Initialize");
48
49 float pt_binvals[n_ptbin+1]={0.0,0.25,0.5,1.0,1.5,2.0,2.5,3.0,3.5,4.0,6.0,8.0,12.0,16.0,20.0,30.0,40.0};
50 float eta_bin_max = 4.0;
51
52
53 SmartIF<ITHistSvc> rootHistSvc{Gaudi::svcLocator()->service("THistSvc")};
54 if (!rootHistSvc) {
55 ATH_MSG_ERROR( "Unable to locate THistSvc" );
56 return StatusCode::FAILURE;
57 }
58
59 std::string histPath = "/FlowOutPut/";
60 char name[100],name1[100];
61 for (int ihar=0;ihar<6;ihar++){
62
63// /*
64 sprintf(name,"hist_Psi_%d_true",ihar+1);
65 sprintf(name1,"Truth Psi_{%d} distribution;%dPsi_{%d} Truth;events",ihar+1,ihar+1,ihar+1);
66 m_hist_Psi_n_true [ihar]=new TH1D (name,name1,1000,-M_PI,M_PI);
67 if ( rootHistSvc->regHist(histPath+m_hist_Psi_n_true[ihar]->GetName(),m_hist_Psi_n_true[ihar]).isFailure() ){
68 msg(MSG::WARNING) << "Can't book "<< histPath+m_hist_Psi_n_true[ihar]->GetName() << endmsg;
69 }
70 m_hist_Psi_n_true[ihar]->GetXaxis()->CenterTitle();
71 m_hist_Psi_n_true[ihar]->GetYaxis()->CenterTitle();
72
73 sprintf(name,"hist_Psi_%d_reco",ihar+1);
74 sprintf(name1,"Reconstructed Psi_{%d} distribution;%dPsi_{%d} Reco;events",ihar+1,ihar+1,ihar+1);
75 m_hist_Psi_n_reco [ihar]=new TH1D (name,name1,1000,-M_PI,M_PI);
76 if ( rootHistSvc->regHist(histPath+m_hist_Psi_n_reco[ihar]->GetName(),m_hist_Psi_n_reco[ihar]).isFailure() ){
77 msg(MSG::WARNING) << "Can't book "<< histPath+m_hist_Psi_n_reco[ihar]->GetName() << endmsg;
78 }
79 m_hist_Psi_n_reco[ihar]->GetXaxis()->CenterTitle();
80 m_hist_Psi_n_reco[ihar]->GetYaxis()->CenterTitle();
81
82
83 for (int ihar2=0;ihar2<6;ihar2++)
84 {
85 int ihar_i=ihar*6+ihar2;
86
87 sprintf(name,"hist_Psi_corr_true_%d_%d",ihar+1,ihar2+1);
88 sprintf(name1,"true Psi_{%d} -Psi_{%d};%dPsi_{%d} -%dPsi_{%d} ;events",ihar+1,ihar2+1,ihar+1,ihar+1,ihar2+1,ihar2+1);
89 m_hist_psi_corr_true [ihar_i]=new TH1D (name,name1,1000,-2*M_PI,2*M_PI);
90 if ( rootHistSvc->regHist(histPath+m_hist_psi_corr_true[ihar_i]->GetName(),m_hist_psi_corr_true[ihar_i]).isFailure() ){
91 msg(MSG::WARNING) << "Can't book "<< histPath+m_hist_psi_corr_true[ihar_i]->GetName() << endmsg;
92 }
93 m_hist_psi_corr_true[ihar_i]->GetXaxis()->CenterTitle();
94 m_hist_psi_corr_true[ihar_i]->GetYaxis()->CenterTitle();
95
96
97 sprintf(name,"hist_Psi_corr_reco_%d_%d",ihar+1,ihar2+1);
98 sprintf(name1,"reco Psi_{%d} -Psi_{%d};%dPsi_{%d} -%dPsi_{%d} ;events",ihar+1,ihar2+1,ihar+1,ihar+1,ihar2+1,ihar2+1);
99 m_hist_psi_corr_reco [ihar_i]=new TH1D (name,name1,1000,-2*M_PI,2*M_PI);
100 if ( rootHistSvc->regHist(histPath+m_hist_psi_corr_reco[ihar_i]->GetName(),m_hist_psi_corr_reco[ihar_i]).isFailure() ){
101 msg(MSG::WARNING) << "Can't book "<< histPath+m_hist_psi_corr_reco[ihar_i]->GetName() << endmsg;
102 }
103 m_hist_psi_corr_reco[ihar_i]->GetXaxis()->CenterTitle();
104 m_hist_psi_corr_reco[ihar_i]->GetYaxis()->CenterTitle();
105 }
106
107
108
109 //integrated vn event by event
110 sprintf(name,"hist_v%d_ebe",ihar+1);
111 sprintf(name1,"v%d;v%d;events",ihar+1,ihar+1);
112 m_hist_vn_ebe [ihar]=new TH1D (name,name1,1000,-0.5,0.5);
113 if ( rootHistSvc->regHist(histPath+m_hist_vn_ebe[ihar]->GetName(),m_hist_vn_ebe[ihar]).isFailure() ){
114 msg(MSG::WARNING) << "Can't book "<< histPath+m_hist_vn_ebe[ihar]->GetName() << endmsg;
115 }
116 m_hist_vn_ebe[ihar]->GetXaxis()->CenterTitle();
117 m_hist_vn_ebe[ihar]->GetYaxis()->CenterTitle();
118
119 sprintf(name ,"hist_Psi%d_ebe",ihar+1);
120 sprintf(name1,"%d#Delta#Psi;%d(#Psi_{reco}-#Psi_{Truth});events",ihar+1,ihar+1);
121 m_hist_Psi_n_ebe [ihar]=new TH1D (name,name1,1000,-M_PI,M_PI);
122 if ( rootHistSvc->regHist(histPath+m_hist_Psi_n_ebe[ihar]->GetName(),m_hist_Psi_n_ebe[ihar]).isFailure() ){
123 msg(MSG::WARNING) << "Can't book "<< histPath+m_hist_Psi_n_ebe[ihar]->GetName() << endmsg;
124 }
125 m_hist_Psi_n_ebe[ihar]->GetXaxis()->CenterTitle();
126 m_hist_Psi_n_ebe[ihar]->GetYaxis()->CenterTitle();
127
128 sprintf(name ,"hist_Psi%d_ebe_pt",ihar+1);
129 sprintf(name1,"%d#Delta#Psi (pT weighted);%d(#Psi_{reco}-#Psi_{Truth});events",ihar+1,ihar+1);
130 m_hist_Psi_n_ebe_pt [ihar]=new TH1D (name,name1,1000,-M_PI,M_PI);
131 if ( rootHistSvc->regHist(histPath+m_hist_Psi_n_ebe_pt[ihar]->GetName(),m_hist_Psi_n_ebe_pt[ihar]).isFailure() ){
132 msg(MSG::WARNING) << "Can't book "<< histPath+m_hist_Psi_n_ebe_pt[ihar]->GetName() << endmsg;
133 }
134 m_hist_Psi_n_ebe_pt[ihar]->GetXaxis()->CenterTitle();
135 m_hist_Psi_n_ebe_pt[ihar]->GetYaxis()->CenterTitle();
136
137
138
139
140
141
142 for(int ieta=0;ieta<n_etabin;ieta++){
143 sprintf(name ,"profile_pt_dep_%d_eta%d" ,ihar+1,ieta);
144 sprintf(name1,"v%d vs pT (eta%d);pT;v%d",ihar+1,ieta,ihar+1);
145 m_profile_pt_dep [ihar][ieta]=new TProfile (name,name1,n_ptbin,pt_binvals);
146 if ( rootHistSvc->regHist(histPath+m_profile_pt_dep[ihar][ieta]->GetName(),m_profile_pt_dep[ihar][ieta]).isFailure() ){
147 msg(MSG::WARNING) << "Can't book "<< histPath+m_profile_pt_dep[ihar][ieta]->GetName() << endmsg;
148 }
149 m_profile_pt_dep [ihar][ieta]->GetXaxis()->CenterTitle();
150 m_profile_pt_dep [ihar][ieta]->GetYaxis()->CenterTitle();
151 }
152
153 for(int ipt=0;ipt<n_ptbin;ipt++){
154 sprintf(name ,"profile_eta_dep_%d_pt%d",ihar+1,ipt);
155 sprintf(name1,"v%d vs #eta; (ipt%d)#eta;v%d",ihar+1,ipt,ihar+1);
156 m_profile_eta_dep [ihar][ipt]=new TProfile (name,name1,2*n_etabin, -eta_bin_max,eta_bin_max);
157 if ( rootHistSvc->regHist(histPath+m_profile_eta_dep[ihar][ipt]->GetName(),m_profile_eta_dep[ihar][ipt]).isFailure() ){
158 msg(MSG::WARNING) << "Can't book "<< histPath+m_profile_eta_dep[ihar][ipt]->GetName() << endmsg;
159 }
160 m_profile_eta_dep [ihar][ipt]->GetXaxis()->CenterTitle();
161 m_profile_eta_dep [ihar][ipt]->GetYaxis()->CenterTitle();
162 }
163
164
165 for(int ieta=0;ieta<n_etabin;ieta++){
166 sprintf(name ,"profile_pt_dep_reco_%d_eta%d",ihar+1,ieta);
167 sprintf(name1,"v%d vs pT (eta%d);pT;v%d",ihar+1,ieta,ihar+1);
168 m_profile_pt_dep_reco [ihar][ieta]=new TProfile (name,name1,n_ptbin,pt_binvals);
169 if ( rootHistSvc->regHist(histPath+m_profile_pt_dep_reco[ihar][ieta]->GetName(),m_profile_pt_dep_reco[ihar][ieta]).isFailure() ){
170 msg(MSG::WARNING) << "Can't book "<< histPath+m_profile_pt_dep_reco[ihar][ieta]->GetName() << endmsg;
171 }
172 m_profile_pt_dep_reco [ihar][ieta]->GetXaxis()->CenterTitle();
173 m_profile_pt_dep_reco [ihar][ieta]->GetYaxis()->CenterTitle();
174 }
175
176
177 for(int ipt=0;ipt<n_ptbin;ipt++){
178 sprintf(name ,"profile_eta_dep_reco_%d_pt%d",ihar+1,ipt);
179 sprintf(name1,"v%d vs #eta (pt%d);#eta;v%d",ihar+1,ipt,ihar+1);
180 m_profile_eta_dep_reco [ihar][ipt]=new TProfile (name,name1,2*n_etabin, -eta_bin_max,eta_bin_max);
181 if ( rootHistSvc->regHist(histPath+m_profile_eta_dep_reco[ihar][ipt]->GetName(),m_profile_eta_dep_reco[ihar][ipt]).isFailure() ){
182 msg(MSG::WARNING) << "Can't book "<< histPath+m_profile_eta_dep_reco[ihar][ipt]->GetName() << endmsg;
183 }
184 m_profile_eta_dep_reco [ihar][ipt]->GetXaxis()->CenterTitle();
185 m_profile_eta_dep_reco [ihar][ipt]->GetYaxis()->CenterTitle();
186 }
187
188
189 }
190 m_profile_resolution=new TProfile("profile_resolution","vn resolution;n;resolution",6, 0.5,6.5);
191 if(rootHistSvc->regHist(histPath+m_profile_resolution->GetName(),m_profile_resolution).isFailure() ){
192 msg(MSG::WARNING) << "Can't book "<< histPath+m_profile_resolution->GetName() << endmsg;
193 }
194
195 msg(MSG::DEBUG) << "Histograms have been booked " << endmsg;
197 return StatusCode::SUCCESS;
198}
199
200
201
203 msg(MSG::INFO) << ">>> CheckFlow_New from execute" << endmsg;
204
205 float pt_binvals[n_ptbin+1]={0.0,0.25,0.5,1.0,1.5,2.0,2.5,3.0,3.5,4.0,6.0,8.0,12.0,16.0,20.0,30.0,40.0};
206 float eta_bin_max = 4.0;
207
209 float b = hijing_pars->get_b();
210 float Psi_n[6],Psi_n_reco[6];
211 float Psi_n_reco_pos[6],Psi_n_reco_neg[6];
212 for(int ihar=0;ihar<6;ihar++){Psi_n[ihar]=hijing_pars->get_psi(ihar+1);}
213 msg(MSG::INFO)<<"SOUMYA "<<hijing_pars->get_psi(1)<<" "<<hijing_pars->get_psi(2)<<" "<<hijing_pars->get_psi(3)
214 <<hijing_pars->get_psi(4)<<" "<<hijing_pars->get_psi(5)<<" "<<hijing_pars->get_psi(6)<<" "<<b << endmsg;
215
216
217 // Check cut on impact parameter b
218 if(b<m_bcut_min || b>m_bcut_max) return StatusCode::SUCCESS;
219
220
221 double ngenerated_pos = 0,ngenerated_pt_pos=0;
222 double cos_n_pos[6],sin_n_pos[6],cos_n_pt_pos[6],sin_n_pt_pos[6];
223 double ngenerated_neg = 0,ngenerated_pt_neg=0;
224 double cos_n_neg[6],sin_n_neg[6],cos_n_pt_neg[6],sin_n_pt_neg[6];
225 for(int ihar=0;ihar<6;ihar++){cos_n_pos[ihar]=0;sin_n_pos[ihar]=0; cos_n_pt_pos[ihar]=0;sin_n_pt_pos[ihar]=0;
226 cos_n_neg[ihar]=0;sin_n_neg[ihar]=0; cos_n_pt_neg[ihar]=0;sin_n_pt_neg[ihar]=0;}
227
228 // Iterate over all MC particles
229 std::vector<HepMC::ConstGenParticlePtr> particles;
230 StatusCode stat = m_tesIO->getMC(particles, false, m_key);
231 if (stat.isFailure()) {
232 msg(MSG::ERROR) << "Could not find " << m_key << endmsg;
233 return stat;
234 }
235 for (auto pitr: particles) {
236 int pid = pitr->pdg_id();
237 double pt = pitr->momentum().perp();
238 double rapid = pitr->momentum().pseudoRapidity();
239 double phi = pitr->momentum().phi();
240 msg(MSG::DEBUG)
241 << " PID = " << pid << " Status = " << pitr->status()
242 << " Eta = " << rapid << " Phi = " << phi<< endmsg;
243
244 if( (std::abs(rapid) >= m_rapcut_min) && (std::abs(rapid) <= m_rapcut_max) &&
245 (std::abs(pt) >= m_ptcut_min) && (std::abs(pt) <= m_ptcut_max) ) {
246
247 for(int ihar=0;ihar<6;ihar++){
248 float temp=(ihar+1)*(phi-Psi_n[ihar]);
249
250 int ieta= (int)(std::abs(rapid)*n_etabin/eta_bin_max);
251 if(ieta>=0 && ieta<n_etabin) m_profile_pt_dep [ihar][ieta]->Fill(pt/1000,cos(temp));
252
253
254 float temp_pt=pt/1000;
255 for(int ipt=0;ipt<n_ptbin;ipt++){
256 if(temp_pt<pt_binvals[ipt+1]){
257 m_profile_eta_dep[ihar][ipt]->Fill(rapid ,cos(temp));
258 break;
259 }
260 }
261 //
262 const auto cosTerm = std::cos( (ihar+1)*phi);
263 const auto sinTerm = std::sin( (ihar+1)*phi);
264 //
265 if( rapid >3.2 && rapid< 4.9){
266 cos_n_pos[ihar]+=cosTerm;
267 sin_n_pos[ihar]+=sinTerm;
268 ngenerated_pos++;
269
270 cos_n_pt_pos[ihar]+=pt*cosTerm;
271 sin_n_pt_pos[ihar]+=pt*sinTerm;
272 ngenerated_pt_pos +=pt;
273 }
274 if( rapid <-3.2 && rapid >-4.9){
275 cos_n_neg[ihar]+=cosTerm;
276 sin_n_neg[ihar]+=sinTerm;
277 ngenerated_neg++;
278
279 cos_n_pt_neg[ihar]+=pt*cosTerm;
280 sin_n_pt_neg[ihar]+=pt*sinTerm;
281 ngenerated_pt_neg +=pt;
282 }
283
284 }
285 }
286 }
287
288
289// Calculate the event by event vn and also the reconstructed Psi_n angles
290// Also make correlation histos between Psi_n_truth and Psi_n_reco
291 float cos_n[6]{},sin_n[6]{},cos_n_pt[6]{},sin_n_pt[6]{};
292 const auto total = ngenerated_pos+ngenerated_neg;
293 const auto nTotalPt = ngenerated_pt_pos+ngenerated_pt_neg;
294 if ((total != 0) and (nTotalPt !=0)) [[likely]]{
295 for(int ihar=0;ihar<6;ihar++){
296 cos_n[ihar] = ( cos_n_pos[ihar]+ cos_n_neg[ihar] ) / total;
297 sin_n[ihar] = ( sin_n_pos[ihar]+ sin_n_neg[ihar] ) / total;
298
299 float psi_reco=std::atan2(sin_n[ihar],cos_n[ihar])/(ihar+1);
300 m_hist_Psi_n_ebe[ihar]->Fill( (ihar+1)*(psi_reco-Psi_n[ihar]) );
301 m_hist_vn_ebe [ihar]->Fill(std::sqrt(cos_n[ihar]*cos_n[ihar] +sin_n[ihar]*sin_n[ihar] ));
302
303 Psi_n_reco_pos[ihar]=std::atan2(sin_n_pos[ihar],cos_n_pos[ihar])/ (ihar+1);
304 Psi_n_reco_neg[ihar]=std::atan2(sin_n_neg[ihar],cos_n_neg[ihar])/ (ihar+1);
305 Psi_n_reco [ihar]=psi_reco;
306
307
308 cos_n_pt[ihar] = ( cos_n_pt_pos[ihar]+ cos_n_pt_neg[ihar] ) / nTotalPt;
309 sin_n_pt[ihar] = ( sin_n_pt_pos[ihar]+ sin_n_pt_neg[ihar] ) / nTotalPt;
310
311 psi_reco=std::atan2(sin_n_pt[ihar],cos_n_pt[ihar])/(ihar+1);
312 m_hist_Psi_n_ebe_pt[ihar]->Fill( (ihar+1)*(psi_reco-Psi_n[ihar]) );
313 }
314 }
315
316// Make the plots for the correlation between Psi_n truth (for different n) (same for Psi_n reco)
317 for(int ihar=0;ihar<6;ihar++){
318 m_hist_Psi_n_true[ihar]->Fill((ihar+1)*Psi_n[ihar]);
319 m_hist_Psi_n_reco[ihar]->Fill((ihar+1)*Psi_n_reco[ihar]);
320
321 float psi1,psi2;
322 for(int ihar2=0;ihar2<6;ihar2++){
323 psi1=(ihar+1)*Psi_n[ihar];psi2=(ihar2+1)*Psi_n[ihar2];
324 m_hist_psi_corr_true[ihar*6+ihar2]->Fill( std::atan2( std::sin(psi1-psi2),std::cos(psi1-psi2) ) );
325
326 psi1=(ihar+1)*Psi_n_reco[ihar];psi2=(ihar2+1)*Psi_n_reco[ihar2];
327 m_hist_psi_corr_reco[ihar*6+ihar2]->Fill( std::atan2( std::sin(psi1-psi2),std::cos(psi1-psi2) ) );
328 }
329 }
330
331
332
333
334
335// calculate the pt and eta dependence using the Psi_reco angles also fill the resolution TProfile
336 for(int ihar=0;ihar<6;ihar++){
337 m_profile_resolution->Fill( ihar+1, cos( (ihar+1) * (Psi_n_reco_pos[ihar] - Psi_n_reco_neg[ihar]) ) );
338 }
339 for (auto pitr: particles) {
340 double pt = pitr->momentum().perp();
341 double rapid = pitr->momentum().pseudoRapidity();
342 double phi = pitr->momentum().phi();
343 if( (std::abs(rapid) >= m_rapcut_min) && (std::abs(rapid) <= m_rapcut_max) &&
344 (std::abs(pt) >= m_ptcut_min) && (std::abs(pt) <= m_ptcut_max) ) {
345
346 for(int ihar=0;ihar<6;ihar++){
347 float temp=(ihar+1)*(phi-Psi_n_reco_pos[ihar]);
348 if(rapid>0) temp=(ihar+1)*(phi-Psi_n_reco_neg[ihar]);
349
350
351 int ieta= (int)(std::abs(rapid)*n_etabin/eta_bin_max);
352 if(ieta>=0 && ieta<n_etabin) m_profile_pt_dep_reco [ihar][ieta]->Fill(pt/1000,std::cos(temp));
353
354 float temp_pt=pt/1000;
355 for(int ipt=0;ipt<n_ptbin;ipt++){
356 if(temp_pt<pt_binvals[ipt+1]){
357 m_profile_eta_dep_reco[ihar][ipt]->Fill(rapid ,cos(temp));
358 break;
359 }
360 }
361 }
362 }
363 }
364
365
366
367
368
369 return StatusCode::SUCCESS;
370}
371
373 msg(MSG::INFO) << ">>> CheckFlow_New from finalize" << endmsg;
374
375 return StatusCode::SUCCESS;
376}
377
#define M_PI
Scalar phi() const
phi method
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
Base class from which all concrete Athena algorithm classes should be derived.
DoubleProperty m_ptcut_min
DoubleProperty m_ptcut_max
TH1D * m_hist_psi_corr_reco[36]
virtual StatusCode finalize() override
TH1D * m_hist_vn_ebe[6]
TProfile * m_profile_pt_dep_reco[6][n_etabin]
StringProperty m_key
SG::ReadHandleKey< HijingEventParams > m_hijingKey
TH1D * m_hist_Psi_n_true[6]
virtual StatusCode execute() override
TProfile * m_profile_eta_dep[6][n_ptbin]
TProfile * m_profile_resolution
DoubleProperty m_rapcut_min
TH1D * m_hist_Psi_n_ebe_pt[6]
virtual StatusCode initialize() override
TProfile * m_profile_eta_dep_reco[6][n_ptbin]
DoubleProperty m_bcut_max
TH1D * m_hist_Psi_n_reco[6]
TH1D * m_hist_Psi_n_ebe[6]
TruthHelper::GenAccessIO * m_tesIO
TProfile * m_profile_pt_dep[6][n_etabin]
DoubleProperty m_rapcut_max
TH1D * m_hist_psi_corr_true[36]
#define likely(x)
Definition dictionary.h:40
MsgStream & msg
Definition testRead.cxx:32