EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
QAG4SimulationEicCalorimeter.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file QAG4SimulationEicCalorimeter.cc
2 
3 #include <qa_modules/QAHistManagerDef.h>
4 
5 #include <g4main/PHG4Hit.h>
7 #include <g4main/PHG4Particle.h>
9 #include <g4main/PHG4VtxPoint.h>
10 
11 #include <calobase/RawCluster.h>
12 #include <calobase/RawClusterContainer.h>
13 #include <calobase/RawTower.h>
14 #include <calobase/RawTowerContainer.h>
15 #include <calobase/RawTowerGeomContainer.h>
16 
17 #include <g4eval/CaloEvalStack.h>
19 
22 #include <fun4all/SubsysReco.h>
23 
24 #include <phool/PHCompositeNode.h>
25 #include <phool/PHNodeIterator.h>
26 #include <phool/getClass.h>
27 #include <phool/phool.h>
28 
29 #include <TAxis.h>
30 #include <TH1.h>
31 #include <TH2.h>
32 #include <TNamed.h>
33 #include <TString.h>
34 #include <TVector3.h>
35 
36 #include <CLHEP/Vector/ThreeVector.h> // for Hep3Vector
37 
38 #include <cassert>
39 #include <cstdlib>
40 #include <iostream>
41 #include <iterator> // for reverse_iterator
42 #include <map>
43 #include <utility>
44 
45 using namespace std;
46 
49  : SubsysReco("QAG4SimulationEicCalorimeter_" + calo_name)
50  , _calo_name(calo_name)
51  , _flags(flags)
52  , _calo_hit_container(nullptr)
53  , _calo_abs_hit_container(nullptr)
54  , _truth_container(nullptr)
55 {
56 }
57 
59 {
60  PHNodeIterator iter(topNode);
61  PHCompositeNode *dstNode = static_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
62  assert(dstNode);
63 
64  if (flag(kProcessG4Hit))
65  {
66  _calo_hit_container = findNode::getClass<PHG4HitContainer>(topNode,
67  "G4HIT_" + _calo_name);
69  {
70  cout << "QAG4SimulationEicCalorimeter::InitRun - Fatal Error - "
71  << "unable to find DST node "
72  << "G4HIT_" + _calo_name << endl;
73  assert(_calo_hit_container);
74  }
75 
76  _calo_abs_hit_container = findNode::getClass<PHG4HitContainer>(topNode,
77  "G4HIT_ABSORBER_" + _calo_name);
79  {
80  cout << "QAG4SimulationEicCalorimeter::InitRun - Fatal Error - "
81  << "unable to find DST node "
82  << "G4HIT_ABSORBER_" + _calo_name
83  << endl;
85  }
86  }
87 
88  _truth_container = findNode::getClass<PHG4TruthInfoContainer>(topNode,
89  "G4TruthInfo");
90  if (!_truth_container)
91  {
92  cout << "QAG4SimulationEicCalorimeter::InitRun - Fatal Error - "
93  << "unable to find DST node "
94  << "G4TruthInfo" << endl;
95  assert(_truth_container);
96  }
97 
98  if (flag(kProcessCluster))
99  {
100  if (!_caloevalstack)
101  {
102  _caloevalstack.reset(new CaloEvalStack(topNode, _calo_name));
103  _caloevalstack->set_strict(true);
104  _caloevalstack->set_verbosity(Verbosity() + 1);
105  }
106  else
107  {
108  _caloevalstack->next_event(topNode);
109  }
110  }
112 }
113 
115 {
117  assert(hm);
118  TH1D *h = new TH1D(TString(get_histo_prefix()) + "_Normalization", //
119  TString(_calo_name) + " Normalization;Items;Count", 10, .5, 10.5);
120  int i = 1;
121  h->GetXaxis()->SetBinLabel(i++, "Event");
122  h->GetXaxis()->SetBinLabel(i++, "G4Hit Active");
123  h->GetXaxis()->SetBinLabel(i++, "G4Hit Absor.");
124  h->GetXaxis()->SetBinLabel(i++, "Tower");
125  h->GetXaxis()->SetBinLabel(i++, "Tower Hit");
126  h->GetXaxis()->SetBinLabel(i++, "Cluster");
127  h->GetXaxis()->LabelsOption("v");
128  hm->registerHisto(h);
129 
130  if (flag(kProcessG4Hit))
131  {
132  if (Verbosity() >= 1)
133  cout << "QAG4SimulationEicCalorimeter::Init - Process sampling fraction"
134  << endl;
135  Init_G4Hit(topNode);
136  }
137  if (flag(kProcessTower))
138  {
139  if (Verbosity() >= 1)
140  cout << "QAG4SimulationEicCalorimeter::Init - Process tower occupancies"
141  << endl;
142  Init_Tower(topNode);
143  }
144  if (flag(kProcessCluster))
145  {
146  if (Verbosity() >= 1)
147  cout << "QAG4SimulationEicCalorimeter::Init - Process tower occupancies"
148  << endl;
149  Init_Cluster(topNode);
150  }
151 
153 }
154 
156 {
157  if (Verbosity() > 2)
158  cout << "QAG4SimulationEicCalorimeter::process_event() entered" << endl;
159 
160  if (_caloevalstack)
161  _caloevalstack->next_event(topNode);
162 
163  if (flag(kProcessG4Hit))
164  {
165  int ret = process_event_G4Hit(topNode);
166 
167  if (ret != Fun4AllReturnCodes::EVENT_OK)
168  return ret;
169  }
170 
171  if (flag(kProcessTower))
172  {
173  int ret = process_event_Tower(topNode);
174 
175  if (ret != Fun4AllReturnCodes::EVENT_OK)
176  return ret;
177  }
178 
179  if (flag(kProcessCluster))
180  {
181  int ret = process_event_Cluster(topNode);
182 
183  if (ret != Fun4AllReturnCodes::EVENT_OK)
184  return ret;
185  }
186 
187  // at the end, count success events
189  assert(hm);
190  TH1D *h_norm = dynamic_cast<TH1D *>(hm->getHisto(
191  get_histo_prefix() + "_Normalization"));
192  assert(h_norm);
193  h_norm->Fill("Event", 1);
194 
196 }
197 
198 string
200 {
201  return "h_QAG4Sim_" + string(_calo_name);
202 }
203 
205 {
207  assert(hm);
208 
209  hm->registerHisto(
210  new TH2F(TString(get_histo_prefix()) + "_G4Hit_RZ", //
211  TString(_calo_name) + " RZ projection;G4 Hit Z (cm);G4 Hit R (cm)", 1200, -300, 300,
212  600, -000, 300));
213 
214  hm->registerHisto(
215  new TH2F(TString(get_histo_prefix()) + "_G4Hit_XY", //
216  TString(_calo_name) + " XY projection;G4 Hit X (cm);G4 Hit Y (cm)", 1200, -300, 300,
217  1200, -300, 300));
218 
219  hm->registerHisto(
220  new TH2F(TString(get_histo_prefix()) + "_G4Hit_LateralTruthProjection", //
221  TString(_calo_name) + " shower lateral projection (last primary);Polar direction (cm);Azimuthal direction (cm)",
222  200, -30, 30, 200, -30, 30));
223 
224  hm->registerHisto(new TH1F(TString(get_histo_prefix()) + "_G4Hit_SF", //
225  TString(_calo_name) + " sampling fraction;Sampling fraction", 1000, 0, .2));
226 
227  hm->registerHisto(
228  new TH1F(TString(get_histo_prefix()) + "_G4Hit_VSF", //
229  TString(_calo_name) + " visible sampling fraction;Visible sampling fraction", 1000, 0,
230  .2));
231 
232  TH1F *h =
233  new TH1F(TString(get_histo_prefix()) + "_G4Hit_HitTime", //
234  TString(_calo_name) + " hit time (edep weighting);Hit time - T0 (ns);Geant4 energy density",
235  1000, 0.5, 10000);
236  QAHistManagerDef::useLogBins(h->GetXaxis());
237  hm->registerHisto(h);
238 
239  hm->registerHisto(
240  new TH1F(TString(get_histo_prefix()) + "_G4Hit_FractionTruthEnergy", //
241  TString(_calo_name) + " fraction truth energy ;G4 edep / particle energy",
242  1000, 0, 1));
243 
244  hm->registerHisto(
245  new TH1F(TString(get_histo_prefix()) + "_G4Hit_FractionEMVisibleEnergy", //
246  TString(_calo_name) + " fraction visible energy from EM; visible energy from e^{#pm} / total visible energy",
247  100, 0, 1));
248 
250 }
251 
253 {
254  if (Verbosity() > 2)
255  cout << "QAG4SimulationEicCalorimeter::process_event_G4Hit() entered" << endl;
256 
257  TH1F *h = nullptr;
258 
260  assert(hm);
261 
262  TH1D *h_norm = dynamic_cast<TH1D *>(hm->getHisto(
263  get_histo_prefix() + "_Normalization"));
264  assert(h_norm);
265 
266  // get primary
267  assert(_truth_container);
268  PHG4TruthInfoContainer::ConstRange primary_range =
270  double total_primary_energy = 1e-9; //make it zero energy epsilon samll so it can be used for denominator
271  for (PHG4TruthInfoContainer::ConstIterator particle_iter = primary_range.first;
272  particle_iter != primary_range.second; ++particle_iter)
273  {
274  const PHG4Particle *particle = particle_iter->second;
275  assert(particle);
276  total_primary_energy += particle->get_e();
277  }
278 
279  assert(not _truth_container->GetMap().empty());
280  const PHG4Particle *last_primary =
281  _truth_container->GetMap().rbegin()->second;
282  assert(last_primary);
283 
284  if (Verbosity() > 2)
285  {
286  cout
287  << "QAG4SimulationEicCalorimeter::process_event_G4Hit() handle this truth particle"
288  << endl;
289  last_primary->identify();
290  }
291  const PHG4VtxPoint *primary_vtx = //
292  _truth_container->GetPrimaryVtx(last_primary->get_vtx_id());
293  assert(primary_vtx);
294  if (Verbosity() > 2)
295  {
296  cout
297  << "QAG4SimulationEicCalorimeter::process_event_G4Hit() handle this vertex"
298  << endl;
299  primary_vtx->identify();
300  }
301 
302  const double t0 = primary_vtx->get_t();
303  const TVector3 vertex(primary_vtx->get_x(), primary_vtx->get_y(),
304  primary_vtx->get_z());
305 
306  // projection axis
307  TVector3 axis_proj(last_primary->get_px(), last_primary->get_py(),
308  last_primary->get_pz());
309  if (axis_proj.Mag() == 0)
310  axis_proj.SetXYZ(0, 0, 1);
311  axis_proj = axis_proj.Unit();
312 
313  // azimuthal direction axis
314  TVector3 axis_azimuth = axis_proj.Cross(TVector3(0, 0, 1));
315  if (axis_azimuth.Mag() == 0)
316  axis_azimuth.SetXYZ(1, 0, 0);
317  axis_azimuth = axis_azimuth.Unit();
318 
319  // polar direction axis
320  TVector3 axis_polar = axis_proj.Cross(axis_azimuth);
321  assert(axis_polar.Mag() > 0);
322  axis_polar = axis_polar.Unit();
323 
324  double e_calo = 0.0; // active energy deposition
325  double ev_calo = 0.0; // visible energy
326  double ea_calo = 0.0; // absorber energy
327  double ev_calo_em = 0.0; // EM visible energy
328 
330  {
331  TH2F *hrz = dynamic_cast<TH2F *>(hm->getHisto(
332  get_histo_prefix() + "_G4Hit_RZ"));
333  assert(hrz);
334  TH2F *hxy = dynamic_cast<TH2F *>(hm->getHisto(
335  get_histo_prefix() + "_G4Hit_XY"));
336  assert(hxy);
337  TH1F *ht = dynamic_cast<TH1F *>(hm->getHisto(
338  get_histo_prefix() + "_G4Hit_HitTime"));
339  assert(ht);
340  TH2F *hlat = dynamic_cast<TH2F *>(hm->getHisto(
341  get_histo_prefix() + "_G4Hit_LateralTruthProjection"));
342  assert(hlat);
343 
344  h_norm->Fill("G4Hit Active", _calo_hit_container->size());
345  PHG4HitContainer::ConstRange calo_hit_range =
347  for (PHG4HitContainer::ConstIterator hit_iter = calo_hit_range.first;
348  hit_iter != calo_hit_range.second; hit_iter++)
349  {
350  PHG4Hit *this_hit = hit_iter->second;
351  assert(this_hit);
352 
353  e_calo += this_hit->get_edep();
354  ev_calo += this_hit->get_light_yield();
355 
356  // EM visible energy that is only associated with electron energy deposition
358  this_hit->get_trkid());
359  if (!particle)
360  {
361  cout << __PRETTY_FUNCTION__ << " - Error - this PHG4hit missing particle: ";
362  this_hit->identify();
363  }
364  assert(particle);
365  if (abs(particle->get_pid()) == 11)
366  ev_calo_em += this_hit->get_light_yield();
367 
368  const TVector3 hit(this_hit->get_avg_x(), this_hit->get_avg_y(),
369  this_hit->get_avg_z());
370 
371  hrz->Fill(hit.Z(), hit.Perp(), this_hit->get_edep());
372  hxy->Fill(hit.X(), hit.Y(), this_hit->get_edep());
373  ht->Fill(this_hit->get_avg_t() - t0, this_hit->get_edep());
374 
375  const double hit_azimuth = axis_azimuth.Dot(hit - vertex);
376  const double hit_polar = axis_polar.Dot(hit - vertex);
377  hlat->Fill(hit_polar, hit_azimuth, this_hit->get_edep());
378  }
379  }
380 
382  {
383  h_norm->Fill("G4Hit Absor.", _calo_abs_hit_container->size());
384 
385  PHG4HitContainer::ConstRange calo_abs_hit_range =
387  for (PHG4HitContainer::ConstIterator hit_iter = calo_abs_hit_range.first;
388  hit_iter != calo_abs_hit_range.second; hit_iter++)
389  {
390  PHG4Hit *this_hit = hit_iter->second;
391  assert(this_hit);
392 
393  ea_calo += this_hit->get_edep();
394  }
395  }
396 
397  if (Verbosity() > 3)
398  cout << "QAG4SimulationEicCalorimeter::process_event_G4Hit::" << _calo_name
399  << " - SF = " << e_calo / (e_calo + ea_calo + 1e-9) << ", VSF = "
400  << ev_calo / (e_calo + ea_calo + 1e-9) << endl;
401 
402  if (e_calo + ea_calo > 0)
403  {
404  h = dynamic_cast<TH1F *>(hm->getHisto(get_histo_prefix() + "_G4Hit_SF"));
405  assert(h);
406  h->Fill(e_calo / (e_calo + ea_calo));
407 
408  h = dynamic_cast<TH1F *>(hm->getHisto(get_histo_prefix() + "_G4Hit_VSF"));
409  assert(h);
410  h->Fill(ev_calo / (e_calo + ea_calo));
411  }
412 
413  h = dynamic_cast<TH1F *>(hm->getHisto(
414  get_histo_prefix() + "_G4Hit_FractionTruthEnergy"));
415  assert(h);
416  h->Fill((e_calo + ea_calo) / total_primary_energy);
417 
418  if (ev_calo > 0)
419  {
420  h = dynamic_cast<TH1F *>(hm->getHisto(
421  get_histo_prefix() + "_G4Hit_FractionEMVisibleEnergy"));
422  assert(h);
423  h->Fill(ev_calo_em / (ev_calo));
424  }
425 
426  if (Verbosity() > 3)
427  cout << "QAG4SimulationEicCalorimeter::process_event_G4Hit::" << _calo_name
428  << " - histogram " << h->GetName() << " Get Sum = " << h->GetSum()
429  << endl;
430 
432 }
433 
435 {
437  assert(hm);
438 
439  TH1F *h = new TH1F(TString(get_histo_prefix()) + "_Tower_1x1", //
440  TString(_calo_name) + " 1x1 tower;1x1 TOWER Energy (GeV)", 100, 9e-4, 100);
441  QAHistManagerDef::useLogBins(h->GetXaxis());
442  hm->registerHisto(h);
443 
444  hm->registerHisto(
445  new TH1F(TString(get_histo_prefix()) + "_Tower_1x1_max", //
446  TString(_calo_name) + " 1x1 tower max per event;1x1 tower max per event (GeV)", 5000,
447  0, 50));
448 
449  h = new TH1F(TString(get_histo_prefix()) + "_Tower_2x2", //
450  TString(_calo_name) + " 2x2 tower;2x2 TOWER Energy (GeV)", 100, 9e-4, 100);
451  QAHistManagerDef::useLogBins(h->GetXaxis());
452  hm->registerHisto(h);
453  hm->registerHisto(
454  new TH1F(TString(get_histo_prefix()) + "_Tower_2x2_max", //
455  TString(_calo_name) + " 2x2 tower max per event;2x2 tower max per event (GeV)", 5000,
456  0, 50));
457 
458  h = new TH1F(TString(get_histo_prefix()) + "_Tower_3x3", //
459  TString(_calo_name) + " 3x3 tower;3x3 TOWER Energy (GeV)", 100, 9e-4, 100);
460  QAHistManagerDef::useLogBins(h->GetXaxis());
461  hm->registerHisto(h);
462  hm->registerHisto(
463  new TH1F(TString(get_histo_prefix()) + "_Tower_3x3_max", //
464  TString(_calo_name) + " 3x3 tower max per event;3x3 tower max per event (GeV)", 5000,
465  0, 50));
466 
467  h = new TH1F(TString(get_histo_prefix()) + "_Tower_4x4", //
468  TString(_calo_name) + " 4x4 tower;4x4 TOWER Energy (GeV)", 100, 9e-4, 100);
469  QAHistManagerDef::useLogBins(h->GetXaxis());
470  hm->registerHisto(h);
471  hm->registerHisto(
472  new TH1F(TString(get_histo_prefix()) + "_Tower_4x4_max", //
473  TString(_calo_name) + " 4x4 tower max per event;4x4 tower max per event (GeV)", 5000,
474  0, 50));
475 
476  h = new TH1F(TString(get_histo_prefix()) + "_Tower_5x5", //
477  TString(_calo_name) + " 5x5 tower;5x5 TOWER Energy (GeV)", 100, 9e-4, 100);
478  QAHistManagerDef::useLogBins(h->GetXaxis());
479  hm->registerHisto(h);
480  hm->registerHisto(
481  new TH1F(TString(get_histo_prefix()) + "_Tower_5x5_max", //
482  TString(_calo_name) + " 5x5 tower max per event;5x5 tower max per event (GeV)", 5000,
483  0, 50));
484 
486 }
487 
489 {
490  const string detector(_calo_name);
491 
492  if (Verbosity() > 2)
493  cout << "QAG4SimulationEicCalorimeter::process_event_Tower() entered" << endl;
494 
496  assert(hm);
497  TH1D *h_norm = dynamic_cast<TH1D *>(hm->getHisto(
498  get_histo_prefix() + "_Normalization"));
499  assert(h_norm);
500 
501  string towernodename = "TOWER_CALIB_" + detector;
502  // Grab the towers
503  RawTowerContainer *towers = findNode::getClass<RawTowerContainer>(topNode,
504  towernodename.c_str());
505  if (!towers)
506  {
507  std::cout << PHWHERE << ": Could not find node " << towernodename.c_str()
508  << std::endl;
510  }
511  string towergeomnodename = "TOWERGEOM_" + detector;
512  RawTowerGeomContainer *towergeom = findNode::getClass<RawTowerGeomContainer>(
513  topNode, towergeomnodename.c_str());
514  if (!towergeom)
515  {
516  cout << PHWHERE << ": Could not find node " << towergeomnodename.c_str()
517  << endl;
519  }
520 
521  static const int max_size = 5;
522  map<int, string> size_label;
523  size_label[1] = "1x1";
524  size_label[2] = "2x2";
525  size_label[3] = "3x3";
526  size_label[4] = "4x4";
527  size_label[5] = "5x5";
528  map<int, double> max_energy;
529  map<int, TH1F *> energy_hist_list;
530  map<int, TH1F *> max_energy_hist_list;
531 
532  for (int size = 1; size <= max_size; ++size)
533  {
534  max_energy[size] = 0;
535 
536  TH1F *h = dynamic_cast<TH1F *>(hm->getHisto(
537  get_histo_prefix() + "_Tower_" + size_label[size]));
538  assert(h);
539  energy_hist_list[size] = h;
540  h = dynamic_cast<TH1F *>(hm->getHisto(
541  get_histo_prefix() + "_Tower_" + size_label[size] + "_max"));
542  assert(h);
543  max_energy_hist_list[size] = h;
544  }
545 
546  h_norm->Fill("Tower", towergeom->size()); // total tower count
547  h_norm->Fill("Tower Hit", towers->size());
548 
549  for (int binphi = 0; binphi < towergeom->get_phibins(); ++binphi)
550  {
551  for (int bineta = 0; bineta < towergeom->get_etabins(); ++bineta)
552  {
553  for (int size = 1; size <= max_size; ++size)
554  {
555  // for 2x2 and 4x4 use slide-2 window as implemented in DAQ
556  if ((size == 2 or size == 4) and ((binphi % 2 != 0) and (bineta % 2 != 0)))
557  continue;
558 
559  double energy = 0;
560 
561  // sliding window made from 2x2 sums
562  for (int iphi = binphi; iphi < binphi + size; ++iphi)
563  {
564  for (int ieta = bineta; ieta < bineta + size; ++ieta)
565  {
566  if (ieta > towergeom->get_etabins())
567  continue;
568 
569  // wrap around
570  int wrapphi = iphi;
571  assert(wrapphi >= 0);
572  if (wrapphi >= towergeom->get_phibins())
573  {
574  wrapphi = wrapphi - towergeom->get_phibins();
575  }
576 
577  RawTower *tower = towers->getTower(ieta, wrapphi);
578 
579  if (tower)
580  {
581  const double e_intput = tower->get_energy();
582 
583  energy += e_intput;
584  }
585  }
586  }
587 
588  energy_hist_list[size]->Fill(energy == 0 ? 9.1e-4 : energy); // trick to fill 0 energy tower to the first bin
589 
590  if (energy > max_energy[size])
591  max_energy[size] = energy;
592 
593  } // for (int size = 1; size <= 4; ++size)
594  }
595  }
596 
597  for (int size = 1; size <= max_size; ++size)
598  {
599  max_energy_hist_list[size]->Fill(max_energy[size]);
600  }
602 }
603 
605 {
607  assert(hm);
608 
609  hm->registerHisto(
610  new TH1F(TString(get_histo_prefix()) + "_Cluster_BestMatchERatio", //
611  TString(_calo_name) + " best matched cluster E/E_{Truth};E_{Cluster}/E_{Truth}", 150,
612  0, 1.5));
613 
614  hm->registerHisto(
615  new TH2F(TString(get_histo_prefix()) + "_Cluster_LateralTruthProjection", //
616  TString(_calo_name) + " best cluster lateral projection (last primary);Polar direction (cm);Azimuthal direction (cm)",
617  200, -15, 15, 200, -15, 15));
618 
620 }
621 
623 {
624  if (Verbosity() > 2)
625  cout << "QAG4SimulationEicCalorimeter::process_event_Cluster() entered"
626  << endl;
627 
629  assert(hm);
630  string towergeomnodename = "TOWERGEOM_" + _calo_name;
631  RawTowerGeomContainer *towergeom = findNode::getClass<RawTowerGeomContainer>(
632  topNode, towergeomnodename.c_str());
633  if (!towergeom)
634  {
635  cout << PHWHERE << ": Could not find node " << towergeomnodename.c_str()
636  << endl;
638  }
639 
640  //get a cluster count
641  TH1D *h_norm = dynamic_cast<TH1D *>(hm->getHisto(
642  get_histo_prefix() + "_Normalization"));
643  assert(h_norm);
644 
645  std::string nodename = "CLUSTER_" + _calo_name;
646  RawClusterContainer *clusters = findNode::getClass<RawClusterContainer>(
647  topNode, nodename.c_str());
648  assert(clusters);
649  h_norm->Fill("Cluster", clusters->size());
650 
651  // get primary
652  assert(_truth_container);
653  assert(not _truth_container->GetMap().empty());
654  PHG4Particle *last_primary = _truth_container->GetMap().rbegin()->second;
655  assert(last_primary);
656 
657  if (Verbosity() > 2)
658  {
659  cout
660  << "QAG4SimulationEicCalorimeter::process_event_Cluster() handle this truth particle"
661  << endl;
662  last_primary->identify();
663  }
664 
665  assert(_caloevalstack);
666  CaloRawClusterEval *clustereval = _caloevalstack->get_rawcluster_eval();
667  assert(clustereval);
668 
669  TH1F *h = dynamic_cast<TH1F *>(hm->getHisto(
670  get_histo_prefix() + "_Cluster_BestMatchERatio"));
671  assert(h);
672 
673  RawCluster *cluster = clustereval->best_cluster_from(last_primary);
674  if (cluster)
675  {
676  // has a cluster matched and best cluster selected
677 
678  if (Verbosity() > 3)
679  cout << "QAG4SimulationEicCalorimeter::process_event_Cluster::"
680  << _calo_name << " - get cluster with energy "
681  << cluster->get_energy() << " VS primary energy "
682  << last_primary->get_e() << endl;
683 
684  h->Fill(cluster->get_energy() / (last_primary->get_e() + 1e-9)); //avoids divide zero
685 
686  // now work on the projection:
687  const CLHEP::Hep3Vector hit(cluster->get_position());
688 
689  const PHG4VtxPoint *primary_vtx = //
690  _truth_container->GetPrimaryVtx(last_primary->get_vtx_id());
691  assert(primary_vtx);
692  if (Verbosity() > 2)
693  {
694  cout
695  << "QAG4SimulationEicCalorimeter::process_event_Cluster() handle this vertex"
696  << endl;
697  primary_vtx->identify();
698  }
699 
700  const CLHEP::Hep3Vector vertex(primary_vtx->get_x(), primary_vtx->get_y(),
701  primary_vtx->get_z());
702 
703  // projection axis
704  CLHEP::Hep3Vector axis_proj(last_primary->get_px(), last_primary->get_py(),
705  last_primary->get_pz());
706  if (axis_proj.mag() == 0)
707  axis_proj.set(0, 0, 1);
708  axis_proj = axis_proj.unit();
709 
710  // azimuthal direction axis
711  CLHEP::Hep3Vector axis_azimuth = axis_proj.cross(CLHEP::Hep3Vector(0, 0, 1));
712  if (axis_azimuth.mag() == 0)
713  axis_azimuth.set(1, 0, 0);
714  axis_azimuth = axis_azimuth.unit();
715 
716  // polar direction axis
717  CLHEP::Hep3Vector axis_polar = axis_proj.cross(axis_azimuth);
718  assert(axis_polar.mag() > 0);
719  axis_polar = axis_polar.unit();
720 
721  TH2F *hlat = dynamic_cast<TH2F *>(hm->getHisto(
722  get_histo_prefix() + "_Cluster_LateralTruthProjection"));
723  assert(hlat);
724 
725  const double hit_azimuth = axis_azimuth.dot(hit - vertex);
726  const double hit_polar = axis_polar.dot(hit - vertex);
727  hlat->Fill(hit_polar, hit_azimuth);
728  }
729  else
730  {
731  if (Verbosity() > 3)
732  cout << "QAG4SimulationEicCalorimeter::process_event_Cluster::"
733  << _calo_name << " - missing cluster !";
734  h->Fill(0); // no cluster matched
735  }
736 
738 }