LArSoft  v10_04_05
Liquid Argon Software toolkit - https://larsoft.org/
SimPhotonCounter_module.cc
Go to the documentation of this file.
1 // \file SimPhotonCounter.h
2 // \author Ben Jones, MIT 2010
3 //
4 // Module to determine how many phots have been detected at each OpDet
5 //
6 // This analyzer takes the SimPhotonsCollection generated by LArG4's sensitive detectors
7 // and fills up to four trees in the histograms file. The four trees are:
8 //
9 // OpDetEvents - count how many phots hit the OpDet face / were detected across all OpDet's per event
10 // OpDets - count how many phots hit the OpDet face / were detected in each OpDet individually for each event
11 // AllPhotons - wavelength information for each phot hitting the OpDet face
12 // DetectedPhotons - wavelength information for each phot detected
13 //
14 // The user may supply a quantum efficiency and sensitive wavelength range for the OpDet's.
15 // with a QE < 1 and a finite wavelength range, a "detected" phot is one which is
16 // in the relevant wavelength range and passes the random sampling condition imposed by
17 // the quantum efficiency of the OpDet
18 //
19 // PARAMETERS REQUIRED:
20 // int32 Verbosity - whether to write to screen a well as to file. levels 0 to 3 specify different levels of detail to display
21 // string InputModule - the module which produced the SimPhotonsCollection
22 // bool MakeAllPhotonsTree - whether to build and store each tree (performance can be enhanced by switching off those not required)
23 // bool MakeDetectedPhotonsTree
24 // bool MakeOpDetsTree
25 // bool MakeOpDetEventsTree
26 // double QantumEfficiency - Quantum efficiency of OpDet
27 // double WavelengthCutLow - Sensitive wavelength range of OpDet
28 // double WavelengthCutHigh
29 
30 // FMWK includes
36 #include "art_root_io/TFileService.h"
37 #include "fhiclcpp/ParameterSet.h"
39 
40 // LArSoft includes
53 
56 
57 // ROOT includes
58 #include "RtypesCore.h"
59 #include "TH1D.h"
60 #include "TLorentzVector.h"
61 #include "TTree.h"
62 #include "TVector3.h"
63 
64 // C++ language includes
65 #include <cassert>
66 #include <cstring>
67 #include <iostream>
68 
69 namespace opdet {
70 
71  class SimPhotonCounter : public art::EDAnalyzer {
72  public:
74 
75  void analyze(art::Event const&);
76 
77  void beginJob();
78  void endJob();
79 
80  private:
82  static constexpr double kVisibleThreshold = 200.0; // nm
83 
85  static constexpr double kVisibleWavelength = 450.0; // nm
86 
88  static constexpr double kVUVWavelength = 128.0; // nm
89 
90  // Trees to output
91 
94  TTree* fTheOpDetTree;
95  TTree* fTheEventTree;
96 
97  // Parameters to read in
98 
99  std::vector<std::string> fInputModule; // Input tag for OpDet collection
100 
101  int fVerbosity; // Level of output to write to std::out
102 
105  bool fMakeOpDetsTree; // Switches to turn on or off each output
108 
109  // float fQE; // Quantum efficiency of tube
110 
111  // float fWavelengthCutLow; // Sensitive wavelength range
112  // float fWavelengthCutHigh; //
113 
116 
117  // Data to store in trees
118 
119  Float_t fWavelength;
120  Float_t fTime;
121  Float_t foriginX;
122  Float_t foriginY;
123  Float_t foriginZ;
124  // Int_t fCount;
128  Float_t fT0_vis;
129 
132  // Int_t fCountEventDetectedwithRefl;
133 
134  Int_t fEventID;
135  Int_t fOpChannel;
136 
137  //for the analysis tree of the light (gamez)
139  std::vector<std::vector<std::vector<double>>> fSignals_vuv;
140  std::vector<std::vector<std::vector<double>>> fSignals_vis;
141 
142  TTree* fLightAnalysisTree = nullptr;
144  double fEnergy, fdEdx;
145  std::vector<double> fPosition0;
146  std::vector<std::vector<double>> fstepPositions;
147  std::vector<double> fstepTimes;
148  std::vector<std::vector<double>> fSignalsvuv;
149  std::vector<std::vector<double>> fSignalsvis;
150  std::string fProcess;
151 
154 
155  bool const fUseLitePhotons;
156 
157  void storeVisibility(int channel,
158  int nDirectPhotons,
159  int nReflectedPhotons,
160  double reflectedT0 = 0.0) const;
161 
163  bool isVisible(double wavelength) const { return fWavelength < kVisibleThreshold; }
164  };
165 }
166 
167 namespace opdet {
168 
170  : EDAnalyzer(pset)
173  {
174  fVerbosity = pset.get<int>("Verbosity");
175  try {
176  fInputModule = pset.get<std::vector<std::string>>("InputModule", {"largeant"});
177  }
178  catch (...) {
179  fInputModule.push_back(pset.get<std::string>("InputModule", "largeant"));
180  }
181  fMakeAllPhotonsTree = pset.get<bool>("MakeAllPhotonsTree");
182  fMakeDetectedPhotonsTree = pset.get<bool>("MakeDetectedPhotonsTree");
183  fSavePhotonOrigins = pset.get<bool>("SavePhotonOrigins", false);
184  fMakeOpDetsTree = pset.get<bool>("MakeOpDetsTree");
185  fMakeOpDetEventsTree = pset.get<bool>("MakeOpDetEventsTree");
186  fMakeLightAnalysisTree = pset.get<bool>("MakeLightAnalysisTree", false);
187  //fQE= pset.get<double>("QuantumEfficiency");
188  //fWavelengthCutLow= pset.get<double>("WavelengthCutLow");
189  //fWavelengthCutHigh= pset.get<double>("WavelengthCutHigh");
190 
193  << "Building a library with reflected light time is not supported when using "
194  "SimPhotonsLite.\n";
195  }
196  }
197 
199  {
200  // Get file service to store trees
203 
204  std::cout << "Optical Channels positions: " << geo->Cryostat().NOpDet() << std::endl;
205  for (int ch = 0; ch != int(geo->Cryostat().NOpDet()); ch++) {
206  auto const OpDetCenter = geo->OpDetGeoFromOpDet(ch).GetCenter();
207  std::cout << ch << " " << OpDetCenter.X() << " " << OpDetCenter.Y() << " "
208  << OpDetCenter.Z() << std::endl;
209  }
210 
211  auto const CryoBounds = geo->Cryostat().Boundaries();
212  std::cout << "Cryo Boundaries" << std::endl;
213  std::cout << "Xmin: " << CryoBounds.MinX() << " Xmax: " << CryoBounds.MaxX()
214  << " Ymin: " << CryoBounds.MinY() << " Ymax: " << CryoBounds.MaxY()
215  << " Zmin: " << CryoBounds.MinZ() << " Zmax: " << CryoBounds.MaxZ() << std::endl;
216 
217  try {
219  }
220  catch (art::Exception const& e) {
221  if (e.categoryCode() != art::errors::ServiceNotFound) throw;
222  mf::LogError("SimPhotonCounter")
223  << "ParticleInventoryService service is not configured!"
224  " Please add it in the job configuration."
225  " In the meanwhile, some checks to particles will be skipped.";
226  }
227 
228  // Create and assign branch addresses to required tree
229  if (fMakeAllPhotonsTree) {
230  fThePhotonTreeAll = tfs->make<TTree>("AllPhotons", "AllPhotons");
231  fThePhotonTreeAll->Branch("EventID", &fEventID, "EventID/I");
232  fThePhotonTreeAll->Branch("Wavelength", &fWavelength, "Wavelength/F");
233  fThePhotonTreeAll->Branch("OpChannel", &fOpChannel, "OpChannel/I");
234  fThePhotonTreeAll->Branch("Time", &fTime, "Time/F");
235  if (fSavePhotonOrigins) {
236  fThePhotonTreeAll->Branch("originX", &foriginX, "originX/F");
237  fThePhotonTreeAll->Branch("originY", &foriginY, "originY/F");
238  fThePhotonTreeAll->Branch("originZ", &foriginZ, "originZ/F");
239  }
240  }
241 
243  fThePhotonTreeDetected = tfs->make<TTree>("DetectedPhotons", "DetectedPhotons");
244  fThePhotonTreeDetected->Branch("EventID", &fEventID, "EventID/I");
245  fThePhotonTreeDetected->Branch("Wavelength", &fWavelength, "Wavelength/F");
246  fThePhotonTreeDetected->Branch("OpChannel", &fOpChannel, "OpChannel/I");
247  fThePhotonTreeDetected->Branch("Time", &fTime, "Time/F");
248  if (fSavePhotonOrigins) {
249  fThePhotonTreeDetected->Branch("originX", &foriginX, "originX/F");
250  fThePhotonTreeDetected->Branch("originY", &foriginY, "originY/F");
251  fThePhotonTreeDetected->Branch("originZ", &foriginZ, "originZ/F");
252  }
253  }
254 
255  if (fMakeOpDetsTree) {
256  fTheOpDetTree = tfs->make<TTree>("OpDets", "OpDets");
257  fTheOpDetTree->Branch("EventID", &fEventID, "EventID/I");
258  fTheOpDetTree->Branch("OpChannel", &fOpChannel, "OpChannel/I");
259  fTheOpDetTree->Branch("CountAll", &fCountOpDetAll, "CountAll/I");
260  fTheOpDetTree->Branch("CountDetected", &fCountOpDetDetected, "CountDetected/I");
261  if (fPVS->StoreReflected())
262  fTheOpDetTree->Branch("CountReflDetected", &fCountOpDetReflDetected, "CountReflDetected/I");
263  fTheOpDetTree->Branch("Time", &fTime, "Time/F");
264  }
265 
266  if (fMakeOpDetEventsTree) {
267  fTheEventTree = tfs->make<TTree>("OpDetEvents", "OpDetEvents");
268  fTheEventTree->Branch("EventID", &fEventID, "EventID/I");
269  fTheEventTree->Branch("CountAll", &fCountEventAll, "CountAll/I");
270  fTheEventTree->Branch("CountDetected", &fCountEventDetected, "CountDetected/I");
271  if (fPVS->StoreReflected())
272  fTheOpDetTree->Branch("CountReflDetected", &fCountOpDetReflDetected, "CountReflDetected/I");
273  }
274 
275  //generating the tree for the light analysis:
277  fLightAnalysisTree = tfs->make<TTree>("LightAnalysis", "LightAnalysis");
278  fLightAnalysisTree->Branch("RunNumber", &fRun);
279  fLightAnalysisTree->Branch("EventID", &fEventID);
280  fLightAnalysisTree->Branch("TrackID", &fTrackID);
281  fLightAnalysisTree->Branch("PdgCode", &fpdg);
282  fLightAnalysisTree->Branch("MotherTrackID", &fmotherTrackID);
283  fLightAnalysisTree->Branch("Energy", &fEnergy);
284  fLightAnalysisTree->Branch("dEdx", &fdEdx);
285  fLightAnalysisTree->Branch("StepPositions", &fstepPositions);
286  fLightAnalysisTree->Branch("StepTimes", &fstepTimes);
287  fLightAnalysisTree->Branch("SignalsVUV", &fSignalsvuv);
288  fLightAnalysisTree->Branch("SignalsVisible", &fSignalsvis);
289  fLightAnalysisTree->Branch("Process", &fProcess);
290  }
291  }
292 
294  {
296  }
297 
299  {
300  // Lookup event ID from event
301  art::EventNumber_t event = evt.id().event();
302  fEventID = Int_t(event);
303 
304  // Service for determining opdet responses
306 
307  // get the geometry to be able to figure out signal types and chan -> plane mappings
308  auto const& wireReadoutGeom = art::ServiceHandle<geo::WireReadout const>()->Get();
309 
310  // GEANT4 info on the particles (only used if making light analysis tree)
311  std::vector<simb::MCParticle> const* mcpartVec = nullptr;
312 
313  //-------------------------initializing light tree vectors------------------------
314  std::vector<double> totalEnergy_track;
315  fstepPositions.clear();
316  fstepTimes.clear();
318  mcpartVec = evt.getHandle<std::vector<simb::MCParticle>>("largeant").product();
319 
320  size_t maxNtracks = 1000U; // mcpartVec->size(); --- { to be fixed soon! ]
321  fSignals_vuv.clear();
322  fSignals_vuv.resize(maxNtracks);
323  fSignals_vis.clear();
324  fSignals_vis.resize(maxNtracks);
325  for (size_t itrack = 0; itrack != maxNtracks; itrack++) {
326  fSignals_vuv[itrack].resize(wireReadoutGeom.NOpChannels());
327  fSignals_vis[itrack].resize(wireReadoutGeom.NOpChannels());
328  }
329  totalEnergy_track.resize(maxNtracks, 0.);
330  //-------------------------stimation of dedx per trackID----------------------
331  //get the list of particles from this event
332  const sim::ParticleList* plist = pi_serv ? &(pi_serv->ParticleList()) : nullptr;
333 
334  // loop over all sim::SimChannels in the event and make sure there are no
335  // sim::IDEs with trackID values that are not in the sim::ParticleList
336  std::vector<const sim::SimChannel*> sccol;
337  for (auto const& mod : fInputModule) {
338  evt.getView(mod, sccol);
339  //loop over the sim channels collection
340  for (size_t sc = 0; sc < sccol.size(); ++sc) {
341  const auto& tdcidemap = sccol[sc]->TDCIDEMap();
342  //loop over all of the tdc IDE map objects
343  for (auto mapitr = tdcidemap.begin(); mapitr != tdcidemap.end(); mapitr++) {
344  const std::vector<sim::IDE> idevec = (*mapitr).second;
345  //go over all of the IDEs in a given simchannel
346  for (size_t iv = 0; iv < idevec.size(); ++iv) {
347  if (plist) {
348  if (plist->find(idevec[iv].trackID) == plist->end() &&
349  idevec[iv].trackID != sim::NoParticleId) {
350  mf::LogWarning("LArG4Ana") << idevec[iv].trackID << " is not in particle list";
351  }
352  }
353  if (idevec[iv].trackID < 0) continue;
354 
355  totalEnergy_track[idevec[iv].trackID] += idevec[iv].energy / 3.;
356  }
357  }
358  }
359  }
360  } //End of if(fMakeLightAnalysisTree)
361 
362  if (!fUseLitePhotons) {
363 
364  //Reset counters
365  fCountEventAll = 0;
367 
368  //Get *ALL* SimPhotonsCollection from Event
369  auto photon_handles = evt.getMany<std::vector<sim::SimPhotons>>();
370  if (photon_handles.size() == 0)
372  << "sim SimPhotons retrieved and you requested them.";
373 
374  for (auto const& mod : fInputModule) {
375  //switching off to add reading in of labelled collections: Andrzej, 02/26/19
376 
377  for (auto const& ph_handle : photon_handles) {
378  // Do some checking before we proceed
379  if (!ph_handle.isValid()) continue;
380  if (ph_handle.provenance()->moduleLabel() != mod)
381  continue; //not the most efficient way of doing this, but preserves the logic of the module. Andrzej
382 
383  bool Reflected = (ph_handle.provenance()->productInstanceName() == "Reflected");
384 
385  if ((*ph_handle).size() > 0) {
387  //resetting the signalt to save in the analysis tree per event
388  const int maxNtracks = 1000;
389  for (size_t itrack = 0; itrack != maxNtracks; itrack++) {
390  for (size_t pmt_i = 0; pmt_i != wireReadoutGeom.NOpChannels(); pmt_i++) {
391  fSignals_vuv[itrack][pmt_i].clear();
392  fSignals_vis[itrack][pmt_i].clear();
393  }
394  }
395  }
396  }
397 
398  if (fVerbosity > 0)
399  std::cout << "Found OpDet hit collection of size " << (*ph_handle).size() << std::endl;
400 
401  if ((*ph_handle).size() > 0) {
402  for (auto const& itOpDet : (*ph_handle)) {
403  //Reset Counters
404  fCountOpDetAll = 0;
407  //Reset t0 for visible light
408  fT0_vis = 999.;
409 
410  //Get data from HitCollection entry
411  fOpChannel = itOpDet.OpChannel();
412  const sim::SimPhotons& TheHit = itOpDet;
413 
414  // Loop through OpDet phots.
415  // Note we make the screen output decision outside the loop
416  // in order to avoid evaluating large numbers of unnecessary
417  // if conditions.
418 
419  for (const sim::OnePhoton& Phot : TheHit) {
420  // Calculate wavelength in nm
421  fWavelength = odresponse->wavelength(Phot.Energy);
422 
423  //Get arrival time from phot
424  fTime = Phot.Time;
425  foriginX = Phot.InitialPosition.X();
426  foriginY = Phot.InitialPosition.Y();
427  foriginZ = Phot.InitialPosition.Z();
428 
429  // special case for LibraryBuildJob: no working "Reflected" handle and all photons stored in single object - must sort using wavelength instead
430  if (fPVS->IsBuildJob() && !Reflected) {
431  // all photons contained in object with Reflected = false flag
432  // Increment per OpDet counters and fill per phot trees
433  fCountOpDetAll++;
434  if (fMakeAllPhotonsTree) {
436  fThePhotonTreeAll->Fill();
437  }
438  }
439 
440  if (odresponse->detected(fOpChannel, Phot)) {
442  //only store direct direct light
444  // reflected and shifted light is in visible range
445  else if (fPVS->StoreReflected()) {
447  // find the first visible arrival time
448  if (fPVS->StoreReflT0() && fTime < fT0_vis) fT0_vis = fTime;
449  }
450  if (fVerbosity > 3)
451  std::cout << "OpDetResponseInterface PerPhoton : Event " << fEventID
452  << " OpChannel " << fOpChannel << " Wavelength " << fWavelength
453  << " Detected 1 " << std::endl;
454  }
455  else {
456  if (fVerbosity > 3)
457  std::cout << "OpDetResponseInterface PerPhoton : Event " << fEventID
458  << " OpChannel " << fOpChannel << " Wavelength " << fWavelength
459  << " Detected 0 " << std::endl;
460  }
461 
462  } // if build library and not reflected
463 
464  else {
465  // store in appropriate trees using "Reflected" handle and fPVS->StoreReflected() flag
466  // Increment per OpDet counters and fill per phot trees
467  fCountOpDetAll++;
468  if (fMakeAllPhotonsTree) {
469  if (!Reflected || (fPVS->StoreReflected() && Reflected)) {
470  fThePhotonTreeAll->Fill();
471  }
472  }
473 
474  if (odresponse->detected(fOpChannel, Phot)) {
476  //only store direct direct light
477  if (!Reflected) fCountOpDetDetected++;
478  // reflected and shifted light is in visible range
479  else if (fPVS->StoreReflected() && Reflected) {
481  // find the first visible arrival time
482  if (fPVS->StoreReflT0() && fTime < fT0_vis) fT0_vis = fTime;
483  }
484  if (fVerbosity > 3)
485  std::cout << "OpDetResponseInterface PerPhoton : Event " << fEventID
486  << " OpChannel " << fOpChannel << " Wavelength " << fWavelength
487  << " Detected 1 " << std::endl;
488  }
489  else {
490  if (fVerbosity > 3)
491  std::cout << "OpDetResponseInterface PerPhoton : Event " << fEventID
492  << " OpChannel " << fOpChannel << " Wavelength " << fWavelength
493  << " Detected 0 " << std::endl;
494  }
495  }
496  } // for each photon in collection
497 
498  // If this is a library building job, fill relevant entry
499  if (
500  fPVS->IsBuildJob() &&
501  !Reflected) // for library build job, both componenents stored in first object with Reflected = false
502  {
504  }
505 
506  // Incremenent per event and fill Per OpDet trees
507  if (fMakeOpDetsTree) fTheOpDetTree->Fill();
510 
511  // Give per OpDet output
512  if (fVerbosity > 2)
513  std::cout << "OpDetResponseInterface PerOpDet : Event " << fEventID << " OpDet "
514  << fOpChannel << " All " << fCountOpDetAll << " Det "
515  << fCountOpDetDetected << std::endl;
516  }
517 
518  // Fill per event tree
519  if (fMakeOpDetEventsTree) fTheEventTree->Fill();
520 
521  // Give per event output
522  if (fVerbosity > 1)
523  std::cout << "OpDetResponseInterface PerEvent : Event " << fEventID << " All "
524  << fCountOpDetAll << " Det " << fCountOpDetDetected << std::endl;
525  }
526  else {
527  // if empty OpDet hit collection,
528  // add an empty record to the per event tree
529  if (fMakeOpDetEventsTree) fTheEventTree->Fill();
530  }
532  assert(mcpartVec);
533  assert(fLightAnalysisTree);
534 
535  std::cout << "Filling the analysis tree" << std::endl;
536  //---------------Filling the analysis tree-----------:
537  fRun = evt.run();
538  std::vector<double> this_xyz;
539 
540  //loop over the particles
541  for (simb::MCParticle const& pPart : *mcpartVec) {
542 
543  if (pPart.Process() == "primary") fEnergy = pPart.E();
544 
545  //resetting the vectors
546  fstepPositions.clear();
547  fstepTimes.clear();
548  fSignalsvuv.clear();
549  fSignalsvis.clear();
550  fdEdx = -1.;
551  //filling the tree fields
552  fTrackID = pPart.TrackId();
553  fpdg = pPart.PdgCode();
554  fmotherTrackID = pPart.Mother();
555  fdEdx = totalEnergy_track[fTrackID];
558  fProcess = pPart.Process();
559  //filling the center positions of each step
560  for (size_t i_s = 1; i_s < pPart.NumberTrajectoryPoints(); i_s++) {
561  this_xyz.clear();
562  this_xyz.resize(3);
563  this_xyz[0] = pPart.Position(i_s).X();
564  this_xyz[1] = pPart.Position(i_s).Y();
565  this_xyz[2] = pPart.Position(i_s).Z();
566  fstepPositions.push_back(this_xyz);
567  fstepTimes.push_back(pPart.Position(i_s).T());
568  }
569  //filling the tree per track
570  fLightAnalysisTree->Fill();
571  }
572  } // if fMakeLightAnalysisTree
573  }
574  }
575  }
576  if (fUseLitePhotons) {
577 
578  //Get *ALL* SimPhotonsCollection from Event
579  auto photon_handles = evt.getMany<std::vector<sim::SimPhotonsLite>>();
580  if (photon_handles.size() == 0)
582  << "sim SimPhotons retrieved and you requested them.";
583 
584  //Get SimPhotonsLite from Event
585  for (auto const& mod : fInputModule) {
586  // Loop over direct/reflected photons
587  for (auto const& ph_handle : photon_handles) {
588  // Do some checking before we proceed
589  if (!ph_handle.isValid()) continue;
590  if (ph_handle.provenance()->moduleLabel() != mod)
591  continue; //not the most efficient way of doing this, but preserves the logic of the module. Andrzej
592 
593  bool Reflected = (ph_handle.provenance()->productInstanceName() == "Reflected");
594 
595  //Reset counters
596  fCountEventAll = 0;
598 
599  if (fVerbosity > 0)
600  std::cout << "Found OpDet hit collection of size " << (*ph_handle).size() << std::endl;
601 
602  if ((*ph_handle).size() > 0) {
603 
604  for (auto const& photon : (*ph_handle)) {
605  //Get data from HitCollection entry
606  fOpChannel = photon.OpChannel;
607  std::map<int, int> PhotonsMap = photon.DetectedPhotons;
608 
609  //Reset Counters
610  fCountOpDetAll = 0;
613 
614  for (auto it = PhotonsMap.begin(); it != PhotonsMap.end(); it++) {
615  // Calculate wavelength in nm
616  if (Reflected) { fWavelength = kVisibleThreshold; }
617  else {
618  fWavelength = kVUVWavelength; // original
619  }
620 
621  //Get arrival time from phot
622  fTime = it->first;
623 
624  for (int i = 0; i < it->second; i++) {
625  // Increment per OpDet counters and fill per phot trees
626  fCountOpDetAll++;
628 
629  if (odresponse->detectedLite(fOpChannel)) {
631  // direct light
632  if (!Reflected) { fCountOpDetDetected++; }
633  else if (Reflected) {
635  }
636  if (fVerbosity > 3)
637  std::cout << "OpDetResponseInterface PerPhoton : Event " << fEventID
638  << " OpChannel " << fOpChannel << " Wavelength " << fWavelength
639  << " Detected 1 " << std::endl;
640  }
641  else if (fVerbosity > 3) {
642  std::cout << "OpDetResponseInterface PerPhoton : Event " << fEventID
643  << " OpChannel " << fOpChannel << " Wavelength " << fWavelength
644  << " Detected 0 " << std::endl;
645  }
646  }
647  }
648 
649  // Incremenent per event and fill Per OpDet trees
650  if (fMakeOpDetsTree) fTheOpDetTree->Fill();
653 
654  if (fPVS->IsBuildJob())
656 
657  // Give per OpDet output
658  if (fVerbosity > 2)
659  std::cout << "OpDetResponseInterface PerOpDet : Event " << fEventID << " OpDet "
660  << fOpChannel << " All " << fCountOpDetAll << " Det "
661  << fCountOpDetDetected << std::endl;
662  }
663  // Fill per event tree
664  if (fMakeOpDetEventsTree) fTheEventTree->Fill();
665 
666  // Give per event output
667  if (fVerbosity > 1)
668  std::cout << "OpDetResponseInterface PerEvent : Event " << fEventID << " All "
669  << fCountOpDetAll << " Det " << fCountOpDetDetected << std::endl;
670  }
671  else {
672  // if empty OpDet hit collection,
673  // add an empty record to the per event tree
674  if (fMakeOpDetEventsTree) fTheEventTree->Fill();
675  }
676  }
677  }
678  }
679  } // SimPhotonCounter::analyze()
680 
681  // ---------------------------------------------------------------------------
683  int nDirectPhotons,
684  int nReflectedPhotons,
685  double reflectedT0 /* = 0.0 */
686  ) const
687  {
689 
690  // ask PhotonVisibilityService which voxel was being served,
691  // and how many photons where there generated (yikes!!);
692  // this value was put there by LightSource (yikes!!)
693  int VoxID;
694  double NProd;
695  fPVS->RetrieveLightProd(VoxID, NProd);
696 
697  pvs.SetLibraryEntry(VoxID, channel, double(nDirectPhotons) / NProd);
698 
699  //store reflected light
700  if (fPVS->StoreReflected()) {
701  pvs.SetLibraryEntry(VoxID, channel, double(nReflectedPhotons) / NProd, true);
702 
703  //store reflected first arrival time
704  if (fPVS->StoreReflT0()) pvs.SetLibraryReflT0Entry(VoxID, channel, reflectedT0);
705 
706  } // if reflected
707 
708  } // SimPhotonCounter::storeVisibility()
709 
710  // ---------------------------------------------------------------------------
711 
712 }
713 
Store parameters for running LArG4.
std::vector< std::vector< double > > fstepPositions
void RetrieveLightProd(int &VoxID, double &N) const
Encapsulate the construction of a single cyostat .
void storeVisibility(int channel, int nDirectPhotons, int nReflectedPhotons, double reflectedT0=0.0) const
Point_t const & GetCenter() const
Definition: OpDetGeo.h:71
double MinX() const
Returns the world x coordinate of the start of the box.
Definition: BoxBoundedGeo.h:88
std::vector< double > fPosition0
All information of a photon entering the sensitive optical detector volume.
Definition: SimPhotons.h:60
phot::PhotonVisibilityService const * fPVS
cout<< "Opened file "<< fin<< " ixs= "<< ixs<< endl;if(ixs==0) hhh=(TH1F *) fff-> Get("h1")
Definition: AddMC.C:8
std::vector< std::vector< double > > fSignalsvuv
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
virtual bool detected(int OpChannel, const sim::OnePhoton &Phot, int &newOpChannel) const
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.cc:6
static constexpr double kVUVWavelength
Value used when a typical ultraviolet light wavelength is needed.
Particle class.
std::vector< double > fstepTimes
iterator find(const key_type &key)
Definition: ParticleList.h:318
std::vector< std::string > fInputModule
Simulation objects for optical detectors.
bool isVisible(double wavelength) const
Returns if we label as "visibile" a photon with specified wavelength [nm].
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
void analyze(art::Event const &)
static const int NoParticleId
Definition: sim.h:21
static constexpr double kVisibleWavelength
Value used when a typical visible light wavelength is needed.
void SetLibraryReflT0Entry(int VoxID, int OpChannel, float value)
T get(std::string const &key) const
Definition: ParameterSet.h:314
Class def header for mctrack data container.
virtual bool detectedLite(int OpChannel, int &newOpChannel) const
CryostatGeo const & Cryostat(CryostatID const &cryoid=details::cryostat_zero) const
Returns the specified cryostat.
const sim::ParticleList & ParticleList() const
std::size_t getView(std::string const &moduleLabel, std::string const &productInstanceName, std::string const &processName, std::vector< ELEMENT const * > &result) const
Encapsulate the geometry of an optical detector.
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
Collection of photons which recorded on one channel.
Definition: SimPhotons.h:127
Handle< PROD > getHandle(SelectorBase const &) const
void SetLibraryEntry(int VoxID, OpDetID_t libOpChannel, float N, bool wantReflected=false)
unsigned int NOpDet() const
Number of optical detectors in this TPC.
Definition: CryostatGeo.h:317
Float_t sc
Definition: plot.C:23
std::vector< std::vector< std::vector< double > > > fSignals_vis
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
object containing MC truth information necessary for making RawDigits and doing back tracking ...
IDNumber_t< Level::Event > EventNumber_t
Definition: IDNumber.h:118
EventNumber_t event() const
Definition: EventID.h:116
static constexpr double kVisibleThreshold
Threshold used to resolve between visible and ultraviolet light.
BoxBoundedGeo const & Boundaries() const
Returns boundaries of the cryostat (in centimetres).
Definition: CryostatGeo.h:113
std::vector< std::vector< double > > fSignalsvis
TCEvent evt
Definition: DataStructs.cxx:8
Particle list in DetSim contains Monte Carlo particle information.
Float_t e
Definition: plot.C:35
RunNumber_t run() const
Definition: Event.cc:29
virtual float wavelength(double energy) const
ROOT libraries.
std::vector< std::vector< std::vector< double > > > fSignals_vuv
Tools and modules for checking out the basics of the Monte Carlo.
EventID id() const
Definition: Event.cc:23
art framework interface to geometry description
OpDetGeo const & OpDetGeoFromOpDet(unsigned int OpDet) const
Returns the geo::OpDetGeo object for the given detector number.
Event finding and building.
std::vector< Handle< PROD > > getMany(SelectorBase const &selector=MatchAllSelector{}) const
cheat::ParticleInventoryService const * pi_serv