LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
PhotonHitConverter_module.cc
Go to the documentation of this file.
1 // Class: PhotonHitConverter
3 // Plugin Type: producer (Unknown Unknown)
4 // File: PhotonHitConverter_module.cc
5 //
6 // Generated at Tue Dec 13 05:30:41 2022 by Alejandro Castillo using cetskelgen
7 // from version .
9 
17 #include "fhiclcpp/ParameterSet.h"
19 
24 
27 
28 #include "artg4tk/pluginDetectors/gdml/PhotonHit.hh"
29 
30 #include <memory>
31 #include <string>
32 #include <vector>
33 
34 namespace sim {
35  class PhotonHitConverter;
36 }
37 
39 public:
40  explicit PhotonHitConverter(fhicl::ParameterSet const& p);
41  // The compiler-generated destructor is fine for non-base
42  // classes without bare pointers or other resource use.
43  // Plugins should not be copied or assigned.
44  PhotonHitConverter(PhotonHitConverter const&) = delete;
48  // Required functions.
49  void produce(art::Event& e) override;
50 
51 private:
52  // Declare member data here.
53  const bool fUseLitePhotons;
54  size_t nOpChannels;
55 };
56 
58  : EDProducer{p}, fUseLitePhotons(p.get<bool>("UseLitePhotons")) // ,
59 // More initializers here.
60 {
61  // Call appropriate produces<>() functions here.
62  if (fUseLitePhotons) {
63  produces<std::vector<sim::SimPhotonsLite>>();
64  produces<std::vector<sim::SimPhotonsLite>>("Reflected");
65  }
66  else {
67  produces<std::vector<sim::SimPhotons>>();
68  produces<std::vector<sim::SimPhotons>>("Reflected");
69  }
70 }
71 
73 {
74  geo::GeometryCore const& geom = *(lar::providerFrom<geo::Geometry>());
75  nOpChannels = geom.NOpDets();
76 
77  //SimPhotonsLite
78  std::unique_ptr<std::vector<sim::SimPhotonsLite>> photLiteCol{
79  new std::vector<sim::SimPhotonsLite>{}};
80  std::unique_ptr<std::vector<sim::SimPhotonsLite>> photLiteCol_ref{
81  new std::vector<sim::SimPhotonsLite>{}};
82  auto& photonLiteCollection{*photLiteCol};
83  auto& photonLiteCollection_ref{*photLiteCol_ref};
84  //SimPhotons
85  std::unique_ptr<std::vector<sim::SimPhotons>> photCol{new std::vector<sim::SimPhotons>{}};
86  std::unique_ptr<std::vector<sim::SimPhotons>> photCol_ref{new std::vector<sim::SimPhotons>{}};
87  auto& photonCollection{*photCol};
88  auto& photonCollection_ref{*photCol_ref};
89 
90  if (fUseLitePhotons) { //SimPhotonsLite
91  photonLiteCollection.resize(nOpChannels);
92  photonLiteCollection_ref.resize(nOpChannels);
93  for (unsigned int i = 0; i < nOpChannels; ++i) {
94  photonLiteCollection[i].OpChannel = i;
95  photonLiteCollection_ref[i].OpChannel = i;
96  }
97  }
98  else { //SimPhotons
99  photonCollection.resize(nOpChannels);
100  photonCollection_ref.resize(nOpChannels);
101  for (unsigned int i = 0; i < nOpChannels; ++i) {
102  photonCollection[i].fOpChannel = i;
103  photonCollection_ref[i].fOpChannel = i;
104  }
105  }
106 
107  typedef std::vector<art::Handle<artg4tk::PhotonHitCollection>> HandleVector;
108  auto allSims = e.getMany<artg4tk::PhotonHitCollection>();
109  for (HandleVector::const_iterator i = allSims.begin(); i != allSims.end(); ++i) {
110  const artg4tk::PhotonHitCollection& sims(**i);
111  for (artg4tk::PhotonHitCollection::const_iterator j = sims.begin(); j != sims.end(); ++j) {
112  const artg4tk::PhotonHit& photonHit = *j;
113  if (fUseLitePhotons) {
114  if (photonHit.GetEdep() > 6.19 * CLHEP::eV) {
115  auto time = static_cast<int>(photonHit.GetTime());
116  auto channel = static_cast<unsigned int>(photonHit.GetID());
117  ++photonLiteCollection[channel].DetectedPhotons[time];
118  }
119  else {
120  auto time = static_cast<int>(photonHit.GetTime());
121  auto channel = static_cast<unsigned int>(photonHit.GetID());
122  ++photonLiteCollection_ref[channel].DetectedPhotons[time];
123  }
124  }
125  else {
126  sim::OnePhoton photon;
127  photon.SetInSD = false;
128  photon.Energy = photonHit.GetEdep();
129  auto time = photonHit.GetTime();
130  photon.Time = time;
131  auto channel = static_cast<unsigned int>(photonHit.GetID());
132  if (photon.Energy > 6.19 * CLHEP::eV) {
133  photonCollection[channel].insert(photonCollection[channel].end(), 1, photon);
134  }
135  else {
136  photonCollection_ref[channel].insert(photonCollection_ref[channel].end(), 1, photon);
137  }
138  }
139  }
140  }
141  if (fUseLitePhotons) {
142  e.put(std::move(photLiteCol));
143  e.put(std::move(photLiteCol_ref), "Reflected");
144  }
145  else {
146  e.put(std::move(photCol));
147  e.put(std::move(photCol_ref), "Reflected");
148  }
149 }
150 
Utilities related to art service access.
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.cc:6
All information of a photon entering the sensitive optical detector volume.
Definition: SimPhotons.h:60
intermediate_table::const_iterator const_iterator
PhotonHitConverter & operator=(PhotonHitConverter const &)=delete
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Event.h:77
decltype(auto) constexpr end(T &&obj)
ADL-aware version of std::end.
Definition: StdUtils.h:77
Simulation objects for optical detectors.
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
Description of geometry of one entire detector.
Definition: GeometryCore.h:119
Monte Carlo Simulation.
unsigned int NOpDets() const
Number of OpDets in the whole detector.
Provides a base class aware of world box coordinates.
Encapsulate the geometry of an optical detector.
PhotonHitConverter(fhicl::ParameterSet const &p)
bool SetInSD
Whether the photon reaches the sensitive detector.
Definition: SimPhotons.h:84
float Energy
Scintillation photon energy [GeV].
Definition: SimPhotons.h:78
Float_t e
Definition: plot.C:35
art framework interface to geometry description
void produce(art::Event &e) override
std::vector< Handle< PROD > > getMany(SelectorBase const &selector=MatchAllSelector{}) const