EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PHG4FullProjSpacalCellReco.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file PHG4FullProjSpacalCellReco.cc
2 
3 #include "PHG4Cell.h" // for PHG4Cell
4 #include "PHG4CylinderGeom.h" // for PHG4CylinderGeom
5 #include "PHG4CylinderGeom_Spacalv1.h" // for PHG4CylinderGeom_Spaca...
9 #include "PHG4CylinderCellGeom.h"
11 
12 #include "PHG4CellContainer.h"
13 #include "PHG4CellDefs.h"
14 #include "PHG4Cellv1.h"
15 
16 #include <phparameter/PHParameterInterface.h> // for PHParameterInterface
17 
18 #include <g4main/PHG4Hit.h>
20 
21 #include <fun4all/Fun4AllBase.h> // for Fun4AllBase::VERBOSITY...
23 #include <fun4all/Fun4AllServer.h>
24 #include <fun4all/SubsysReco.h> // for SubsysReco
25 
26 #include <phool/PHCompositeNode.h>
27 #include <phool/PHIODataNode.h>
28 #include <phool/PHNode.h> // for PHNode
29 #include <phool/PHNodeIterator.h>
30 #include <phool/PHObject.h> // for PHObject
31 #include <phool/getClass.h>
32 #include <phool/phool.h> // for PHWHERE
33 
34 #include <TAxis.h> // for TAxis
35 #include <TFile.h>
36 #include <TH1.h>
37 #include <TH2.h>
38 #include <TObject.h> // for TObject
39 
40 #include <boost/foreach.hpp>
41 
42 #include <cassert>
43 #include <cmath>
44 #include <cstdlib>
45 #include <exception>
46 #include <iostream>
47 #include <sstream>
48 #include <utility> // for pair
49 
50 
51 using namespace std;
52 
54  SubsysReco(name),
56  sum_energy_g4hit(0),
57  chkenergyconservation(0),
58  tmin(NAN),
59  tmax(NAN),
60  light_collection_model()
61 {
63 }
64 
65 int
67 {
68  sum_energy_g4hit = 0.;
70 }
71 
72 int
74 {
75  PHNodeIterator iter(topNode);
76 
77  // Looking for the DST node
78  PHCompositeNode *dstNode;
79  dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
80  if (!dstNode)
81  {
82  std::cout << PHWHERE << "DST Node missing, doing nothing." << std::endl;
83  exit(1);
84  }
85  hitnodename = "G4HIT_" + detector;
86  PHG4HitContainer *g4hit = findNode::getClass<PHG4HitContainer>(topNode, hitnodename.c_str());
87  if (!g4hit)
88  {
89  cout << "PHG4FullProjSpacalCellReco::InitRun - Could not locate g4 hit node "
90  << hitnodename << endl;
91  topNode->print();
92  exit(1);
93  }
94  cellnodename = "G4CELL_" + detector;
95  PHG4CellContainer *cells = findNode::getClass<PHG4CellContainer>(topNode, cellnodename);
96  if (!cells)
97  {
98  PHNodeIterator dstiter(dstNode);
99  PHCompositeNode *DetNode = dynamic_cast<PHCompositeNode*>(dstiter.findFirst("PHCompositeNode", detector));
100  if (!DetNode)
101  {
102  DetNode = new PHCompositeNode(detector);
103  dstNode->addNode(DetNode);
104  }
105  cells = new PHG4CellContainer();
106  PHIODataNode<PHObject> *newNode = new PHIODataNode<PHObject>(cells, cellnodename.c_str(), "PHObject");
107  DetNode->addNode(newNode);
108  }
109 
110  geonodename = "CYLINDERGEOM_" + detector;
112  findNode::getClass<PHG4CylinderGeomContainer>(topNode, geonodename.c_str());
113  if (!geo)
114  {
115  cout << "PHG4FullProjSpacalCellReco::InitRun - Could not locate geometry node "
116  << geonodename << endl;
117  topNode->print();
118  exit(1);
119  }
120  if (Verbosity() > 0)
121  {
122  cout << "PHG4FullProjSpacalCellReco::InitRun - incoming geometry:"
123  << endl;
124  geo->identify();
125  }
126  seggeonodename = "CYLINDERCELLGEOM_" + detector;
127  PHG4CylinderCellGeomContainer *seggeo = findNode::getClass<PHG4CylinderCellGeomContainer>(topNode, seggeonodename.c_str());
128  if (!seggeo)
129  {
130  seggeo = new PHG4CylinderCellGeomContainer();
131  PHCompositeNode *runNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "RUN"));
132  PHIODataNode<PHObject> *newNode = new PHIODataNode<PHObject>(seggeo, seggeonodename.c_str(), "PHObject");
133  runNode->addNode(newNode);
134  }
135 
136  const PHG4CylinderGeom *layergeom_raw = geo->GetFirstLayerGeom();
137  assert(layergeom_raw);
138 
139  // a special implimentation of PHG4CylinderGeom is required here.
140  const PHG4CylinderGeom_Spacalv3 *layergeom =
141  dynamic_cast<const PHG4CylinderGeom_Spacalv3 *>(layergeom_raw);
142 
143  if (!layergeom)
144  {
145  cout << "PHG4FullProjSpacalCellReco::InitRun - Fatal Error -"
146  << " require to work with a version of SPACAL geometry of PHG4CylinderGeom_Spacalv3 or higher. "
147  << "However the incoming geometry has version "
148  << layergeom_raw->ClassName() << endl;
149  exit(1);
150  }
151  if (Verbosity() > 1)
152  {
153  layergeom->identify();
154  }
155 
156  layergeom->subtower_consistency_check();
157 
158  // int layer = layergeom->get_layer();
159 
160  // create geo object and fill with variables common to all binning methods
161  PHG4CylinderCellGeom_Spacalv1 *layerseggeo =
163  layerseggeo->set_layer(layergeom->get_layer());
164  layerseggeo->set_radius(layergeom->get_radius());
165  layerseggeo->set_thickness(layergeom->get_thickness());
167 
168  // construct a map to convert tower_ID into the older eta bins.
169 
170  const PHG4CylinderGeom_Spacalv3::tower_map_t & tower_map =
171  layergeom->get_sector_tower_map();
172  const PHG4CylinderGeom_Spacalv3::sector_map_t & sector_map =
173  layergeom->get_sector_map();
174  const int nphibin = layergeom->get_azimuthal_n_sec() // sector
175  * layergeom->get_max_phi_bin_in_sec() // blocks per sector
176  * layergeom->get_n_subtower_phi(); // subtower per block
177  const double deltaphi = 2. * M_PI / nphibin;
178 
179  typedef map<double, int> map_z_tower_z_ID_t;
180  map_z_tower_z_ID_t map_z_tower_z_ID;
181  double phi_min = NAN;
182 
183  BOOST_FOREACH(const PHG4CylinderGeom_Spacalv3::tower_map_t::value_type& tower_pair, tower_map)
184  {
185 
186  const int & tower_ID = tower_pair.first;
188  tower_pair.second;
189 
190  // inspect index in sector 0
191  pair<int, int> tower_z_phi_ID = layergeom->get_tower_z_phi_ID(tower_ID, 0);
192 
193  const int & tower_ID_z = tower_z_phi_ID.first;
194  const int & tower_ID_phi = tower_z_phi_ID.second;
195 
196  if (tower_ID_phi == 0)
197  {
198  //assign phi min according phi bin 0
199  phi_min = M_PI_2 - deltaphi *(layergeom->get_max_phi_bin_in_sec()* layergeom->get_n_subtower_phi()/2) // shift of first tower in sector
200  + sector_map.begin()->second;
201  }
202 
203  if (tower_ID_phi == layergeom->get_max_phi_bin_in_sec() / 2)
204  {
205  //assign eta min according phi bin 0
206  map_z_tower_z_ID[tower.centralZ] = tower_ID_z;
207  }
208  // ...
209  } // BOOST_FOREACH(const PHG4CylinderGeom_Spacalv3::tower_map_t::value_type& tower_pair, tower_map)
210 
211  assert(! std::isnan(phi_min));
212  layerseggeo->set_phimin(phi_min);
213  layerseggeo->set_phistep(deltaphi);
214  layerseggeo->set_phibins(nphibin);
215 
217  int eta_bin = 0;
218  BOOST_FOREACH( map_z_tower_z_ID_t::value_type& z_tower_z_ID, map_z_tower_z_ID)
219  {
220  tower_z_ID_eta_bin_map[z_tower_z_ID.second] = eta_bin;
221  eta_bin++;
222  }
223  layerseggeo->set_tower_z_ID_eta_bin_map(tower_z_ID_eta_bin_map);
224  layerseggeo->set_etabins(eta_bin * layergeom->get_n_subtower_eta());
225  layerseggeo->set_etamin(NAN);
226  layerseggeo->set_etastep(NAN);
227 
228  //build eta bin maps
229  BOOST_FOREACH(const PHG4CylinderGeom_Spacalv3::tower_map_t::value_type& tower_pair, tower_map)
230  {
231 
232  const int & tower_ID = tower_pair.first;
234  tower_pair.second;
235 
236  // inspect index in sector 0
237  std::pair<int, int> tower_z_phi_ID = layergeom->get_tower_z_phi_ID(tower_ID, 0);
238  const int & tower_ID_z = tower_z_phi_ID.first;
239  const int & tower_ID_phi = tower_z_phi_ID.second;
240  const int & etabin = tower_z_ID_eta_bin_map[tower_ID_z];
241 
242  if (tower_ID_phi == layergeom->get_max_phi_bin_in_sec() / 2)
243  {
244  // half z-range
245  const double dz = fabs(0.5 * (tower.pDy1 + tower.pDy2) / sin(tower.pRotationAngleX));
246  const double tower_radial = layergeom->get_tower_radial_position(tower);
247 
248  auto z_to_eta = [&tower_radial](const double &z){return -log(tan(0.5 * atan2(tower_radial, z)));};
249 
250  const double eta_central = z_to_eta(tower.centralZ);
251  // half eta-range
252  const double deta = (z_to_eta( tower.centralZ + dz) - z_to_eta( tower.centralZ - dz))/2;
253  assert(deta > 0);
254 
255  for (int sub_tower_ID_y = 0; sub_tower_ID_y < tower.NSubtowerY;
256  ++sub_tower_ID_y)
257  {
258  assert(tower.NSubtowerY <=layergeom->get_n_subtower_eta());
259  // do not overlap to the next bin.
260  const int sub_tower_etabin = etabin * layergeom->get_n_subtower_eta() + sub_tower_ID_y;
261 
262  const pair<double, double>etabounds (eta_central - deta + sub_tower_ID_y * 2 * deta / tower.NSubtowerY,
263  eta_central - deta + (sub_tower_ID_y + 1) * 2 * deta / tower.NSubtowerY);
264 
265  const pair<double, double>zbounds (tower.centralZ - dz + sub_tower_ID_y * 2 * dz / tower.NSubtowerY,
266  tower.centralZ - dz + (sub_tower_ID_y + 1) * 2 * dz / tower.NSubtowerY);
267 
268  layerseggeo->set_etabounds(sub_tower_etabin,etabounds);
269  layerseggeo->set_zbounds(sub_tower_etabin,zbounds);
270 
271  if (Verbosity() >= VERBOSITY_SOME)
272  {
273  cout << "PHG4FullProjSpacalCellReco::InitRun::" << Name()
274  << "\t tower_ID_z = " << tower_ID_z
275  << "\t tower_ID_phi = " << tower_ID_phi
276  << "\t sub_tower_ID_y = " << sub_tower_ID_y
277  << "\t sub_tower_etabin = " << sub_tower_etabin
278  << "\t dz = " << dz
279  << "\t tower_radial = " << tower_radial
280  << "\t eta_central = " << eta_central
281  << "\t deta = " << deta
282  << "\t etabounds = [" <<etabounds.first << ", " << etabounds.second<<"]"
283  << "\t zbounds = [" <<zbounds.first << ", " << zbounds.second<<"]"
284  <<endl;
285  }
286 
287  }
288 
289  }
290  // ...
291  } // BOOST_FOREACH(const PHG4CylinderGeom_Spacalv3::tower_map_t::value_type& tower_pair, tower_map)
292 
293  // add geo object filled by different binning methods
294  seggeo->AddLayerCellGeom(layerseggeo);
295  if (Verbosity() >= VERBOSITY_SOME)
296  {
297  cout << "PHG4FullProjSpacalCellReco::InitRun::" << Name()
298  << " - Done layer" << (layergeom->get_layer())
299  << ". Print out the cell geometry:" << endl;
300  layerseggeo->identify();
301  }
303  // save this to the run wise tree to store on DST
304  PHCompositeNode *runNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "RUN" ));
305  PHCompositeNode *parNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "PAR" ));
306  PHNodeIterator runIter(runNode);
307  PHCompositeNode *RunDetNode = dynamic_cast<PHCompositeNode*>(runIter.findFirst("PHCompositeNode",detector));
308  if (! RunDetNode)
309  {
310  RunDetNode = new PHCompositeNode(detector);
311  runNode->addNode(RunDetNode);
312  }
313  string paramnodename = "G4CELLPARAM_" + detector;
314  SaveToNodeTree(RunDetNode,paramnodename);
315  // save this to the parNode for use
316  PHNodeIterator parIter(parNode);
317  PHCompositeNode *ParDetNode = dynamic_cast<PHCompositeNode*>(parIter.findFirst("PHCompositeNode",detector));
318  if (! ParDetNode)
319  {
320  ParDetNode = new PHCompositeNode(detector);
321  parNode->addNode(ParDetNode);
322  }
323  string geonodename = "G4CELLGEO_" + detector;
324  PutOnParNode(ParDetNode,geonodename);
325  tmin = get_double_param("tmin");
326  tmax = get_double_param("tmax");
327 
329 }
330 
331 int
333 {
334  PHG4HitContainer *g4hit = findNode::getClass<PHG4HitContainer>(topNode, hitnodename.c_str());
335  if (!g4hit)
336  {
337  cout << "PHG4FullProjSpacalCellReco::process_event - Fatal Error - Could not locate g4 hit node "
338  << hitnodename << endl;
339  exit(1);
340  }
341  PHG4CellContainer *cells = findNode::getClass<PHG4CellContainer>(topNode, cellnodename);
342  if (!cells)
343  {
344  cout << "PHG4FullProjSpacalCellReco::process_event - Fatal Error - could not locate cell node "
345  << cellnodename << endl;
346  exit(1);
347  }
348 
349  PHG4CylinderCellGeomContainer *seggeo = findNode::getClass<PHG4CylinderCellGeomContainer>(topNode, seggeonodename.c_str());
350  if (!seggeo)
351  {
352  cout << "PHG4FullProjSpacalCellReco::process_event - Fatal Error - could not locate cell geometry node "
353  << seggeonodename << endl;
354  exit(1);
355  }
356 
357  PHG4CylinderGeomContainer *layergeo = findNode::getClass<PHG4CylinderGeomContainer>(topNode, geonodename.c_str());
358  if (!layergeo)
359  {
360  cout << "PHG4FullProjSpacalCellReco::process_event - Fatal Error - Could not locate sim geometry node "
361  << geonodename << endl;
362  exit(1);
363  }
364 
365 
367  PHG4HitContainer::ConstRange hit_begin_end = g4hit->getHits();
368 
369  PHG4CylinderCellGeom *geo_raw = seggeo->GetFirstLayerCellGeom();
370  PHG4CylinderCellGeom_Spacalv1 * geo = dynamic_cast<PHG4CylinderCellGeom_Spacalv1 *>(geo_raw);
371  assert(geo);
372  const PHG4CylinderGeom *layergeom_raw = layergeo->GetFirstLayerGeom();
373  assert(layergeom_raw);
374  // a special implimentation of PHG4CylinderGeom is required here.
375  const PHG4CylinderGeom_Spacalv3 *layergeom =
376  dynamic_cast<const PHG4CylinderGeom_Spacalv3 *>(layergeom_raw);
377  assert(layergeom);
378 
379  for (hiter = hit_begin_end.first; hiter != hit_begin_end.second; ++hiter)
380  {
381  // checking ADC timing integration window cut
382  if (hiter->second->get_t(0)>tmax) continue;
383  if (hiter->second->get_t(1)<tmin) continue;
384 
385  sum_energy_g4hit += hiter->second->get_edep();
386 
387  // hit loop
388  int scint_id = hiter->second->get_scint_id();
389 
390  // decode scint_id
392 
393  // convert to z_ID, phi_ID
394  std::pair<int, int> tower_z_phi_ID = layergeom->get_tower_z_phi_ID(decoder.tower_ID, decoder.sector_ID);
395  const int & tower_ID_z = tower_z_phi_ID.first;
396  const int & tower_ID_phi = tower_z_phi_ID.second;
397 
398  PHG4CylinderGeom_Spacalv3::tower_map_t::const_iterator it_tower =
399  layergeom->get_sector_tower_map().find(decoder.tower_ID);
400  assert(it_tower != layergeom->get_sector_tower_map().end());
401 
402  unsigned int key = static_cast<unsigned int>(scint_id);
403  PHG4Cell *cell = nullptr;
404  map<unsigned int, PHG4Cell *>::iterator it = celllist.find(key);
405  if (it != celllist.end())
406  {
407  cell = it->second;
408  }
409  else
410  {
411 
412  // convert tower_ID_z to to eta bin number
413  int etabin = -1;
414  try
415  {
416  etabin = geo->get_etabin_block(tower_ID_z); // block eta bin
417  }
418  catch (exception & e)
419  {
420  cout << "Print cell geometry:" << endl;
421  geo->identify();
422  cout << "Print scint_id_coder:" << endl;
423  decoder.identify();
424  cout << "Print the hit:" << endl;
425  hiter->second->print();
426  cout << "PHG4FullProjSpacalCellReco::process_event::"
427  << Name() << " - Fatal Error - " << e.what() << endl;
428  exit(1);
429  }
430 
431  const int sub_tower_ID_x = it_tower->second.get_sub_tower_ID_x(decoder.fiber_ID);
432  const int sub_tower_ID_y = it_tower->second.get_sub_tower_ID_y(decoder.fiber_ID);
433  unsigned short fiber_ID = decoder.fiber_ID;
434  unsigned short etabinshort = etabin * layergeom->get_n_subtower_eta() + sub_tower_ID_y;
435  unsigned short phibin = tower_ID_phi * layergeom->get_n_subtower_phi() + sub_tower_ID_x;
436  PHG4CellDefs::keytype cellkey = PHG4CellDefs::SpacalBinning::genkey(etabinshort,phibin,fiber_ID);
437  cell = new PHG4Cellv1(cellkey);
438  celllist[key] = cell;
439  }
440 
441  double light_yield = hiter->second->get_light_yield();
442 
443  // light yield correction from fiber attenuation:
445  {
446  const double z = 0.5
447  * (hiter->second->get_local_z(0)
448  + hiter->second->get_local_z(1));
449  assert(not std::isnan(z));
450 
452  }
453 
454  // light yield correction from light guide collection efficiency:
456  {
457  const double x = it_tower->second.get_position_fraction_x_in_sub_tower(decoder.fiber_ID);
458  const double y = it_tower->second.get_position_fraction_y_in_sub_tower(decoder.fiber_ID);
459 
461  }
462 
463  cell->add_edep(hiter->first, hiter->second->get_edep());
464  cell->add_edep(hiter->second->get_edep());
465  cell->add_light_yield(light_yield);
466  cell->add_shower_edep(hiter->second->get_shower_id(), hiter->second->get_edep());
467 
468  } // end loop over g4hits
469  int numcells = 0;
470  for (map<unsigned int, PHG4Cell *>::const_iterator mapiter =
471  celllist.begin(); mapiter != celllist.end(); ++mapiter)
472  {
473  cells->AddCell(mapiter->second);
474  numcells++;
475  if (Verbosity() > 1)
476  {
477  cout << "PHG4FullProjSpacalCellReco::process_event::" << Name()
478  << " - " << "Adding cell in bin eta "
479  << PHG4CellDefs::SpacalBinning::get_etabin(mapiter->second->get_cellid())
480  << " phi "
481  << PHG4CellDefs::SpacalBinning::get_phibin(mapiter->second->get_cellid())
482  << " fiber "
483  << PHG4CellDefs::SpacalBinning::get_fiberid(mapiter->second->get_cellid())
484  << ", energy dep: "
485  << mapiter->second->get_edep() << ", light yield: "
486  << mapiter->second->get_light_yield() << endl;
487  }
488  }
489  celllist.clear();
490  if (Verbosity() > 0)
491  {
492  cout << "PHG4FullProjSpacalCellReco::process_event::" << Name()
493  << " - " << " found " << numcells
494  << " fibers with energy deposition" << endl;
495  }
496 
497 
498  if (chkenergyconservation || Verbosity() > 4)
499  {
500  CheckEnergy(topNode);
501  }
503 }
504 
505 int
507 {
508  PHG4CellContainer *cells = findNode::getClass<PHG4CellContainer>(topNode, cellnodename);
509  double sum_energy_cells = 0.;
510  PHG4CellContainer::ConstRange cell_begin_end = cells->getCells();
512  for (citer = cell_begin_end.first; citer != cell_begin_end.second; ++citer)
513  {
514  sum_energy_cells += citer->second->get_edep();
515 
516  }
517  // the fractional eloss for particles traversing eta bins leads to minute rounding errors
518  if (fabs(sum_energy_cells - sum_energy_g4hit) / sum_energy_g4hit > 1e-6)
519  {
520  cout
521  << "PHG4FullProjSpacalCellReco::CheckEnergy - energy mismatch between cells: "
522  << sum_energy_cells << " and hits: " << sum_energy_g4hit
523  << " diff sum(cells) - sum(hits): "
524  << sum_energy_cells - sum_energy_g4hit << endl;
525  return -1;
526  }
527  else
528  {
529  if (Verbosity() > 0)
530  {
531  cout << "PHG4FullProjSpacalCellReco::CheckEnergy::" << Name()
532  << " - total energy for this event: " << sum_energy_g4hit
533  << " GeV. Passed CheckEnergy" << endl;
534  }
535  }
536  return 0;
537 }
538 
540  data_grid_light_guide_efficiency(nullptr),
541  data_grid_fiber_trans(nullptr)
542 {
543 
544  data_grid_light_guide_efficiency_verify = new TH2F("data_grid_light_guide_efficiency_verify",
545  "light collection efficiency as used in PHG4FullProjSpacalCellReco::LightCollectionModel;x positio fraction;y position fraction", //
546  100, 0., 1., 100, 0., 1.);
547 
548  data_grid_fiber_trans_verify = new TH1F("data_grid_fiber_trans",
549  "SCSF-78 Fiber Transmission as used in PHG4FullProjSpacalCellReco::LightCollectionModel;position in fiber (cm);Effective transmission",
550  100, -15, 15);
551 
552  // register histograms
554 
557 
558 }
559 
561 {
562  delete data_grid_light_guide_efficiency;
563  delete data_grid_fiber_trans;
564 }
565 
566 void
568  const std::string & input_file,
569  const std::string & histogram_light_guide_model,
570  const std::string & histogram_fiber_model)
571 {
572  TFile * fin = TFile::Open(input_file.c_str());
573 
574  assert(fin);
575  assert(fin->IsOpen());
576 
577  if (data_grid_light_guide_efficiency) delete data_grid_light_guide_efficiency;
578  data_grid_light_guide_efficiency = dynamic_cast<TH2 *>(fin->Get(histogram_light_guide_model.c_str()));
579  assert(data_grid_light_guide_efficiency);
580  data_grid_light_guide_efficiency->SetDirectory(nullptr);
581 
582  if (data_grid_fiber_trans) delete data_grid_fiber_trans;
583  data_grid_fiber_trans = dynamic_cast<TH1 *>(fin->Get(histogram_fiber_model.c_str()));
584  assert(data_grid_fiber_trans);
585  data_grid_fiber_trans->SetDirectory(nullptr);
586 
587  delete fin;
588 }
589 
590 double
592 {
593  assert(data_grid_light_guide_efficiency);
594  assert(x_fraction >= 0);
595  assert(x_fraction <= 1);
596  assert(y_fraction >= 0);
597  assert(y_fraction <= 1);
598 
599  const double eff = data_grid_light_guide_efficiency->Interpolate(x_fraction,
600  y_fraction);
601 
602  data_grid_light_guide_efficiency_verify->SetBinContent( //
603  data_grid_light_guide_efficiency_verify->GetXaxis()->FindBin(x_fraction), //
604  data_grid_light_guide_efficiency_verify->GetYaxis()->FindBin(y_fraction), //
605  eff //
606  );
607 
608  return eff;
609 }
610 
611 double
613 {
614  assert(data_grid_fiber_trans);
615 
616  const double eff = data_grid_fiber_trans->Interpolate(z_distance);
617 
618  data_grid_fiber_trans_verify->SetBinContent( //
619  data_grid_fiber_trans_verify->GetXaxis()->FindBin(z_distance), //
620  eff //
621  );
622 
623  return eff;
624 }
625 
626 void
628 {
629  set_default_double_param("tmax",60.0);
630  set_default_double_param("tmin",-20.0); // collision has a timing spread around the triggered event. Accepting negative time too.
631  return;
632 }
633 
634 void
635 PHG4FullProjSpacalCellReco::set_timing_window(const double tmi, const double tma)
636 {
637  set_double_param("tmin",tmi);
638  set_double_param("tmax",tma);
639 }