LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
OpFastScintillation.cxx
Go to the documentation of this file.
1 // Class adapted for LArSoft by Ben Jones, MIT 10/10/12
2 //
3 // This class is a physics process based on the standard Geant4
4 // scintillation process.
5 //
6 // It has been stripped down and adapted to form the backbone of
7 // the LArG4 fast optical simulation. Photons, instead of being
8 // produced and added to the geant4 particle stack, are logged
9 // and used to predict the visibility of this step to each PMT in
10 // the detector.
11 //
12 // The photonvisibilityservice looks up the visibility of the relevant
13 // xyz point, and if a photon is detected at a given PMT, one OnePhoton
14 // object is logged in the OpDetPhotonTable
15 //
16 // At the end of the event, the OpDetPhotonTable is read out
17 // by LArG4, and detected photons are stored in the event.
18 //
19 // This process can be used alongside the standard Cerenkov process,
20 // which does step geant4 opticalphotons. Both the fast scintillation
21 // table and the geant4 sensitive detectors are read out by LArG4 to
22 // produce a combined SimPhoton collection.
23 //
24 // Added disclaimer : This code is gross. Thats basically because
25 // it adheres to the original, gross Geant4 implementation.
26 //
27 //
28 // ********************************************************************
29 // * License and Disclaimer *
30 // * *
31 // * The Geant4 software is copyright of the Copyright Holders of *
32 // * the Geant4 Collaboration. It is provided under the terms and *
33 // * conditions of the Geant4 Software License, included in the file *
34 // * LICENSE and available at http://cern.ch/geant4/license . These *
35 // * include a list of copyright holders. *
36 // * *
37 // * Neither the authors of this software system, nor their employing *
38 // * institutes,nor the agencies providing financial support for this *
39 // * work make any representation or warranty, express or implied, *
40 // * regarding this software system or assume any liability for its *
41 // * use. Please see the license in the file LICENSE and URL above *
42 // * for the full disclaimer and the limitation of liability. *
43 // * *
44 // * This code implementation is the result of the scientific and *
45 // * technical work of the GEANT4 collaboration. *
46 // * By using, copying, modifying or distributing the software (or *
47 // * any work based on the software) you agree to acknowledge its *
48 // * use in resulting scientific publications, and indicate your *
49 // * acceptance of all terms of the Geant4 Software license. *
50 // ********************************************************************
51 //
52 //
53 // GEANT4 tag $Name: not supported by cvs2svn $
54 //
56 // Scintillation Light Class Implementation
58 //
59 // File: OpFastScintillation.cc
60 // Description: RestDiscrete Process - Generation of Scintillation Photons
61 // Version: 1.0
62 // Created: 1998-11-07
63 // Author: Peter Gumplinger
64 // Updated: 2010-10-20 Allow the scintillation yield to be a function
65 // of energy deposited by particle type
66 // Thanks to Zach Hartwig (Department of Nuclear
67 // Science and Engineeering - MIT)
68 // 2010-09-22 by Peter Gumplinger
69 // > scintillation rise time included, thanks to
70 // > Martin Goettlich/DESY
71 // 2005-08-17 by Peter Gumplinger
72 // > change variable name MeanNumPhotons -> MeanNumberOfPhotons
73 // 2005-07-28 by Peter Gumplinger
74 // > add G4ProcessType to constructor
75 // 2004-08-05 by Peter Gumplinger
76 // > changed StronglyForced back to Forced in GetMeanLifeTime
77 // 2002-11-21 by Peter Gumplinger
78 // > change to use G4Poisson for small MeanNumberOfPhotons
79 // 2002-11-07 by Peter Gumplinger
80 // > now allow for fast and slow scintillation component
81 // 2002-11-05 by Peter Gumplinger
82 // > now use scintillation constants from G4Material
83 // 2002-05-09 by Peter Gumplinger
84 // > use only the PostStepPoint location for the origin of
85 // scintillation photons when energy is lost to the medium
86 // by a neutral particle
87 // 2000-09-18 by Peter Gumplinger
88 // > change: aSecondaryPosition=x0+rand*aStep.GetDeltaPosition();
89 // aSecondaryTrack->SetTouchable(0);
90 // 2001-09-17, migration of Materials to pure STL (mma)
91 // 2003-06-03, V.Ivanchenko fix compilation warnings
92 //
93 // mail: gum@triumf.ca
94 //
96 
97 #include "Geant4/G4Alpha.hh"
98 #include "Geant4/G4DynamicParticle.hh"
99 #include "Geant4/G4Electron.hh"
100 #include "Geant4/G4ExceptionSeverity.hh"
101 #include "Geant4/G4Gamma.hh"
102 #include "Geant4/G4KaonMinus.hh"
103 #include "Geant4/G4KaonPlus.hh"
104 #include "Geant4/G4Material.hh"
105 #include "Geant4/G4MaterialPropertiesTable.hh"
106 #include "Geant4/G4MaterialPropertyVector.hh"
107 #include "Geant4/G4MaterialTable.hh"
108 #include "Geant4/G4MuonMinus.hh"
109 #include "Geant4/G4MuonPlus.hh"
110 #include "Geant4/G4ParticleChange.hh"
111 #include "Geant4/G4PhysicsVector.hh"
112 #include "Geant4/G4PionMinus.hh"
113 #include "Geant4/G4PionPlus.hh"
114 #include "Geant4/G4Poisson.hh"
115 #include "Geant4/G4Proton.hh"
116 #include "Geant4/G4Step.hh"
117 #include "Geant4/G4StepPoint.hh"
118 #include "Geant4/G4Track.hh"
119 #include "Geant4/G4VPhysicalVolume.hh"
120 #include "Geant4/G4ios.hh"
121 #include "Geant4/Randomize.hh"
122 #include "Geant4/globals.hh"
123 
129 #include "larcorealg/Geometry/geo_vectors_utils.h" // geo::vect::fillCoords()
130 #include "larcorealg/Geometry/geo_vectors_utils_TVector.h" // geo::vect::toTVector3()
140 
141 // support libraries
142 #include "cetlib_except/exception.h"
144 
145 #include "TGeoSphere.h"
146 #include "TMath.h"
147 #include "TRandom3.h"
148 
149 #include <cassert>
150 #include <cmath>
151 #include <limits>
152 
153 #include "boost/math/special_functions/ellint_1.hpp"
154 #include "boost/math/special_functions/ellint_3.hpp"
155 
156 // Define a new policy *not* internally promoting RealType to double:
157 typedef boost::math::policies::policy<
158  // boost::math::policies::digits10<8>,
159  boost::math::policies::promote_double<false>>
161 
162 namespace larg4 {
163 
165  // Class Implementation
167 
169  // Operators
171 
172  // OpFastScintillation::operator=(const OpFastScintillation &right)
173  // {
174  // }
175 
177  // Constructors
179 
180  OpFastScintillation::OpFastScintillation(const G4String& processName, G4ProcessType type)
181  : G4VRestDiscreteProcess(processName, type)
182  , fActiveVolumes{extractActiveVolumes(*(lar::providerFrom<geo::Geometry>()))}
183  , bPropagate(!(art::ServiceHandle<sim::LArG4Parameters const>()->NoPhotonPropagation()))
186  // for now, limit to the active volume only if semi-analytic model is used
188  {
189  SetProcessSubType(25); // TODO: unhardcode
190  fTrackSecondariesFirst = false;
191  fFiniteRiseTime = false;
192  YieldFactor = 1.0;
193  ExcitationRatio = 1.0;
194 
195  const detinfo::LArProperties* larp = lar::providerFrom<detinfo::LArPropertiesService>();
196 
198 
199  if (verboseLevel > 0) { G4cout << GetProcessName() << " is created " << G4endl; }
200 
202  emSaturation = NULL;
203 
204  if (bPropagate) {
205  assert(fPVS);
206 
207  // Loading the position of each optical channel, neccessary for the parametrizatiuons of Nhits and prop-time
208  geo::GeometryCore const& geom = *(lar::providerFrom<geo::Geometry>());
209 
210  {
211  auto log = mf::LogTrace("OpFastScintillation")
212  << "OpFastScintillation: active volume boundaries from " << fActiveVolumes.size()
213  << " volumes:";
214  for (auto const& [iCryo, box] : util::enumerate(fActiveVolumes)) {
215  log << "\n - C:" << iCryo << ": " << box.Min() << " -- " << box.Max() << " cm";
216  } // for
217  log << "\n (scintillation photons are propagated "
218  << (fOnlyActiveVolume ? "only from active volumes" : "from anywhere") << ")";
219  } // local scope
220 
221  if (usesSemiAnalyticModel() && (geom.Ncryostats() > 1U)) {
222  if (fOnlyOneCryostat) {
223  mf::LogWarning("OpFastScintillation")
224  << std::string(80, '=') << "\nA detector with " << geom.Ncryostats()
225  << " cryostats is configured"
226  << " , and semi-analytic model is requested for scintillation photon propagation."
227  << " THIS CONFIGURATION IS NOT SUPPORTED and it is open to bugs"
228  << " (e.g. scintillation may be detected only in cryostat #0)."
229  << "\nThis would be normally a fatal error, but it has been forcibly overridden."
230  << "\n"
231  << std::string(80, '=');
232  }
233  else {
235  << "Photon propagation via semi-analytic model is not supported yet"
236  << " on detectors with more than one cryostat.";
237  }
238  } // if
239 
240  geo::Point_t const Cathode_centre{geom.TPC().GetCathodeCenter().X(),
241  fActiveVolumes[0].CenterY(),
242  fActiveVolumes[0].CenterZ()};
243  mf::LogTrace("OpFastScintillation") << "Cathode_centre: " << Cathode_centre << " cm";
244 
245  // std::cout << "\nInitialize acos_arr with " << acos_bins+1
246  // << " hence with a resolution of " << 1./acos_bins << std::endl;
247  // for(size_t i=0; i<=acos_bins; ++i){
248  // acos_arr[i] = std::acos(i/double(acos_bins));
249  // }
250 
251  for (size_t const i : util::counter(fPVS->NOpChannels())) {
252  geo::OpDetGeo const& opDet = geom.OpDetGeoFromOpDet(i);
253  fOpDetCenter.push_back(opDet.GetCenter());
254 
255  if (dynamic_cast<TGeoSphere const*>(opDet.Shape()) != nullptr) { // sphere/dome
256  fOpDetType.push_back(1); // Dome PMTs
257  fOpDetLength.push_back(-1);
258  fOpDetHeight.push_back(-1);
259  }
260  else if (opDet.isBar()) { // box
261  fOpDetType.push_back(0); //Arapucas
262  fOpDetLength.push_back(opDet.Length());
263  fOpDetHeight.push_back(opDet.Height());
264  }
265  else { // disk
266  fOpDetType.push_back(2); // Disk PMTs
267  // std::cout<<"Radio: "<<geom.OpDetGeoFromOpDet(i).RMax()<<std::endl;
268  fOpDetLength.push_back(-1);
269  fOpDetHeight.push_back(-1);
270  }
271  }
272 
273  if (fPVS->IncludePropTime()) {
274  std::cout << "Using parameterisation of timings." << std::endl;
275  //OLD VUV time parapetrization (to be removed soon)
276  //fPVS->SetDirectLightPropFunctions(functions_vuv, fd_break, fd_max, ftf1_sampling_factor);
277  //fPVS->SetReflectedCOLightPropFunctions(functions_vis, ft0_max, ft0_break_point);
278  //New VUV time parapetrization
280  fstep_size,
281  fmax_d,
282  fmin_d,
287 
288  // create vector of empty TF1s that will be replaces with the parameterisations that are generated as they are required
289  // default TF1() constructor gives function with 0 dimensions, can then check numDim to qucikly see if a parameterisation has been generated
290  const size_t num_params =
291  (fmax_d - fmin_d) /
292  fstep_size; // for d < fmin_d, no parameterisaton, a delta function is used instead
293  size_t num_angles = std::round(90 / fangle_bin_timing_vuv);
294  VUV_timing = std::vector(num_angles, std::vector(num_params, TF1()));
295 
296  // initialise vectors to contain range parameterisations sampled to in each case
297  // when using TF1->GetRandom(xmin,xmax), must be in same range otherwise sampling is regenerated, this is the slow part!
298  VUV_max = std::vector(num_angles, std::vector(num_params, 0.0));
299  VUV_min = std::vector(num_angles, std::vector(num_params, 0.0));
300 
301  // VIS time parameterisation
302  if (fPVS->StoreReflected()) {
303  // load parameters
307  ftau_pars,
308  fvis_vmean,
310  }
311  }
312  if (usesSemiAnalyticModel()) {
313  mf::LogVerbatim("OpFastScintillation")
314  << "OpFastScintillation: using semi-analytic model for number of hits";
315 
316  // LAr absorption length in cm
317  std::map<double, double> abs_length_spectrum =
318  lar::providerFrom<detinfo::LArPropertiesService>()->AbsLengthSpectrum();
319  std::vector<double> x_v, y_v;
320  for (auto elem : abs_length_spectrum) {
321  x_v.push_back(elem.first);
322  y_v.push_back(elem.second);
323  }
324  fL_abs_vuv =
325  std::round(interpolate(x_v, y_v, 9.7, false)); // 9.7 eV: peak of VUV emission spectrum
326 
327  // Load Gaisser-Hillas corrections for VUV semi-analytic hits
328  std::cout << "Loading the GH corrections" << std::endl;
331  if (!fIsFlatPDCorr && !fIsDomePDCorr) {
332  throw cet::exception("OpFastScintillation")
333  << "Both isFlatPDCorr and isDomePDCorr parameters are false, at least one type of "
334  "parameterisation is required for the semi-analytic light simulation."
335  << "\n";
336  }
337  if (fIsFlatPDCorr) {
339  }
340  if (fIsDomePDCorr) {
342  }
343  // cathode center coordinates required for corrections
345  fcathode_centre.SetY(fActiveVolumes[0].CenterY());
346  fcathode_centre.SetZ(
347  fActiveVolumes[0].CenterZ()); // to get full cathode dimension rather than just single tpc
348 
349  if (fPVS->StoreReflected()) {
350  fStoreReflected = true;
351  // Load corrections for VIS semi-anlytic hits
352  std::cout << "Loading visible light corrections" << std::endl;
354  if (fIsFlatPDCorr) {
356  }
357  if (fIsDomePDCorr) {
359  }
360 
361  // cathode dimensions
362  fcathode_ydimension = fActiveVolumes[0].SizeY();
363  fcathode_zdimension = fActiveVolumes[0].SizeZ();
364 
365  // set cathode plane struct for solid angle function
369  }
370  else
371  fStoreReflected = false;
372  }
373  }
374  tpbemission = lar::providerFrom<detinfo::LArPropertiesService>()->TpbEm();
375  std::vector<double> parent;
376  parent.reserve(tpbemission.size());
377  for (auto iter = tpbemission.begin(); iter != tpbemission.end(); ++iter) {
378  parent.push_back(iter->second);
379  }
380  fTPBEm = std::make_unique<CLHEP::RandGeneral>(parent.data(), parent.size());
381  }
382 
384  // Destructors
387  {
388  if (theFastIntegralTable) theFastIntegralTable->clearAndDestroy();
389  if (theSlowIntegralTable) theSlowIntegralTable->clearAndDestroy();
390  }
391 
393  // Methods
395 
396  // AtRestDoIt
397  // ----------
398  //
399  G4VParticleChange* OpFastScintillation::AtRestDoIt(const G4Track& aTrack, const G4Step& aStep)
400  // This routine simply calls the equivalent PostStepDoIt since all the
401  // necessary information resides in aStep.GetTotalEnergyDeposit()
402  {
403  return OpFastScintillation::PostStepDoIt(aTrack, aStep);
404  }
405 
406  // PostStepDoIt
407  // -------------
408  //
409  G4VParticleChange* OpFastScintillation::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
410  // This routine is called for each tracking step of a charged particle
411  // in a scintillator. A Poisson/Gauss-distributed number of photons is
412  // generated according to the scintillation yield formula, distributed
413  // evenly along the track segment and uniformly into 4pi.
414  {
415  aParticleChange.Initialize(aTrack);
416  // Check that we are in a material with a properties table, if not
417  // just return
418  const G4Material* aMaterial = aTrack.GetMaterial();
419  G4MaterialPropertiesTable* aMaterialPropertiesTable = aMaterial->GetMaterialPropertiesTable();
420  if (!aMaterialPropertiesTable) return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
421 
422  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
423  G4ThreeVector x0 = pPreStepPoint->GetPosition();
424  G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
425 
427  // This is the old G4 way - but we do things differently - Ben J, Oct Nov 2012.
429  //
430  // if (MeanNumberOfPhotons > 10.)
431  // {
432  // G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons);
433  // NumPhotons = G4int(G4RandGauss::shoot(MeanNumberOfPhotons,sigma)+0.5);
434  // }
435  // else
436  // {
437  // NumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
438  // }
439  //
440  //
441  //
442  // if (NumPhotons <= 0)
443  // {
444  // // return unchanged particle and no secondaries
445  //
446  // aParticleChange.SetNumberOfSecondaries(0);
447  //
448  // return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
449  // }
450  //
451  //
452  // aParticleChange.SetNumberOfSecondaries(NumPhotons);
453  //
454  //
455  // if (fTrackSecondariesFirst) {
456  // if (aTrack.GetTrackStatus() == fAlive )
457  // aParticleChange.ProposeTrackStatus(fSuspend);
458  // }
459  //
460  //
461  //
462  //
464  //
465 
467  // The fast sim way - Ben J, Nov 2012
469  //
470  //
471  // We don't want to produce any trackable G4 secondaries
472  aParticleChange.SetNumberOfSecondaries(0);
473 
474  // Retrieve the Scintillation Integral for this material
475  // new G4PhysicsOrderedFreeVector allocated to hold CII's
476 
477  // Some explanation for later improvements to scint yield code:
478  //
479  // What does G4 do here?
480  // It produces light in 2 steps, fast (scnt=1) then slow (scnt=2)
481  //
482  // The ratio of slow photons to fast photons is related by the yieldratio
483  // parameter. G4's poisson fluctuating scheme is a bit different to ours
484  // - we should check that they are equivalent.
485  //
486  // G4 poisson fluctuates the number of initial photons then divides them
487  // with a constant factor between fast + slow, whereas we poisson
488  // fluctuate separateyly the fast and slow detection numbers.
489  //
490 
491  // get the number of photons produced from the IonizationAndScintillation
492  // singleton
494  double MeanNumberOfPhotons =
496  // double stepEnergy = larg4::IonizationAndScintillation::Instance()->VisibleEnergyDeposit()/CLHEP::MeV;
497  RecordPhotonsProduced(aStep, MeanNumberOfPhotons); //, stepEnergy);
498  if (verboseLevel > 0) {
499  G4cout << "\n Exiting from OpFastScintillation::DoIt -- NumberOfSecondaries = "
500  << aParticleChange.GetNumberOfSecondaries() << G4endl;
501  }
502  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
503  }
504 
505  void OpFastScintillation::ProcessStep(const G4Step& step)
506  {
507  if (step.GetTotalEnergyDeposit() <= 0) return;
508 
510  -1,
511  -1,
512  1.0, //scintillation yield
513  (double)(step.GetTotalEnergyDeposit() / CLHEP::MeV), //energy in MeV
514  (float)(step.GetPreStepPoint()->GetPosition().x() / CLHEP::cm),
515  (float)(step.GetPreStepPoint()->GetPosition().y() / CLHEP::cm),
516  (float)(step.GetPreStepPoint()->GetPosition().z() / CLHEP::cm),
517  (float)(step.GetPostStepPoint()->GetPosition().x() / CLHEP::cm),
518  (float)(step.GetPostStepPoint()->GetPosition().y() / CLHEP::cm),
519  (float)(step.GetPostStepPoint()->GetPosition().z() / CLHEP::cm),
520  (double)(step.GetPreStepPoint()->GetGlobalTime()),
521  (double)(step.GetPostStepPoint()->GetGlobalTime()),
522  //step.GetTrack()->GetTrackID(),
524  step.GetTrack()->GetParticleDefinition()->GetPDGEncoding(),
526  step.GetPreStepPoint()->GetPhysicalVolume()->GetName());
527  }
528 
530  double MeanNumberOfPhotons) //, double stepEnergy)
531  {
532  // make sure that whatever happens afterwards, the energy deposition is stored
534  if (lgp->FillSimEnergyDeposits()) ProcessStep(aStep);
535 
536  // Get the pointer to the fast scintillation table
538 
539  const G4Track* aTrack = aStep.GetTrack();
540 
541  G4StepPoint const* pPreStepPoint = aStep.GetPreStepPoint();
542  // unused G4StepPoint const* pPostStepPoint = aStep.GetPostStepPoint();
543 
544  const G4DynamicParticle* aParticle = aTrack->GetDynamicParticle();
545  const G4Material* aMaterial = aTrack->GetMaterial();
546 
547  G4int materialIndex = aMaterial->GetIndex();
548  G4int tracknumber = aStep.GetTrack()->GetTrackID();
549 
550  G4ThreeVector x0 = pPreStepPoint->GetPosition();
551  G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
552  //G4double t0 = pPreStepPoint->GetGlobalTime() - fGlobalTimeOffset;
553  G4double t0 = pPreStepPoint->GetGlobalTime();
554 
555  G4MaterialPropertiesTable* aMaterialPropertiesTable = aMaterial->GetMaterialPropertiesTable();
556 
557  G4MaterialPropertyVector* Fast_Intensity =
558  aMaterialPropertiesTable->GetProperty("FASTCOMPONENT");
559  G4MaterialPropertyVector* Slow_Intensity =
560  aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT");
561 
562  if (!Fast_Intensity && !Slow_Intensity) return 1;
563 
564  if (!bPropagate) return 0;
565 
566  // Get the visibility vector for this point
567  assert(fPVS);
568 
569  G4int nscnt = 1;
570  if (Fast_Intensity && Slow_Intensity) nscnt = 2;
571 
572  double Num = 0;
573  double YieldRatio = 0;
574 
576  // The scintillation response is a function of the energy
577  // deposited by particle types.
578 
579  // Get the definition of the current particle
580  G4ParticleDefinition* pDef = aParticle->GetDefinition();
581 
582  // Obtain the G4MaterialPropertyVectory containing the
583  // scintillation light yield as a function of the deposited
584  // energy for the current particle type
585 
586  // Protons
587  if (pDef == G4Proton::ProtonDefinition()) {
588  YieldRatio = aMaterialPropertiesTable->GetConstProperty("PROTONYIELDRATIO");
589  }
590  // Muons
591  else if (pDef == G4MuonPlus::MuonPlusDefinition() ||
592  pDef == G4MuonMinus::MuonMinusDefinition()) {
593  YieldRatio = aMaterialPropertiesTable->GetConstProperty("MUONYIELDRATIO");
594  }
595  // Pions
596  else if (pDef == G4PionPlus::PionPlusDefinition() ||
597  pDef == G4PionMinus::PionMinusDefinition()) {
598  YieldRatio = aMaterialPropertiesTable->GetConstProperty("PIONYIELDRATIO");
599  }
600  // Kaons
601  else if (pDef == G4KaonPlus::KaonPlusDefinition() ||
602  pDef == G4KaonMinus::KaonMinusDefinition()) {
603  YieldRatio = aMaterialPropertiesTable->GetConstProperty("KAONYIELDRATIO");
604  }
605  // Alphas
606  else if (pDef == G4Alpha::AlphaDefinition()) {
607  YieldRatio = aMaterialPropertiesTable->GetConstProperty("ALPHAYIELDRATIO");
608  }
609  // Electrons (must also account for shell-binding energy
610  // attributed to gamma from standard PhotoElectricEffect)
611  else if (pDef == G4Electron::ElectronDefinition() || pDef == G4Gamma::GammaDefinition()) {
612  YieldRatio = aMaterialPropertiesTable->GetConstProperty("ELECTRONYIELDRATIO");
613  }
614  // Default for particles not enumerated/listed above
615  else {
616  YieldRatio = aMaterialPropertiesTable->GetConstProperty("ELECTRONYIELDRATIO");
617  }
618  // If the user has not specified yields for (p,d,t,a,carbon)
619  // then these unspecified particles will default to the
620  // electron's scintillation yield
621  if (YieldRatio == 0) {
622  YieldRatio = aMaterialPropertiesTable->GetConstProperty("ELECTRONYIELDRATIO");
623  }
624  }
625 
626  geo::Point_t const ScintPoint = {x0[0] / CLHEP::cm, x0[1] / CLHEP::cm, x0[2] / CLHEP::cm};
627  if (fOnlyActiveVolume && !isScintInActiveVolume(ScintPoint)) return 0;
628  const phot::MappedCounts_t& Visibilities = fPVS->GetAllVisibilities(ScintPoint);
629 
630  phot::MappedCounts_t ReflVisibilities;
631 
632  // Store timing information in the object for use in propagationTime method
633  if (fPVS->StoreReflected()) {
634  ReflVisibilities = fPVS->GetAllVisibilities(ScintPoint, true);
635  if (fPVS->StoreReflT0()) ReflT0s = fPVS->GetReflT0s(ScintPoint);
636  }
637  if (fPVS->IncludeParPropTime()) { ParPropTimeTF1 = fPVS->GetTimingTF1(ScintPoint); }
638 
639  /*
640  // For Kazu to debug # photons generated using csv file, by default should be commented out
641  // but don't remove as it's useful. Separated portion of code relevant to this block
642  // is labeled as "CASE-DEBUG DO NOT REMOVE THIS COMMENT"
643  // (2017 May 2nd ... "kazuhiro at nevis dot columbia dot edu")
644  //
645  static bool first_time=false;
646  if(!first_time) {
647  std::cout<<"LOGMEid,pdg,mass,ke,te,de,x,y,z,t,dr,mean_npe,gen_npe,det_npe"<<std::endl;
648  first_time=true;
649  }
650 
651  std::cout<<"LOGME"
652  <<aStep.GetTrack()->GetTrackID()<<","
653  <<aParticle->GetDefinition()->GetPDGEncoding()<<","
654  <<aParticle->GetDefinition()->GetPDGMass()/CLHEP::MeV<<","
655  <<pPreStepPoint->GetKineticEnergy()<<","
656  <<pPreStepPoint->GetTotalEnergy()<<","
657  <<aStep.GetTotalEnergyDeposit()<<","
658  <<ScintPoint<<","<<t0<<","
659  <<aStep.GetDeltaPosition().mag()<<","
660  <<MeanNumberOfPhotons<<","<<std::flush;
661 
662  double gen_photon_ctr=0;
663  double det_photon_ctr=0;
664  */
665  for (G4int scnt = 1; scnt <= nscnt; scnt++) {
666  G4double ScintillationTime = 0. * CLHEP::ns;
667  G4double ScintillationRiseTime = 0. * CLHEP::ns;
668  G4PhysicsOrderedFreeVector* ScintillationIntegral = NULL;
669  if (scnt == 1) {
670  if (nscnt == 1) {
671  if (Fast_Intensity) {
672  ScintillationTime = aMaterialPropertiesTable->GetConstProperty("FASTTIMECONSTANT");
673  if (fFiniteRiseTime) {
674  ScintillationRiseTime =
675  aMaterialPropertiesTable->GetConstProperty("FASTSCINTILLATIONRISETIME");
676  }
677  ScintillationIntegral =
678  (G4PhysicsOrderedFreeVector*)((*theFastIntegralTable)(materialIndex));
679  }
680  if (Slow_Intensity) {
681  ScintillationTime = aMaterialPropertiesTable->GetConstProperty("SLOWTIMECONSTANT");
682  if (fFiniteRiseTime) {
683  ScintillationRiseTime =
684  aMaterialPropertiesTable->GetConstProperty("SLOWSCINTILLATIONRISETIME");
685  }
686  ScintillationIntegral =
687  (G4PhysicsOrderedFreeVector*)((*theSlowIntegralTable)(materialIndex));
688  }
689  } //endif nscnt=1
690  else {
691  if (YieldRatio == 0)
692  YieldRatio = aMaterialPropertiesTable->GetConstProperty("YIELDRATIO");
693 
694  if (ExcitationRatio == 1.0) { Num = std::min(YieldRatio, 1.0) * MeanNumberOfPhotons; }
695  else {
696  Num = std::min(ExcitationRatio, 1.0) * MeanNumberOfPhotons;
697  }
698  ScintillationTime = aMaterialPropertiesTable->GetConstProperty("FASTTIMECONSTANT");
699  if (fFiniteRiseTime) {
700  ScintillationRiseTime =
701  aMaterialPropertiesTable->GetConstProperty("FASTSCINTILLATIONRISETIME");
702  }
703  ScintillationIntegral =
704  (G4PhysicsOrderedFreeVector*)((*theFastIntegralTable)(materialIndex));
705  } //endif nscnt!=1
706  } //end scnt=1
707 
708  else {
709  Num = MeanNumberOfPhotons - Num;
710  ScintillationTime = aMaterialPropertiesTable->GetConstProperty("SLOWTIMECONSTANT");
711  if (fFiniteRiseTime) {
712  ScintillationRiseTime =
713  aMaterialPropertiesTable->GetConstProperty("SLOWSCINTILLATIONRISETIME");
714  }
715  ScintillationIntegral =
716  (G4PhysicsOrderedFreeVector*)((*theSlowIntegralTable)(materialIndex));
717  }
718 
719  if (!ScintillationIntegral) continue;
720  //gen_photon_ctr += Num; // CASE-DEBUG DO NOT REMOVE THIS COMMENT
721  // Max Scintillation Integral
722  // G4double CIImax = ScintillationIntegral->GetMaxValue();
723  //std::cout << "++++++++++++" << Num << "++++++++++" << std::endl;
724 
725  // here we go: now if visibilities are invalid, we are in trouble
726  //if (!Visibilities && (fPVS->NOpChannels() > 0)) {
727  // throw cet::exception("OpFastScintillator")
728  // << "Photon library does not cover point " << ScintPoint << " cm.\n";
729  //}
730 
731  if (!Visibilities && !usesSemiAnalyticModel()) continue;
732 
733  // detected photons from direct light
734  std::map<size_t, int> DetectedNum;
735  if (Visibilities && !usesSemiAnalyticModel()) {
736  for (size_t const OpDet : util::counter(fPVS->NOpChannels())) {
737  if (fOpaqueCathode && !isOpDetInSameTPC(ScintPoint, fOpDetCenter.at(OpDet))) continue;
738  int const DetThis = std::round(G4Poisson(Visibilities[OpDet] * Num));
739  if (DetThis > 0) DetectedNum[OpDet] = DetThis;
740  }
741  }
742  else {
743  detectedDirectHits(DetectedNum, Num, ScintPoint);
744  }
745 
746  // detected photons from reflected light
747  std::map<size_t, int> ReflDetectedNum;
748  if (fPVS->StoreReflected()) {
749  if (!usesSemiAnalyticModel()) {
750  for (size_t const OpDet : util::counter(fPVS->NOpChannels())) {
751  if (fOpaqueCathode && !isOpDetInSameTPC(ScintPoint, fOpDetCenter.at(OpDet))) continue;
752  int const ReflDetThis = std::round(G4Poisson(ReflVisibilities[OpDet] * Num));
753  if (ReflDetThis > 0) ReflDetectedNum[OpDet] = ReflDetThis;
754  }
755  }
756  else {
757  detectedReflecHits(ReflDetectedNum, Num, ScintPoint);
758  }
759  }
760 
761  std::vector<double> arrival_time_dist;
762  // Now we run through each PMT figuring out num of detected photons
763  for (size_t Reflected = 0; Reflected <= 1; ++Reflected) {
764  // Only do the reflected loop if we have reflected visibilities
765  if (Reflected && !fPVS->StoreReflected()) continue;
766 
769  if (Reflected) {
770  itstart = ReflDetectedNum.begin();
771  itend = ReflDetectedNum.end();
772  }
773  else {
774  itstart = DetectedNum.begin();
775  itend = DetectedNum.end();
776  }
777  for (auto itdetphot = itstart; itdetphot != itend; ++itdetphot) {
778  const size_t OpChannel = itdetphot->first;
779  const int NPhotons = itdetphot->second;
780 
781  // Set up the OpDetBTR information
782  sim::OpDetBacktrackerRecord tmpOpDetBTRecord(OpChannel);
783  int thisG4TrackID = ParticleListAction::GetCurrentTrackID();
784  double xyzPos[3];
785  average_position(aStep, xyzPos);
786  double Edeposited = 0;
788  //We use this when it is the only sensical information. It may be of limited use to end users.
789  Edeposited = aStep.GetTotalEnergyDeposit();
790  }
791  else if (emSaturation) {
792  //If Birk Coefficient used, log VisibleEnergies.
793  Edeposited =
795  }
796  else {
797  //We use this when it is the only sensical information. It may be of limited use to end users.
798  Edeposited = aStep.GetTotalEnergyDeposit();
799  }
800 
801  // Get the transport time distribution
802  arrival_time_dist.resize(NPhotons);
803  propagationTime(arrival_time_dist, x0, OpChannel, Reflected);
804 
805  //We need to split the energy up by the number of photons so that we never try to write a 0 energy.
806  Edeposited = Edeposited / double(NPhotons);
807 
808  // Loop through the photons
809  for (G4int i = 0; i < NPhotons; ++i) {
810  //std::cout<<"VUV time correction: "<<arrival_time_dist[i]<<std::endl;
811  G4double Time = t0 + scint_time(aStep, ScintillationTime, ScintillationRiseTime) +
812  arrival_time_dist[i] * CLHEP::ns;
813 
814  // Always store the BTR
815  tmpOpDetBTRecord.AddScintillationPhotons(thisG4TrackID, Time, 1, xyzPos, Edeposited);
816 
817  // Store as lite photon or as OnePhoton
818  if (lgp->UseLitePhotons()) {
819  fst->AddLitePhoton(OpChannel, static_cast<int>(Time), 1, Reflected);
820  }
821  else {
822  // The sim photon in this case stores its production point and time
823  TVector3 PhotonPosition(x0[0], x0[1], x0[2]);
824 
825  float PhotonEnergy = 0;
826  if (Reflected)
827  PhotonEnergy = reemission_energy() * CLHEP::eV;
828  else
829  PhotonEnergy = 9.7 * CLHEP::eV; // 9.7 eV peak of VUV emission spectrum
830 
831  // Make a photon object for the collection
832  sim::OnePhoton PhotToAdd;
833  PhotToAdd.InitialPosition = PhotonPosition;
834  PhotToAdd.Energy = PhotonEnergy;
835  PhotToAdd.Time = Time;
836  PhotToAdd.SetInSD = false;
837  PhotToAdd.MotherTrackID = tracknumber;
838 
839  fst->AddPhoton(OpChannel, std::move(PhotToAdd), Reflected);
840  }
841  }
842  fst->AddOpDetBacktrackerRecord(tmpOpDetBTRecord, Reflected);
843  }
844  }
845  }
846  //std::cout<<gen_photon_ctr<<","<<det_photon_ctr<<std::endl; // CASE-DEBUG DO NOT REMOVE THIS COMMENT
847  return 0;
848  }
849 
850  // BuildThePhysicsTable for the scintillation process
851  // --------------------------------------------------
852  //
854  {
856 
857  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
858  G4int numOfMaterials = G4Material::GetNumberOfMaterials();
859 
860  // create new physics table
862  theFastIntegralTable = std::make_unique<G4PhysicsTable>(numOfMaterials);
864  theSlowIntegralTable = std::make_unique<G4PhysicsTable>(numOfMaterials);
865 
866  // loop for materials
867  for (G4int i = 0; i < numOfMaterials; i++) {
868  G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector = new G4PhysicsOrderedFreeVector();
869  G4PhysicsOrderedFreeVector* bPhysicsOrderedFreeVector = new G4PhysicsOrderedFreeVector();
870 
871  // Retrieve vector of scintillation wavelength intensity for
872  // the material from the material's optical properties table.
873  G4Material* aMaterial = (*theMaterialTable)[i];
874 
875  G4MaterialPropertiesTable* aMaterialPropertiesTable = aMaterial->GetMaterialPropertiesTable();
876 
877  if (aMaterialPropertiesTable) {
878 
879  G4MaterialPropertyVector* theFastLightVector =
880  aMaterialPropertiesTable->GetProperty("FASTCOMPONENT");
881 
882  if (theFastLightVector) {
883  // Retrieve the first intensity point in vector
884  // of (photon energy, intensity) pairs
885  G4double currentIN = (*theFastLightVector)[0];
886  if (currentIN >= 0.0) {
887  // Create first (photon energy, Scintillation
888  // Integral pair
889  G4double currentPM = theFastLightVector->Energy(0);
890  G4double currentCII = 0.0;
891 
892  aPhysicsOrderedFreeVector->InsertValues(currentPM, currentCII);
893 
894  // Set previous values to current ones prior to loop
895  G4double prevPM = currentPM;
896  G4double prevCII = currentCII;
897  G4double prevIN = currentIN;
898 
899  // loop over all (photon energy, intensity)
900  // pairs stored for this material
901  for (size_t i = 1; i < theFastLightVector->GetVectorLength(); i++) {
902  currentPM = theFastLightVector->Energy(i);
903  currentIN = (*theFastLightVector)[i];
904 
905  currentCII = 0.5 * (prevIN + currentIN);
906 
907  currentCII = prevCII + (currentPM - prevPM) * currentCII;
908 
909  aPhysicsOrderedFreeVector->InsertValues(currentPM, currentCII);
910 
911  prevPM = currentPM;
912  prevCII = currentCII;
913  prevIN = currentIN;
914  }
915  }
916  }
917 
918  G4MaterialPropertyVector* theSlowLightVector =
919  aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT");
920  if (theSlowLightVector) {
921  // Retrieve the first intensity point in vector
922  // of (photon energy, intensity) pairs
923  G4double currentIN = (*theSlowLightVector)[0];
924  if (currentIN >= 0.0) {
925  // Create first (photon energy, Scintillation
926  // Integral pair
927  G4double currentPM = theSlowLightVector->Energy(0);
928  G4double currentCII = 0.0;
929 
930  bPhysicsOrderedFreeVector->InsertValues(currentPM, currentCII);
931 
932  // Set previous values to current ones prior to loop
933  G4double prevPM = currentPM;
934  G4double prevCII = currentCII;
935  G4double prevIN = currentIN;
936 
937  // loop over all (photon energy, intensity)
938  // pairs stored for this material
939  for (size_t i = 1; i < theSlowLightVector->GetVectorLength(); i++) {
940  currentPM = theSlowLightVector->Energy(i);
941  currentIN = (*theSlowLightVector)[i];
942 
943  currentCII = 0.5 * (prevIN + currentIN);
944 
945  currentCII = prevCII + (currentPM - prevPM) * currentCII;
946 
947  bPhysicsOrderedFreeVector->InsertValues(currentPM, currentCII);
948 
949  prevPM = currentPM;
950  prevCII = currentCII;
951  prevIN = currentIN;
952  }
953  }
954  }
955  }
956  // The scintillation integral(s) for a given material
957  // will be inserted in the table(s) according to the
958  // position of the material in the material table.
959  theFastIntegralTable->insertAt(i, aPhysicsOrderedFreeVector);
960  theSlowIntegralTable->insertAt(i, bPhysicsOrderedFreeVector);
961  }
962  }
963 
964  // Called by the user to set the scintillation yield as a function
965  // of energy deposited by particle type
967  {
968  if (emSaturation) {
969  G4Exception("OpFastScintillation::SetScintillationByParticleType",
970  "Scint02",
971  JustWarning,
972  "Redefinition: Birks Saturation is replaced by ScintillationByParticleType!");
974  }
975  scintillationByParticleType = scintType;
976  }
977 
978  G4double OpFastScintillation::GetMeanFreePath(const G4Track&,
979  G4double,
980  G4ForceCondition* condition)
981  {
982  *condition = StronglyForced;
983  return DBL_MAX;
984  }
985 
986  G4double OpFastScintillation::GetMeanLifeTime(const G4Track&, G4ForceCondition* condition)
987  {
988  *condition = Forced;
989  return DBL_MAX;
990  }
991 
992  G4double OpFastScintillation::scint_time(const G4Step& aStep,
993  G4double ScintillationTime,
994  G4double ScintillationRiseTime) const
995  {
996  G4StepPoint const* pPreStepPoint = aStep.GetPreStepPoint();
997  G4StepPoint const* pPostStepPoint = aStep.GetPostStepPoint();
998  G4double avgVelocity = (pPreStepPoint->GetVelocity() + pPostStepPoint->GetVelocity()) / 2.;
999  G4double deltaTime = aStep.GetStepLength() / avgVelocity;
1000  if (ScintillationRiseTime == 0.0) {
1001  deltaTime = deltaTime - ScintillationTime * std::log(G4UniformRand());
1002  }
1003  else {
1004  deltaTime = deltaTime + sample_time(ScintillationRiseTime, ScintillationTime);
1005  }
1006  return deltaTime;
1007  }
1008 
1009  void OpFastScintillation::propagationTime(std::vector<double>& arrival_time_dist,
1010  G4ThreeVector x0,
1011  const size_t OpChannel,
1012  bool Reflected) //const
1013  {
1014  assert(fPVS);
1016  throw cet::exception("OpFastScintillation")
1017  << "Cannot have both propagation time models simultaneously.";
1018  }
1019  else if (fPVS->IncludeParPropTime() &&
1020  !(ParPropTimeTF1 && (ParPropTimeTF1[OpChannel].GetNdim() == 1))) {
1021  //Warning: TF1::GetNdim()==1 will tell us if the TF1 is really defined or it is the default one.
1022  //This will fix a segfault when using timing and interpolation.
1023  G4cout << "WARNING: Requested parameterized timing, but no function found. Not applying "
1024  "propagation time."
1025  << G4endl;
1026  }
1027  else if (fPVS->IncludeParPropTime()) {
1028  if (Reflected)
1029  throw cet::exception("OpFastScintillation")
1030  << "No parameterized propagation time for reflected light";
1031  for (size_t i = 0; i < arrival_time_dist.size(); ++i) {
1032  arrival_time_dist[i] = ParPropTimeTF1[OpChannel].GetRandom();
1033  }
1034  }
1035  else if (fPVS->IncludePropTime()) {
1036  // Get VUV photons arrival time distribution from the parametrization
1037  geo::Point_t const& opDetCenter = fOpDetCenter.at(OpChannel);
1038  if (!Reflected) {
1039  const G4ThreeVector OpDetPoint(
1040  opDetCenter.X() * CLHEP::cm, opDetCenter.Y() * CLHEP::cm, opDetCenter.Z() * CLHEP::cm);
1041  double distance_in_cm = (x0 - OpDetPoint).mag() / CLHEP::cm; // this must be in CENTIMETERS!
1042  double cosine = std::abs(x0[0] * CLHEP::cm - OpDetPoint[0] * CLHEP::cm) / distance_in_cm;
1043  double theta = fast_acos(cosine) * 180. / CLHEP::pi;
1044  int angle_bin = theta / fangle_bin_timing_vuv;
1045  getVUVTimes(arrival_time_dist, distance_in_cm, angle_bin); // in ns
1046  }
1047  else {
1048  TVector3 const ScintPoint(x0[0] / CLHEP::cm, x0[1] / CLHEP::cm, x0[2] / CLHEP::cm); // in cm
1049  getVISTimes(arrival_time_dist, ScintPoint, geo::vect::toTVector3(opDetCenter)); // in ns
1050  }
1051  }
1052  }
1053 
1054  G4double OpFastScintillation::sample_time(const G4double tau1, const G4double tau2) const
1055  {
1056  // tau1: rise time and tau2: decay time
1057  while (1) {
1058  // two random numbers
1059  G4double ran1 = G4UniformRand();
1060  G4double ran2 = G4UniformRand();
1061  //
1062  // exponential distribution as envelope function: very efficient
1063  //
1064  G4double d = (tau1 + tau2) / tau2;
1065  // make sure the envelope function is
1066  // always larger than the bi-exponential
1067  G4double t = -1.0 * tau2 * std::log(1 - ran1);
1068  G4double g = d * single_exp(t, tau2);
1069  if (ran2 <= bi_exp(t, tau1, tau2) / g) return t;
1070  }
1071  return -1.0;
1072  }
1073 
1075  {
1076  return fTPBEm->fire() * ((*(--tpbemission.end())).first - (*tpbemission.begin()).first) +
1077  (*tpbemission.begin()).first;
1078  }
1079 
1080  void OpFastScintillation::average_position(G4Step const& aStep, double* xyzPos) const
1081  {
1082  CLHEP::Hep3Vector prePoint = (aStep.GetPreStepPoint())->GetPosition();
1083  CLHEP::Hep3Vector postPoint = (aStep.GetPostStepPoint())->GetPosition();
1084  xyzPos[0] = (((prePoint.getX() + postPoint.getX()) / 2.0) / CLHEP::cm);
1085  xyzPos[1] = (((prePoint.getY() + postPoint.getY()) / 2.0) / CLHEP::cm);
1086  xyzPos[2] = (((prePoint.getZ() + postPoint.getZ()) / 2.0) / CLHEP::cm);
1087  }
1088 
1089  // ======TIMING PARAMETRIZATION===== //
1090  /*
1091  // Parametrization of the VUV light timing (result from direct transport + Rayleigh scattering ONLY)
1092  // using a landau + expo function.The function below returns the arrival time distribution given the
1093  // distance IN CENTIMETERS between the scintillation/ionization point and the optical detectotr.
1094  std::vector<double> OpFastScintillation::GetVUVTime(double distance, int number_photons) {
1095 
1096  //-----Distances in cm and times in ns-----//
1097  //gRandom->SetSeed(0);
1098  std::vector<double> arrival_time_distrb;
1099  arrival_time_distrb.clear();
1100 
1101  // Parametrization data:
1102  if(distance < 10 || distance > fd_max) {
1103  G4cout<<"WARNING: Parametrization of Direct Light not fully reliable"<<G4endl;
1104  G4cout<<"Too close/far to the PMT -> set 0 VUV photons(?)!!!!!!"<<G4endl;
1105  return arrival_time_distrb;
1106  }
1107  //signals (remember this is transportation) no longer than 1us
1108  const double signal_t_range = 1000.;
1109  const double vuv_vgroup = 10.13;//cm/ns
1110  double t_direct = distance/vuv_vgroup;
1111  // Defining the two functions (Landau + Exponential) describing the timing vs distance
1112  double pars_landau[3] = {functions_vuv[1]->Eval(distance), functions_vuv[2]->Eval(distance),
1113  std::pow(10.,functions_vuv[0]->Eval(distance))};
1114  if(distance > fd_break) {
1115  pars_landau[0]=functions_vuv[6]->Eval(distance);
1116  pars_landau[1]=functions_vuv[2]->Eval(fd_break);
1117  pars_landau[2]=std::pow(10.,functions_vuv[5]->Eval(distance));
1118  }
1119  TF1 *flandau = new TF1("flandau","[2]*TMath::Landau(x,[0],[1])",0,signal_t_range/2);
1120  flandau->SetParameters(pars_landau);
1121  double pars_expo[2] = {functions_vuv[3]->Eval(distance), functions_vuv[4]->Eval(distance)};
1122  if(distance > (fd_break - 50.)) {
1123  pars_expo[0] = functions_vuv[7]->Eval(distance);
1124  pars_expo[1] = functions_vuv[4]->Eval(fd_break - 50.);
1125  }
1126  TF1 *fexpo = new TF1("fexpo","expo",0, signal_t_range/2);
1127  fexpo->SetParameters(pars_expo);
1128  //this is to find the intersection point between the two functions:
1129  TF1 *fint = new TF1("fint","finter_d",flandau->GetMaximumX(),3*t_direct,5);
1130  double parsInt[5] = {pars_landau[0], pars_landau[1], pars_landau[2], pars_expo[0], pars_expo[1]};
1131  fint->SetParameters(parsInt);
1132  double t_int = fint->GetMinimumX();
1133  double minVal = fint->Eval(t_int);
1134  //the functions must intersect!!!
1135  if(minVal>0.015)
1136  G4cout<<"WARNING: Parametrization of Direct Light discontinuous (landau + expo)!!!!!!"<<G4endl;
1137 
1138  TF1 *fVUVTiming = new TF1("fTiming","LandauPlusExpoFinal",0,signal_t_range,6);
1139  double parsfinal[6] = {t_int, pars_landau[0], pars_landau[1], pars_landau[2], pars_expo[0], pars_expo[1]};
1140  fVUVTiming->SetParameters(parsfinal);
1141  // Set the number of points used to sample the function
1142 
1143  int f_sampling = 1000;
1144  if(distance < 50)
1145  f_sampling *= ftf1_sampling_factor;
1146  fVUVTiming->SetNpx(f_sampling);
1147 
1148  for(int i=0; i<number_photons; i++)
1149  arrival_time_distrb.push_back(fVUVTiming->GetRandom());
1150 
1151  //deleting ...
1152  delete flandau;
1153  delete fexpo;
1154  delete fint;
1155  delete fVUVTiming;
1156  //G4cout<<"BEAMAUS timing distribution hecha"<<G4endl;
1157  return arrival_time_distrb;
1158  }
1159 
1160  // Parametrization of the Visible light timing (result from direct transport + Rayleigh scattering ONLY)
1161  // using a landau + exponential function. The function below returns the arrival time distribution given the
1162  // time of the first visible photon in the PMT. The light generated has been reflected by the cathode ONLY.
1163  std::vector<double> OpFastScintillation::GetVisibleTimeOnlyCathode(double t0, int number_photons) {
1164  //-----Distances in cm and times in ns-----//
1165  //gRandom->SetSeed(0);
1166 
1167  std::vector<double> arrival_time_distrb;
1168  arrival_time_distrb.clear();
1169  // Parametrization data:
1170 
1171  if(t0 < 8 || t0 > ft0_max) {
1172  G4cout<<"WARNING: Parametrization of Cathode-Only reflected Light not fully reliable"<<G4endl;
1173  G4cout<<"Too close/far to the PMT -> set 0 Visible photons(?)!!!!!!"<<G4endl;
1174  return arrival_time_distrb;
1175  }
1176  //signals (remember this is transportation) no longer than 1us
1177  const double signal_t_range = 1000.;
1178  double pars_landau[3] = {functions_vis[1]->Eval(t0), functions_vis[2]->Eval(t0),
1179  std::pow(10.,functions_vis[0]->Eval(t0))};
1180  double pars_expo[2] = {functions_vis[3]->Eval(t0), functions_vis[4]->Eval(t0)};
1181  if(t0 > ft0_break_point) {
1182  pars_landau[0] = -0.798934 + 1.06216*t0;
1183  pars_landau[1] = functions_vis[2]->Eval(ft0_break_point);
1184  pars_landau[2] = std::pow(10.,functions_vis[0]->Eval(ft0_break_point));
1185  pars_expo[0] = functions_vis[3]->Eval(ft0_break_point);
1186  pars_expo[1] = functions_vis[4]->Eval(ft0_break_point);
1187  }
1188 
1189  // Defining the two functions (Landau + Exponential) describing the timing vs t0
1190  TF1 *flandau = new TF1("flandau","[2]*TMath::Landau(x,[0],[1])",0,signal_t_range/2);
1191  flandau->SetParameters(pars_landau);
1192  TF1 *fexpo = new TF1("fexpo","expo",0, signal_t_range/2);
1193  fexpo->SetParameters(pars_expo);
1194  //this is to find the intersection point between the two functions:
1195  TF1 *fint = new TF1("fint",finter_d,flandau->GetMaximumX(),2*t0,5);
1196  double parsInt[5] = {pars_landau[0], pars_landau[1], pars_landau[2], pars_expo[0], pars_expo[1]};
1197  fint->SetParameters(parsInt);
1198  double t_int = fint->GetMinimumX();
1199  double minVal = fint->Eval(t_int);
1200  //the functions must intersect!!!
1201  if(minVal>0.015)
1202  G4cout<<"WARNING: Parametrization of Direct Light discontinuous (landau + expo)!!!!!!"<<G4endl;
1203 
1204  TF1 *fVisTiming = new TF1("fTiming",LandauPlusExpoFinal,0,signal_t_range,6);
1205  double parsfinal[6] = {t_int, pars_landau[0], pars_landau[1], pars_landau[2], pars_expo[0], pars_expo[1]};
1206  fVisTiming->SetParameters(parsfinal);
1207  // Set the number of points used to sample the function
1208 
1209  int f_sampling = 1000;
1210  if(t0 < 20)
1211  f_sampling *= ftf1_sampling_factor;
1212  fVisTiming->SetNpx(f_sampling);
1213 
1214  for(int i=0; i<number_photons; i++)
1215  arrival_time_distrb.push_back(fVisTiming->GetRandom());
1216 
1217  //deleting ...
1218 
1219  delete flandau;
1220  delete fexpo;
1221  delete fint;
1222  delete fVisTiming;
1223  //G4cout<<"Timing distribution BEAMAUS!"<<G4endl;
1224  return arrival_time_distrb;
1225  }*/
1226 
1227  // New Parametrization code
1228  // parameterisation generation function
1229  void OpFastScintillation::generateParam(const size_t index, const size_t angle_bin)
1230  {
1231  // get distance
1232  double distance_in_cm = (index * fstep_size) + fmin_d;
1233 
1234  // time range
1235  const double signal_t_range = 5000.; // TODO: unhardcode
1236 
1237  // parameterisation TF1
1238  TF1 fVUVTiming;
1239 
1240  // For very short distances the time correction is just a shift
1241  double t_direct_mean = distance_in_cm / fvuv_vgroup_mean;
1242  double t_direct_min = distance_in_cm / fvuv_vgroup_max;
1243 
1244  // Defining the model function(s) describing the photon transportation timing vs distance
1245  // Getting the landau parameters from the time parametrization
1246  std::array<double, 3> pars_landau;
1247  interpolate3(pars_landau,
1248  fparameters[0][0],
1249  fparameters[2][angle_bin],
1250  fparameters[3][angle_bin],
1251  fparameters[1][angle_bin],
1252  distance_in_cm,
1253  true);
1254  // Deciding which time model to use (depends on the distance)
1255  // defining useful times for the VUV arrival time shapes
1256  if (distance_in_cm >= finflexion_point_distance) {
1257  double pars_far[4] = {t_direct_min, pars_landau[0], pars_landau[1], pars_landau[2]};
1258  // Set model: Landau
1259  fVUVTiming = TF1("fVUVTiming", model_far, 0, signal_t_range, 4);
1260  fVUVTiming.SetParameters(pars_far);
1261  }
1262  else {
1263  // Set model: Landau + Exponential
1264  fVUVTiming = TF1("fVUVTiming", model_close, 0, signal_t_range, 7);
1265  // Exponential parameters
1266  double pars_expo[2];
1267  // Getting the exponential parameters from the time parametrization
1268  pars_expo[1] =
1269  interpolate(fparameters[4][0], fparameters[5][angle_bin], distance_in_cm, true);
1270  pars_expo[0] =
1271  interpolate(fparameters[4][0], fparameters[6][angle_bin], distance_in_cm, true);
1272  pars_expo[0] *= pars_landau[2];
1273  pars_expo[0] = std::log(pars_expo[0]);
1274  // this is to find the intersection point between the two functions:
1275  TF1 fint = TF1("fint", finter_d, pars_landau[0], 4 * t_direct_mean, 5);
1276  double parsInt[5] = {
1277  pars_landau[0], pars_landau[1], pars_landau[2], pars_expo[0], pars_expo[1]};
1278  fint.SetParameters(parsInt);
1279  double t_int = fint.GetMinimumX();
1280  double minVal = fint.Eval(t_int);
1281  // the functions must intersect - output warning if they don't
1282  if (minVal > 0.015) {
1283  std::cout << "WARNING: Parametrization of VUV light discontinuous for distance = "
1284  << distance_in_cm << std::endl;
1285  std::cout << "WARNING: This shouldn't be happening " << std::endl;
1286  }
1287  double parsfinal[7] = {t_int,
1288  pars_landau[0],
1289  pars_landau[1],
1290  pars_landau[2],
1291  pars_expo[0],
1292  pars_expo[1],
1293  t_direct_min};
1294  fVUVTiming.SetParameters(parsfinal);
1295  }
1296 
1297  // set the number of points used to sample parameterisation
1298  // for shorter distances, peak is sharper so more sensitive sampling required
1299  int fsampling; // TODO: unhardcode
1300  if (distance_in_cm < 50) { fsampling = 10000; }
1301  else if (distance_in_cm < 100) {
1302  fsampling = 5000;
1303  }
1304  else {
1305  fsampling = 1000;
1306  }
1307  fVUVTiming.SetNpx(fsampling);
1308 
1309  // calculate max and min distance relevant to sample parameterisation
1310  // max
1311  const size_t nq_max = 1;
1312  double xq_max[nq_max];
1313  double yq_max[nq_max];
1314  xq_max[0] = 0.995; // include 99.5%
1315  fVUVTiming.GetQuantiles(nq_max, yq_max, xq_max);
1316  double max = yq_max[0];
1317  // min
1318  double min = t_direct_min;
1319 
1320  // store TF1 and min/max, this allows identical TF1 to be used every time sampling
1321  // the first call of GetRandom generates the timing sampling and stores it in the TF1 object, this is the slow part
1322  // all subsequent calls check if it has been generated previously and are ~100+ times quicker
1323  VUV_timing[angle_bin][index] = fVUVTiming;
1324  VUV_max[angle_bin][index] = max;
1325  VUV_min[angle_bin][index] = min;
1326  }
1327 
1328  // VUV arrival times calculation function
1329  void OpFastScintillation::getVUVTimes(std::vector<double>& arrivalTimes,
1330  const double distance,
1331  const size_t angle_bin)
1332  {
1333  if (distance < fmin_d) {
1334  // times are fixed shift i.e. direct path only
1335  double t_prop_correction = distance / fvuv_vgroup_mean;
1336  for (size_t i = 0; i < arrivalTimes.size(); ++i) {
1337  arrivalTimes[i] = t_prop_correction;
1338  }
1339  }
1340  else { // distance >= fmin_d
1341  // determine nearest parameterisation in discretisation
1342  int index = std::round((distance - fmin_d) / fstep_size);
1343  // check whether required parameterisation has been generated, generating if not
1344  if (VUV_timing[angle_bin][index].GetNdim() == 0) { generateParam(index, angle_bin); }
1345  // randomly sample parameterisation for each photon
1346  for (size_t i = 0; i < arrivalTimes.size(); ++i) {
1347  arrivalTimes[i] = VUV_timing[angle_bin][index].GetRandom(VUV_min[angle_bin][index],
1348  VUV_max[angle_bin][index]);
1349  }
1350  }
1351  }
1352 
1353  // VIS arrival times calculation functions
1354  void OpFastScintillation::getVISTimes(std::vector<double>& arrivalTimes,
1355  const TVector3& ScintPoint,
1356  const TVector3& OpDetPoint)
1357  {
1358  // *************************************************************************************************
1359  // Calculation of earliest arrival times and corresponding unsmeared distribution
1360  // *************************************************************************************************
1361 
1362  // set plane_depth for correct TPC:
1363  double plane_depth;
1364  if (ScintPoint[0] < 0) { plane_depth = -fplane_depth; }
1365  else {
1366  plane_depth = fplane_depth;
1367  }
1368 
1369  // calculate point of reflection for shortest path
1370  TVector3 bounce_point(plane_depth, ScintPoint[1], ScintPoint[2]);
1371 
1372  // calculate distance travelled by VUV light and by vis light
1373  double VUVdist = (bounce_point - ScintPoint).Mag();
1374  double Visdist = (OpDetPoint - bounce_point).Mag();
1375 
1376  // calculate times taken by VUV part of path
1377  int angle_bin_vuv = 0; // on-axis by definition
1378  getVUVTimes(arrivalTimes, VUVdist, angle_bin_vuv);
1379 
1380  // add visible direct path transport time
1381  for (size_t i = 0; i < arrivalTimes.size(); ++i) {
1382  arrivalTimes[i] += Visdist / fvis_vmean;
1383  }
1384 
1385  // *************************************************************************************************
1386  // Smearing of arrival time distribution
1387  // *************************************************************************************************
1388  // calculate fastest time possible
1389  // vis part
1390  double vis_time = Visdist / fvis_vmean;
1391  // vuv part
1392  double vuv_time;
1393  if (VUVdist < fmin_d) { vuv_time = VUVdist / fvuv_vgroup_max; }
1394  else {
1395  // find index of required parameterisation
1396  const size_t index = std::round((VUVdist - fmin_d) / fstep_size);
1397  // find shortest time
1398  vuv_time = VUV_min[angle_bin_vuv][index];
1399  }
1400  // sum
1401  double fastest_time = vis_time + vuv_time;
1402 
1403  // calculate angle theta between bound_point and optical detector
1404  double cosine_theta = std::abs(OpDetPoint[0] - bounce_point[0]) / Visdist;
1405  double theta = fast_acos(cosine_theta) * 180. / CLHEP::pi;
1406 
1407  // determine smearing parameters using interpolation of generated points:
1408  // 1). tau = exponential smearing factor, varies with distance and angle
1409  // 2). cutoff = largest smeared time allowed, preventing excessively large
1410  // times caused by exponential distance to cathode
1411  double distance_cathode_plane = std::abs(plane_depth - ScintPoint[0]);
1412  // angular bin
1413  size_t theta_bin = theta / fangle_bin_timing_vis;
1414  // radial distance from centre of TPC (y,z plane)
1415  double r = std::sqrt(std::pow(ScintPoint.Y() - fcathode_centre.Z(), 2) +
1416  std::pow(ScintPoint.Y() - fcathode_centre.Z(), 2));
1417 
1418  // cut-off and tau
1419  // cut-off
1420  // interpolate in d_c for each r bin
1421  std::vector<double> interp_vals(fcut_off_pars[theta_bin].size(), 0.0);
1422  for (size_t i = 0; i < fcut_off_pars[theta_bin].size(); i++) {
1423  interp_vals[i] =
1424  interpolate(fdistances_refl, fcut_off_pars[theta_bin][i], distance_cathode_plane, true);
1425  }
1426  // interpolate in r
1427  double cutoff = interpolate(fradial_distances_refl, interp_vals, r, true);
1428 
1429  // tau
1430  // interpolate in x for each r bin
1431  std::vector<double> interp_vals_tau(ftau_pars[theta_bin].size(), 0.0);
1432  for (size_t i = 0; i < ftau_pars[theta_bin].size(); i++) {
1433  interp_vals_tau[i] =
1434  interpolate(fdistances_refl, ftau_pars[theta_bin][i], distance_cathode_plane, true);
1435  }
1436  // interpolate in r
1437  double tau = interpolate(fradial_distances_refl, interp_vals_tau, r, true);
1438 
1439  if (tau < 0) { tau = 0; } // failsafe if tau extrapolate goes wrong
1440 
1441  // apply smearing:
1442  for (size_t i = 0; i < arrivalTimes.size(); ++i) {
1443  double arrival_time_smeared;
1444  // if time is already greater than cutoff, do not apply smearing
1445  if (arrivalTimes[i] >= cutoff) { continue; }
1446  // otherwise smear
1447  else {
1448  unsigned int counter = 0;
1449  // loop until time generated is within cutoff limit
1450  // most are within single attempt, very few take more than two
1451  do {
1452  // don't attempt smearings too many times
1453  if (counter >= 10) { // TODO: unhardcode
1454  arrival_time_smeared = arrivalTimes[i]; // don't smear
1455  break;
1456  }
1457  else {
1458  // generate random number in appropriate range
1459  double x = gRandom->Uniform(0.5, 1.0); // TODO: unhardcode
1460  // apply the exponential smearing
1461  arrival_time_smeared =
1462  arrivalTimes[i] + (arrivalTimes[i] - fastest_time) * (std::pow(x, -tau) - 1);
1463  }
1464  counter++;
1465  } while (arrival_time_smeared > cutoff);
1466  }
1467  arrivalTimes[i] = arrival_time_smeared;
1468  }
1469  }
1470 
1471  // ---------------------------------------------------------------------------
1473  {
1474  return fUseNhitsModel;
1475  } // OpFastScintillation::usesSemiAnalyticModel()
1476 
1477  // ---------------------------------------------------------------------------
1478  void OpFastScintillation::detectedDirectHits(std::map<size_t, int>& DetectedNum,
1479  const double Num,
1480  geo::Point_t const& ScintPoint) const
1481  {
1482  for (size_t const OpDet : util::counter(fPVS->NOpChannels())) {
1483  if (!isOpDetInSameTPC(ScintPoint, fOpDetCenter.at(OpDet))) continue;
1484 
1485  // set detector struct for solid angle function
1487  fOpDetLength.at(OpDet),
1488  fOpDetCenter.at(OpDet),
1489  fOpDetType.at(OpDet)};
1490  const int DetThis = VUVHits(Num, ScintPoint, op);
1491  if (DetThis > 0) {
1492  DetectedNum[OpDet] = DetThis;
1493  // mf::LogInfo("OpFastScintillation") << "FastScint: " <<
1494  // // it->second<<" " << Num << " " << DetThisPMT;
1495  //det_photon_ctr += DetThisPMT; // CASE-DEBUG DO NOT REMOVE THIS COMMENT
1496  }
1497  }
1498  }
1499 
1500  void OpFastScintillation::detectedReflecHits(std::map<size_t, int>& ReflDetectedNum,
1501  const double Num,
1502  geo::Point_t const& ScintPoint) const
1503  {
1504  // 1). calculate total number of hits of VUV photons on
1505  // reflective foils via solid angle + Gaisser-Hillas
1506  // corrections:
1507 
1508  // set plane_depth for correct TPC:
1509  double plane_depth;
1510  if (ScintPoint.X() < 0.) { plane_depth = -fplane_depth; }
1511  else {
1512  plane_depth = fplane_depth;
1513  }
1514 
1515  // get scintpoint coords relative to centre of cathode plane
1516  std::array<double, 3> ScintPoint_relative = {std::abs(ScintPoint.X() - plane_depth),
1517  std::abs(ScintPoint.Y() - fcathode_centre.Y()),
1518  std::abs(ScintPoint.Z() - fcathode_centre.Z())};
1519  // calculate solid angle of cathode from the scintillation point
1520  double solid_angle_cathode = Rectangle_SolidAngle(fcathode_plane, ScintPoint_relative);
1521 
1522  // calculate distance and angle between ScintPoint and hotspot
1523  // vast majority of hits in hotspot region directly infront of scintpoint,
1524  // therefore consider attenuation for this distance and on axis GH instead of for the centre coordinate
1525  double distance_cathode = std::abs(plane_depth - ScintPoint.X());
1526  // calculate hits on cathode plane via geometric acceptance
1527  double cathode_hits_geo = std::exp(-1. * distance_cathode / fL_abs_vuv) *
1528  (solid_angle_cathode / (4. * CLHEP::pi)) * Num;
1529  // determine Gaisser-Hillas correction including border effects
1530  // use flat correction
1531  double r = std::sqrt(std::pow(ScintPoint.Y() - fcathode_centre.Y(), 2) +
1532  std::pow(ScintPoint.Z() - fcathode_centre.Z(), 2));
1533  double pars_ini[4] = {0, 0, 0, 0};
1534  double s1 = 0;
1535  double s2 = 0;
1536  double s3 = 0;
1537  if (fIsFlatPDCorr) {
1538  pars_ini[0] = fGHvuvpars_flat[0][0];
1539  pars_ini[1] = fGHvuvpars_flat[1][0];
1540  pars_ini[2] = fGHvuvpars_flat[2][0];
1541  pars_ini[3] = fGHvuvpars_flat[3][0];
1545  }
1546  else
1547  std::cout
1548  << "Error: flat optical detector VUV correction required for reflected semi-analytic hits."
1549  << std::endl;
1550 
1551  // add border correction
1552  pars_ini[0] = pars_ini[0] + s1 * r;
1553  pars_ini[1] = pars_ini[1] + s2 * r;
1554  pars_ini[2] = pars_ini[2] + s3 * r;
1555  pars_ini[3] = pars_ini[3];
1556 
1557  // calculate corrected number of hits
1558  double GH_correction = Gaisser_Hillas(distance_cathode, pars_ini);
1559  const double cathode_hits_rec = GH_correction * cathode_hits_geo;
1560 
1561  // detemine hits on each PD
1562  const std::array<double, 3> hotspot = {plane_depth, ScintPoint.Y(), ScintPoint.Z()};
1563  for (size_t const OpDet : util::counter(fPVS->NOpChannels())) {
1564  if (!isOpDetInSameTPC(ScintPoint, fOpDetCenter.at(OpDet))) continue;
1565 
1566  // set detector struct for solid angle function
1568  fOpDetLength.at(OpDet),
1569  fOpDetCenter.at(OpDet),
1570  fOpDetType.at(OpDet)};
1571 
1572  int const ReflDetThis = VISHits(ScintPoint, op, cathode_hits_rec, hotspot);
1573  if (ReflDetThis > 0) { ReflDetectedNum[OpDet] = ReflDetThis; }
1574  }
1575  }
1576 
1577  // VUV semi-analytic hits calculation
1578  int OpFastScintillation::VUVHits(const double Nphotons_created,
1579  geo::Point_t const& ScintPoint_v,
1580  OpticalDetector const& opDet) const
1581  {
1582  // the interface has been converted into geo::Point_t, the implementation not yet
1583  std::array<double, 3U> ScintPoint;
1584  std::array<double, 3U> OpDetPoint;
1585  geo::vect::fillCoords(ScintPoint, ScintPoint_v);
1586  geo::vect::fillCoords(OpDetPoint, opDet.OpDetPoint);
1587 
1588  // distance and angle between ScintPoint and OpDetPoint
1589  double distance = dist(&ScintPoint[0], &OpDetPoint[0], 3);
1590  double cosine = std::abs(ScintPoint[0] - OpDetPoint[0]) / distance;
1591  // double theta = std::acos(cosine) * 180. / CLHEP::pi;
1592  double theta = fast_acos(cosine) * 180. / CLHEP::pi;
1593 
1594  // calculate solid angle:
1595  double solid_angle = 0;
1596  // Arapucas/Bars (rectangle)
1597  if (opDet.type == 0) {
1598  // get scintillation point coordinates relative to arapuca window centre
1599  std::array<double, 3> ScintPoint_rel = {std::abs(ScintPoint[0] - OpDetPoint[0]),
1600  std::abs(ScintPoint[1] - OpDetPoint[1]),
1601  std::abs(ScintPoint[2] - OpDetPoint[2])};
1602  // calculate solid angle
1603  solid_angle = Rectangle_SolidAngle(Dims{opDet.h, opDet.w}, ScintPoint_rel);
1604  }
1605  // PMTs (dome)
1606  else if (opDet.type == 1) {
1607  solid_angle = Omega_Dome_Model(distance, theta);
1608  }
1609  // PMTs (disk approximation)
1610  else if (opDet.type == 2) {
1611  // offset in z-y plane
1612  double d = dist(&ScintPoint[1], &OpDetPoint[1], 2);
1613  // drift distance (in x)
1614  double h = std::abs(ScintPoint[0] - OpDetPoint[0]);
1615  // Solid angle of a disk
1616  solid_angle = Disk_SolidAngle(d, h, fradius);
1617  }
1618  else {
1619  std::cout << "Error: Invalid optical detector type. 0 = rectangular, 1 = dome, 2 = disk"
1620  << std::endl;
1621  }
1622 
1623  // calculate number of photons hits by geometric acceptance: accounting for solid angle and LAr absorbtion length
1624  double hits_geo =
1625  std::exp(-1. * distance / fL_abs_vuv) * (solid_angle / (4 * CLHEP::pi)) * Nphotons_created;
1626 
1627  // apply Gaisser-Hillas correction for Rayleigh scattering distance
1628  // and angular dependence offset angle bin
1629  const size_t j = (theta / fdelta_angulo_vuv);
1630 
1631  // determine GH parameters, accounting for border effects
1632  // radial distance from centre of detector (Y-Z)
1633  // FIXME (KJK): Original implementation only accounts for the Y component
1634  double r = std::abs(ScintPoint[1] - fcathode_centre.Y());
1635  double pars_ini[4] = {0, 0, 0, 0};
1636  double s1 = 0;
1637  double s2 = 0;
1638  double s3 = 0;
1639  // flat PDs
1640  if ((opDet.type == 0 || opDet.type == 2) && fIsFlatPDCorr) {
1641  pars_ini[0] = fGHvuvpars_flat[0][j];
1642  pars_ini[1] = fGHvuvpars_flat[1][j];
1643  pars_ini[2] = fGHvuvpars_flat[2][j];
1644  pars_ini[3] = fGHvuvpars_flat[3][j];
1648  }
1649  // dome PDs
1650  else if (opDet.type == 1 && fIsDomePDCorr) {
1651  pars_ini[0] = fGHvuvpars_dome[0][j];
1652  pars_ini[1] = fGHvuvpars_dome[1][j];
1653  pars_ini[2] = fGHvuvpars_dome[2][j];
1654  pars_ini[3] = fGHvuvpars_dome[3][j];
1658  }
1659  else
1660  std::cout << "Error: Invalid optical detector type. 0 = rectangular, 1 = dome, 2 = disk. Or "
1661  "corrections for chosen optical detector type missing."
1662  << std::endl;
1663 
1664  // add border correction
1665  pars_ini[0] = pars_ini[0] + s1 * r;
1666  pars_ini[1] = pars_ini[1] + s2 * r;
1667  pars_ini[2] = pars_ini[2] + s3 * r;
1668  pars_ini[3] = pars_ini[3];
1669 
1670  // calculate correction
1671  double GH_correction = Gaisser_Hillas(distance, pars_ini);
1672 
1673  // calculate number photons
1674  double hits_rec = GH_correction * hits_geo / cosine;
1675  int hits_vuv = std::round(G4Poisson(hits_rec));
1676 
1677  return hits_vuv;
1678  }
1679 
1680  // VIS hits semi-analytic model calculation
1682  OpticalDetector const& opDet,
1683  const double cathode_hits_rec,
1684  const std::array<double, 3> hotspot) const
1685  {
1686  // 1). calculate total number of hits of VUV photons on reflective
1687  // foils via solid angle + Gaisser-Hillas corrections.
1688  // Done outside as it doesn't depend on OpDetPoint
1689 
1690  // set plane_depth for correct TPC:
1691  double plane_depth;
1692  if (ScintPoint_v.X() < 0.) { plane_depth = -fplane_depth; }
1693  else {
1694  plane_depth = fplane_depth;
1695  }
1696 
1697  // 2). calculate number of these hits which reach the optical
1698  // detector from the hotspot via solid angle
1699 
1700  // the interface has been converted into geo::Point_t, the implementation not yet
1701  std::array<double, 3U> ScintPoint;
1702  std::array<double, 3U> OpDetPoint;
1703  geo::vect::fillCoords(ScintPoint, ScintPoint_v);
1704  geo::vect::fillCoords(OpDetPoint, opDet.OpDetPoint);
1705 
1706  // calculate distances and angles for application of corrections
1707  // distance from hotspot to optical detector
1708  double distance_vis = dist(&hotspot[0], &OpDetPoint[0], 3);
1709  // angle between hotspot and optical detector
1710  double cosine_vis = std::abs(hotspot[0] - OpDetPoint[0]) / distance_vis;
1711  // double theta_vis = std::acos(cosine_vis) * 180. / CLHEP::pi;
1712  double theta_vis = fast_acos(cosine_vis) * 180. / CLHEP::pi;
1713 
1714  // calculate solid angle of optical detector
1715  double solid_angle_detector = 0;
1716  // rectangular aperture
1717  if (opDet.type == 0) {
1718  // get hotspot coordinates relative to opDet
1719  std::array<double, 3> emission_relative = {std::abs(hotspot[0] - OpDetPoint[0]),
1720  std::abs(hotspot[1] - OpDetPoint[1]),
1721  std::abs(hotspot[2] - OpDetPoint[2])};
1722  // calculate solid angle
1723  solid_angle_detector = Rectangle_SolidAngle(Dims{opDet.h, opDet.w}, emission_relative);
1724  }
1725  // dome aperture
1726  else if (opDet.type == 1) {
1727  solid_angle_detector = Omega_Dome_Model(distance_vis, theta_vis);
1728  }
1729  // disk aperture
1730  else if (opDet.type == 2) {
1731  // offset in z-y plane
1732  double d = dist(&hotspot[1], &OpDetPoint[1], 2);
1733  // drift distance (in x)
1734  double h = std::abs(hotspot[0] - OpDetPoint[0]);
1735  // calculate solid angle
1736  solid_angle_detector = Disk_SolidAngle(d, h, fradius);
1737  }
1738  else {
1739  std::cout << "Error: Invalid optical detector type. 0 = rectangular, 1 = dome, 2 = disk"
1740  << std::endl;
1741  }
1742 
1743  // calculate number of hits via geometeric acceptance
1744  double hits_geo = (solid_angle_detector / (2. * CLHEP::pi)) *
1745  cathode_hits_rec; // 2*pi due to presence of reflective foils
1746 
1747  // determine correction factor, depending on PD type
1748  const size_t k = (theta_vis / fdelta_angulo_vis); // off-set angle bin
1749  // radial distance from centre of detector (Y-Z)
1750  // FIXME (KJK): Original implementation only accounts for the Y component
1751  double r = std::abs(ScintPoint[1] - fcathode_centre.Y());
1752  double d_c = std::abs(ScintPoint[0] - plane_depth); // distance to cathode
1753  double border_correction = 0;
1754  // flat PDs
1755  if ((opDet.type == 0 || opDet.type == 2) && fIsFlatPDCorr) {
1756  // interpolate in d_c for each r bin
1757  const size_t nbins_r = fvispars_flat[k].size();
1758  std::vector<double> interp_vals(nbins_r, 0.0);
1759  {
1760  size_t idx = 0;
1761  size_t size = fvis_distances_x_flat.size();
1762  if (d_c >= fvis_distances_x_flat[size - 2])
1763  idx = size - 2;
1764  else {
1765  while (d_c > fvis_distances_x_flat[idx + 1])
1766  idx++;
1767  }
1768  for (size_t i = 0; i < nbins_r; ++i) {
1769  interp_vals[i] = interpolate(fvis_distances_x_flat, fvispars_flat[k][i], d_c, false, idx);
1770  }
1771  }
1772  // interpolate in r
1773  border_correction = interpolate(fvis_distances_r_flat, interp_vals, r, false);
1774  }
1775  // dome PDs
1776  else if (opDet.type == 1 && fIsDomePDCorr) {
1777  // interpolate in d_c for each r bin
1778  const size_t nbins_r = fvispars_dome[k].size();
1779  std::vector<double> interp_vals(nbins_r, 0.0);
1780  {
1781  size_t idx = 0;
1782  size_t size = fvis_distances_x_dome.size();
1783  if (d_c >= fvis_distances_x_dome[size - 2])
1784  idx = size - 2;
1785  else {
1786  while (d_c > fvis_distances_x_dome[idx + 1])
1787  idx++;
1788  }
1789  for (size_t i = 0; i < nbins_r; ++i) {
1790  interp_vals[i] = interpolate(fvis_distances_x_dome, fvispars_dome[k][i], d_c, false, idx);
1791  }
1792  }
1793  // interpolate in r
1794  border_correction = interpolate(fvis_distances_r_dome, interp_vals, r, false);
1795  }
1796  else {
1797  std::cout << "Error: Invalid optical detector type. 0 = rectangular, 1 = dome, 2 = disk. Or "
1798  "corrections for chosen optical detector type missing."
1799  << std::endl;
1800  }
1801 
1802  // apply correction factor
1803  double hits_rec = border_correction * hits_geo / cosine_vis;
1804  double hits_vis = std::round(G4Poisson(hits_rec));
1805 
1806  return hits_vis;
1807  }
1808 
1810  geo::Point_t const& OpDetPoint) const
1811  {
1812  // check optical channel is in same TPC as scintillation light, if not return 0 hits
1813  // temporary method working for SBND, uBooNE, DUNE 1x2x6; to be replaced to work in full DUNE geometry
1814  // check x coordinate has same sign or is close to zero, otherwise return 0 hits
1815  if (((ScintPoint.X() < 0.) != (OpDetPoint.X() < 0.)) &&
1816  std::abs(OpDetPoint.X()) > 10.) { // TODO: unhardcode
1817  return false;
1818  }
1819  return true;
1820  }
1821 
1823  {
1824  //semi-analytic approach only works in the active volume
1825  return fActiveVolumes[0].ContainsPosition(ScintPoint);
1826  }
1827 
1828  G4double OpFastScintillation::single_exp(const G4double t, const G4double tau2) const
1829  {
1830  return std::exp(-1.0 * t / tau2) / tau2;
1831  }
1832 
1833  G4double OpFastScintillation::bi_exp(const G4double t,
1834  const G4double tau1,
1835  const G4double tau2) const
1836  { // TODO: what's up with this? ... / tau2 / tau2 ...
1837  return std::exp(-1.0 * t / tau2) * (1 - std::exp(-1.0 * t / tau1)) / tau2 / tau2 *
1838  (tau1 + tau2);
1839  }
1840 
1841  G4double OpFastScintillation::Gaisser_Hillas(const double x, const double* par) const
1842  {
1843  double X_mu_0 = par[3];
1844  double Normalization = par[0];
1845  double Diff = par[1] - X_mu_0;
1846  double Term = std::pow((x - X_mu_0) / Diff, Diff / par[2]);
1847  double Exponential = std::exp((par[1] - x) / par[2]);
1848  return (Normalization * Term * Exponential);
1849  }
1850 
1851  double finter_d(double* x, double* par)
1852  {
1853  double y1 = par[2] * TMath::Landau(x[0], par[0], par[1]);
1854  double y2 = TMath::Exp(par[3] + x[0] * par[4]);
1855  return TMath::Abs(y1 - y2);
1856  }
1857 
1858  double LandauPlusExpoFinal(double* x, double* par)
1859  {
1860  // par0 = joining point
1861  // par1 = Landau MPV
1862  // par2 = Landau widt
1863  // par3 = normalization
1864  // par4 = Expo cte
1865  // par5 = Expo tau
1866  double y1 = par[3] * TMath::Landau(x[0], par[1], par[2]);
1867  double y2 = TMath::Exp(par[4] + x[0] * par[5]);
1868  if (x[0] > par[0]) y1 = 0.;
1869  if (x[0] < par[0]) y2 = 0.;
1870  return (y1 + y2);
1871  }
1872 
1873  double finter_r(double* x, double* par)
1874  {
1875  double y1 = par[2] * TMath::Landau(x[0], par[0], par[1]);
1876  double y2 = par[5] * TMath::Landau(x[0], par[3], par[4]);
1877  return TMath::Abs(y1 - y2);
1878  }
1879 
1880  double model_close(double* x, double* par)
1881  {
1882  // par0 = joining point
1883  // par1 = Landau MPV
1884  // par2 = Landau width
1885  // par3 = normalization
1886  // par4 = Expo cte
1887  // par5 = Expo tau
1888  // par6 = t_min
1889  double y1 = par[3] * TMath::Landau(x[0], par[1], par[2]);
1890  double y2 = TMath::Exp(par[4] + x[0] * par[5]);
1891  if (x[0] <= par[6] || x[0] > par[0]) y1 = 0.;
1892  if (x[0] < par[0]) y2 = 0.;
1893  return (y1 + y2);
1894  }
1895 
1896  double model_far(double* x, double* par)
1897  {
1898  // par1 = Landau MPV
1899  // par2 = Landau width
1900  // par3 = normalization
1901  // par0 = t_min
1902  double y = par[3] * TMath::Landau(x[0], par[1], par[2]);
1903  if (x[0] <= par[0]) y = 0.;
1904  return y;
1905  }
1906 
1907  //======================================================================
1908  // Returns interpolated value at x from parallel arrays ( xData, yData )
1909  // Assumes that xData has at least two elements, is sorted and is strictly monotonic increasing
1910  // boolean argument extrapolate determines behaviour beyond ends of array (if needed)
1911  double OpFastScintillation::interpolate(const std::vector<double>& xData,
1912  const std::vector<double>& yData,
1913  double x,
1914  bool extrapolate,
1915  size_t i) const
1916  {
1917  if (i == 0) {
1918  size_t size = xData.size();
1919  if (x >= xData[size - 2]) { // special case: beyond right end
1920  i = size - 2;
1921  }
1922  else {
1923  while (x > xData[i + 1])
1924  i++;
1925  }
1926  }
1927  double xL = xData[i];
1928  double xR = xData[i + 1];
1929  double yL = yData[i];
1930  double yR = yData[i + 1]; // points on either side (unless beyond ends)
1931  if (!extrapolate) { // if beyond ends of array and not extrapolating
1932  if (x < xL) return yL;
1933  if (x > xR) return yL;
1934  }
1935  const double dydx = (yR - yL) / (xR - xL); // gradient
1936  return yL + dydx * (x - xL); // linear interpolation
1937  }
1938 
1939  void OpFastScintillation::interpolate3(std::array<double, 3>& inter,
1940  const std::vector<double>& xData,
1941  const std::vector<double>& yData1,
1942  const std::vector<double>& yData2,
1943  const std::vector<double>& yData3,
1944  double x,
1945  bool extrapolate) const
1946  {
1947  size_t size = xData.size();
1948  size_t i = 0; // find left end of interval for interpolation
1949  if (x >= xData[size - 2]) { // special case: beyond right end
1950  i = size - 2;
1951  }
1952  else {
1953  while (x > xData[i + 1])
1954  i++;
1955  }
1956  double xL = xData[i];
1957  double xR = xData[i + 1]; // points on either side (unless beyond ends)
1958  double yL1 = yData1[i];
1959  double yR1 = yData1[i + 1];
1960  double yL2 = yData2[i];
1961  double yR2 = yData2[i + 1];
1962  double yL3 = yData3[i];
1963  double yR3 = yData3[i + 1];
1964 
1965  if (!extrapolate) { // if beyond ends of array and not extrapolating
1966  if (x < xL) {
1967  inter[0] = yL1;
1968  inter[1] = yL2;
1969  inter[2] = yL3;
1970  return;
1971  }
1972  if (x > xR) {
1973  inter[0] = yL1;
1974  inter[1] = yL2;
1975  inter[2] = yL3;
1976  return;
1977  }
1978  }
1979  const double m = (x - xL) / (xR - xL);
1980  inter[0] = m * (yR1 - yL1) + yL1;
1981  inter[1] = m * (yR2 - yL2) + yL2;
1982  inter[2] = m * (yR3 - yL3) + yL3;
1983  }
1984 
1985  // solid angle of circular aperture
1986  // TODO: allow greater tolerance in comparisons, by default its using:
1987  // std::numeric_limits<double>::epsilon(): 2.22045e-16
1988  // that's an unrealistic small number, better setting
1989  // constexpr double tolerance = 0.0000001; // 1 nm
1990  double OpFastScintillation::Disk_SolidAngle(const double d, const double h, const double b) const
1991  {
1992  if (b <= 0. || d < 0. || h <= 0.) return 0.;
1993  const double leg2 = (b + d) * (b + d);
1994  const double aa = std::sqrt(h * h / (h * h + leg2));
1995  if (isApproximatelyZero(d)) { return 2. * CLHEP::pi * (1. - aa); }
1996  double bb = 2. * std::sqrt(b * d / (h * h + leg2));
1997  double cc = 4. * b * d / leg2;
1998 
1999  if (isDefinitelyGreaterThan(d, b)) {
2000  try {
2001  return 2. * aa *
2002  (std::sqrt(1. - cc) * boost::math::ellint_3(bb, cc, noLDoublePromote()) -
2003  boost::math::ellint_1(bb, noLDoublePromote()));
2004  }
2005  catch (std::domain_error& e) {
2006  if (isApproximatelyEqual(d, b, 1e-9)) {
2007  mf::LogWarning("OpFastScintillation")
2008  << "Elliptic Integral in Disk_SolidAngle() given parameters outside domain."
2009  << "\nbb: " << bb << "\ncc: " << cc << "\nException message: " << e.what()
2010  << "\nRelax condition and carry on.";
2011  return CLHEP::pi - 2. * aa * boost::math::ellint_1(bb, noLDoublePromote());
2012  }
2013  else {
2014  mf::LogError("OpFastScintillation")
2015  << "Elliptic Integral inside Disk_SolidAngle() given parameters outside domain.\n"
2016  << "\nbb: " << bb << "\ncc: " << cc << "Exception message: " << e.what();
2017  return 0.;
2018  }
2019  }
2020  }
2021  if (isDefinitelyLessThan(d, b)) {
2022  try {
2023  return 2. * CLHEP::pi -
2024  2. * aa *
2025  (boost::math::ellint_1(bb, noLDoublePromote()) +
2026  std::sqrt(1. - cc) * boost::math::ellint_3(bb, cc, noLDoublePromote()));
2027  }
2028  catch (std::domain_error& e) {
2029  if (isApproximatelyEqual(d, b, 1e-9)) {
2030  mf::LogWarning("OpFastScintillation")
2031  << "Elliptic Integral in Disk_SolidAngle() given parameters outside domain."
2032  << "\nbb: " << bb << "\ncc: " << cc << "\nException message: " << e.what()
2033  << "\nRelax condition and carry on.";
2034  return CLHEP::pi - 2. * aa * boost::math::ellint_1(bb, noLDoublePromote());
2035  }
2036  else {
2037  mf::LogError("OpFastScintillation")
2038  << "Elliptic Integral inside Disk_SolidAngle() given parameters outside domain.\n"
2039  << "\nbb: " << bb << "\ncc: " << cc << "Exception message: " << e.what();
2040  return 0.;
2041  }
2042  }
2043  }
2044  if (isApproximatelyEqual(d, b)) {
2045  return CLHEP::pi - 2. * aa * boost::math::ellint_1(bb, noLDoublePromote());
2046  }
2047  return 0.;
2048  }
2049 
2050  // solid angle of rectangular aperture
2052  const double b,
2053  const double d) const
2054  {
2055  double aa = a / (2. * d);
2056  double bb = b / (2. * d);
2057  double aux = (1. + aa * aa + bb * bb) / ((1. + aa * aa) * (1. + bb * bb));
2058  // return 4 * std::acos(std::sqrt(aux));
2059  return 4. * fast_acos(std::sqrt(aux));
2060  }
2061 
2062  // TODO: allow greater tolerance in comparisons, see note above on Disk_SolidAngle()
2064  const std::array<double, 3> v) const
2065  {
2066  // v is the position of the track segment with respect to
2067  // the center position of the arapuca window
2068 
2069  // arapuca plane fixed in x direction
2070  if (isApproximatelyZero(v[1]) && isApproximatelyZero(v[2])) {
2071  return Rectangle_SolidAngle(o.h, o.w, v[0]);
2072  }
2073  if (isDefinitelyGreaterThan(v[1], o.h * .5) && isDefinitelyGreaterThan(v[2], o.w * .5)) {
2074  double A = v[1] - o.h * .5;
2075  double B = v[2] - o.w * .5;
2076  double to_return = (Rectangle_SolidAngle(2. * (A + o.h), 2. * (B + o.w), v[0]) -
2077  Rectangle_SolidAngle(2. * A, 2. * (B + o.w), v[0]) -
2078  Rectangle_SolidAngle(2. * (A + o.h), 2. * B, v[0]) +
2079  Rectangle_SolidAngle(2. * A, 2. * B, v[0])) *
2080  .25;
2081  return to_return;
2082  }
2083  if ((v[1] <= o.h * .5) && (v[2] <= o.w * .5)) {
2084  double A = -v[1] + o.h * .5;
2085  double B = -v[2] + o.w * .5;
2086  double to_return = (Rectangle_SolidAngle(2. * (o.h - A), 2. * (o.w - B), v[0]) +
2087  Rectangle_SolidAngle(2. * A, 2. * (o.w - B), v[0]) +
2088  Rectangle_SolidAngle(2. * (o.h - A), 2. * B, v[0]) +
2089  Rectangle_SolidAngle(2. * A, 2. * B, v[0])) *
2090  .25;
2091  return to_return;
2092  }
2093  if (isDefinitelyGreaterThan(v[1], o.h * .5) && (v[2] <= o.w * .5)) {
2094  double A = v[1] - o.h * .5;
2095  double B = -v[2] + o.w * .5;
2096  double to_return = (Rectangle_SolidAngle(2. * (A + o.h), 2. * (o.w - B), v[0]) -
2097  Rectangle_SolidAngle(2. * A, 2. * (o.w - B), v[0]) +
2098  Rectangle_SolidAngle(2. * (A + o.h), 2. * B, v[0]) -
2099  Rectangle_SolidAngle(2. * A, 2. * B, v[0])) *
2100  .25;
2101  return to_return;
2102  }
2103  if ((v[1] <= o.h * .5) && isDefinitelyGreaterThan(v[2], o.w * .5)) {
2104  double A = -v[1] + o.h * .5;
2105  double B = v[2] - o.w * .5;
2106  double to_return = (Rectangle_SolidAngle(2. * (o.h - A), 2. * (B + o.w), v[0]) -
2107  Rectangle_SolidAngle(2. * (o.h - A), 2. * B, v[0]) +
2108  Rectangle_SolidAngle(2. * A, 2. * (B + o.w), v[0]) -
2109  Rectangle_SolidAngle(2. * A, 2. * B, v[0])) *
2110  .25;
2111  return to_return;
2112  }
2113  // error message if none of these cases, i.e. something has gone wrong!
2114  // std::cout << "Warning: invalid solid angle call." << std::endl;
2115  return 0.;
2116  }
2117 
2118  double OpFastScintillation::Omega_Dome_Model(const double distance, const double theta) const
2119  {
2120  // this function calculates the solid angle of a semi-sphere of radius b,
2121  // as a correction to the analytic formula of the on-axix solid angle,
2122  // as we move off-axis an angle theta. We have used 9-angular bins
2123  // with delta_theta width.
2124 
2125  // par0 = Radius correction close
2126  // par1 = Radius correction far
2127  // par2 = breaking distance betwween "close" and "far"
2128 
2129  double par0[9] = {0., 0., 0., 0., 0., 0.597542, 1.00872, 1.46993, 2.04221};
2130  double par1[9] = {0, 0, 0.19569, 0.300449, 0.555598, 0.854939, 1.39166, 2.19141, 2.57732};
2131  const double delta_theta = 10.;
2132  int j = int(theta / delta_theta);
2133  // PMT radius
2134  const double b = fradius; // cm
2135  // distance form which the model parameters break (empirical value)
2136  const double d_break = 5 * b; //par2
2137 
2138  if (distance >= d_break) {
2139  double R_apparent_far = b - par1[j];
2140  return (2 * CLHEP::pi * (1 - std::sqrt(1 - std::pow(R_apparent_far / distance, 2))));
2141  }
2142  else {
2143  double R_apparent_close = b - par0[j];
2144  return (2 * CLHEP::pi * (1 - std::sqrt(1 - std::pow(R_apparent_close / distance, 2))));
2145  }
2146  }
2147 
2148  // ---------------------------------------------------------------------------
2149  std::vector<geo::BoxBoundedGeo> OpFastScintillation::extractActiveVolumes(
2150  geo::GeometryCore const& geom)
2151  {
2152  std::vector<geo::BoxBoundedGeo> activeVolumes;
2153  activeVolumes.reserve(geom.Ncryostats());
2154 
2155  for (geo::CryostatGeo const& cryo : geom.Iterate<geo::CryostatGeo>()) {
2156 
2157  // can't use it default-constructed since it would always include origin
2158  geo::BoxBoundedGeo box{cryo.TPC(0).ActiveBoundingBox()};
2159 
2160  for (geo::TPCGeo const& TPC : cryo.IterateTPCs())
2161  box.ExtendToInclude(TPC.ActiveBoundingBox());
2162 
2163  activeVolumes.push_back(std::move(box));
2164 
2165  } // for cryostats
2166 
2167  return activeVolumes;
2168  } // OpFastScintillation::extractActiveVolumes()
2169 
2170  // ---------------------------------------------------------------------------
2171 
2172  double fast_acos(double x)
2173  {
2174  double negate = double(x < 0);
2175  x = std::abs(x);
2176  x -= double(x > 1.0) * (x - 1.0); // <- equivalent to min(1.0,x), but faster
2177  double ret = -0.0187293;
2178  ret = ret * x;
2179  ret = ret + 0.0742610;
2180  ret = ret * x;
2181  ret = ret - 0.2121144;
2182  ret = ret * x;
2183  ret = ret + 1.5707288;
2184  ret = ret * std::sqrt(1.0 - x);
2185  ret = ret - 2. * negate * ret;
2186  return negate * 3.14159265358979 + ret;
2187  }
2188 
2189 }
Float_t x
Definition: compare.C:6
TRandom r
Definition: spectrum.C:23
code to link reconstructed objects back to the MC truth information
void LoadTimingsForVISPar(std::vector< double > &distances, std::vector< double > &radial_distances, std::vector< std::vector< std::vector< double >>> &cut_off, std::vector< std::vector< std::vector< double >>> &tau, double &vis_vmean, double &angle_bin_timing_vis) const
details::range_type< T > Iterate() const
Initializes the specified ID with the ID of the first cryostat.
Definition: GeometryCore.h:541
Store parameters for running LArG4.
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
Point_t GetCathodeCenter() const
Definition: TPCGeo.h:254
std::vector< std::vector< std::vector< double > > > fvispars_dome
Specializations of geo_vectors_utils.h for ROOT old vector types.
void detectedDirectHits(std::map< size_t, int > &DetectedNum, const double Num, geo::Point_t const &ScintPoint) const
OpDetGeo const & OpDetGeoFromOpDet(unsigned int OpDet) const
Returns the geo::OpDetGeo object for the given detector number.
void detectedReflecHits(std::map< size_t, int > &ReflDetectedNum, const double Num, geo::Point_t const &ScintPoint) const
void average_position(G4Step const &aStep, double *xzyPos) const
Utilities related to art service access.
bool isScintInActiveVolume(geo::Point_t const &ScintPoint)
std::vector< geo::Point_t > fOpDetCenter
G4double bi_exp(const G4double t, const G4double tau1, const G4double tau2) const
Float_t y1[n_points_granero]
Definition: compare.C:5
void LoadTimingsForVUVPar(std::vector< std::vector< double >>(&v)[7], double &step_size, double &max_d, double &min_d, double &vuv_vgroup_mean, double &vuv_vgroup_max, double &inflexion_point_distance, double &angle_bin_timing_vuv) const
Definition of util::enumerate().
G4double Gaisser_Hillas(const double x, const double *par) const
std::vector< double > fOpDetHeight
static std::vector< geo::BoxBoundedGeo > extractActiveVolumes(geo::GeometryCore const &geom)
void propagationTime(std::vector< double > &arrival_time_dist, G4ThreeVector x0, const size_t OpChannel, bool Reflected=false)
bool isBar() const
Returns whether the detector shape is a bar (TGeoBBox).
Definition: OpDetGeo.h:180
Float_t y
Definition: compare.C:6
bool const fOnlyActiveVolume
Whether photon propagation is performed only from active volumes.
int VISHits(geo::Point_t const &ScintPoint, OpticalDetector const &opDet, const double cathode_hits_rec, const std::array< double, 3 > hotspot) const
std::map< double, double > tpbemission
Geometry information for a single TPC.
Definition: TPCGeo.h:36
All information of a photon entering the sensitive optical detector volume.
Definition: SimPhotons.h:60
std::vector< std::vector< std::vector< double > > > fvispars_flat
constexpr auto abs(T v)
Returns the absolute value of the argument.
bool isOpDetInSameTPC(geo::Point_t const &ScintPoint, geo::Point_t const &OpDetPoint) const
Geant4 interface.
void LoadVisParsFlat(std::vector< double > &vis_distances_x_flat, std::vector< double > &vis_distances_r_flat, std::vector< std::vector< std::vector< double >>> &vispars_flat) const
MappedFunctions_t GetTimingTF1(Point const &p) const
std::vector< double > fvis_distances_r_flat
void ProcessStep(const G4Step &step)
double finter_d(double *x, double *par)
std::vector< double > fdistances_refl
std::unique_ptr< G4PhysicsTable > theSlowIntegralTable
intermediate_table::const_iterator const_iterator
double model_far(double *x, double *par)
std::vector< std::vector< TF1 > > VUV_timing
std::vector< std::vector< double > > fparameters[7]
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
Geometry information for a single cryostat.
Definition: CryostatGeo.h:43
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
Definition: GeometryCore.h:430
double model_close(double *x, double *par)
void LoadGHDome(std::vector< std::vector< double >> &GHvuvpars_dome, std::vector< double > &border_corr_angulo_dome, std::vector< std::vector< double >> &border_corr_dome) const
std::vector< std::vector< double > > fGHvuvpars_flat
double finter_r(double *x, double *par)
std::vector< std::vector< double > > VUV_min
G4double scint_time(const G4Step &aStep, G4double ScintillationTime, G4double ScintillationRiseTime) const
Energy deposited on a readout Optical Detector by simulated tracks.
void getVISTimes(std::vector< double > &arrivalTimes, const TVector3 &ScintPoint, const TVector3 &OpDetPoint)
geo::Point_t InitialPosition
Scintillation position in world coordinates [cm].
Definition: SimPhotons.h:63
OpFastScintillation(const G4String &processName="Scintillation", G4ProcessType type=fElectromagnetic)
std::vector< std::vector< double > > fGHvuvpars_dome
auto enumerate(Iterables &&...iterables)
Range-for loop helper tracking the number of iteration.
Definition: enumerate.h:65
void LoadVUVSemiAnalyticProperties(bool &isFlatPDCorr, bool &isDomePDCorr, double &delta_angulo_vuv, double &radius) const
Float_t y2[n_points_geant4]
Definition: compare.C:26
int VUVHits(const double Nphotons_created, geo::Point_t const &ScintPoint, OpticalDetector const &opDet) const
unsigned int fillCoords(Coords &dest, Vector const &src)
Fills a coordinate array with the coordinates of a vector.
TPCGeo const & TPC(TPCID const &tpcid=tpc_zero) const
Returns the specified TPC.
Definition: GeometryCore.h:722
std::vector< double > fborder_corr_angulo_dome
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:101
virtual G4VParticleChange * AtRestDoIt(const G4Track &aTrack, const G4Step &aStep)
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
bool const bPropagate
Whether propagation of photons is enabled.
std::vector< std::vector< double > > fborder_corr_dome
std::vector< std::vector< double > > fborder_corr_flat
Simulation objects for optical detectors.
std::vector< double > fvis_distances_x_flat
void LoadGHFlat(std::vector< std::vector< double >> &GHvuvpars_flat, std::vector< double > &border_corr_angulo_flat, std::vector< std::vector< double >> &border_corr_flat) const
auto counter(T begin, T end)
Returns an object to iterate values from begin to end in a range-for loop.
Definition: counter.h:295
Use Geant4&#39;s user "hooks" to maintain a list of particles generated by Geant4.
bool usesSemiAnalyticModel() const
Returns whether the semi-analytic visibility parametrization is being used.
MappedT0s_t GetReflT0s(Point const &p) const
void LoadVisParsDome(std::vector< double > &vis_distances_x_dome, std::vector< double > &vis_distances_r_dome, std::vector< std::vector< std::vector< double >>> &vispars_dome) const
std::vector< double > fborder_corr_angulo_flat
void AddPhoton(size_t opchannel, sim::OnePhoton &&photon, bool Reflected=false)
std::vector< double > fOpDetLength
void SetScintillationByParticleType(const G4bool)
bool RecordPhotonsProduced(const G4Step &aStep, double N)
boost::math::policies::policy< boost::math::policies::promote_double< false > > noLDoublePromote
std::unique_ptr< CLHEP::RandGeneral > fTPBEm
std::vector< geo::BoxBoundedGeo > const fActiveVolumes
std::vector< std::vector< std::vector< double > > > ftau_pars
double Length() const
Definition: OpDetGeo.h:78
static IonizationAndScintillation * Instance()
Utilities to extend the interface of geometry vectors.
bool FillSimEnergyDeposits() const
Float_t d
Definition: plot.C:235
void AddOpDetBacktrackerRecord(sim::OpDetBacktrackerRecord soc, bool Reflected=false)
Test of util::counter and support utilities.
void AddLitePhoton(int opchannel, int time, int nphotons, bool Reflected=false)
Description of geometry of one entire detector.
Definition: GeometryCore.h:119
static constexpr bool isDefinitelyGreaterThan(TReal a, TReal b, TReal tolerance=std::numeric_limits< TReal >::epsilon())
MappedCounts_t GetAllVisibilities(Point const &p, bool wantReflected=false) const
G4double sample_time(const G4double tau1, const G4double tau2) const
void AddScintillationPhotons(TrackID_t trackID, timePDclock_t timePDclock, double numberPhotons, double const *xyz, double energy)
Add scintillation photons and energy to this OpticalDetector.
std::vector< double > fradial_distances_refl
Encapsulate the geometry of an optical detector.
geo::Point_t const & GetCenter() const
Definition: OpDetGeo.h:72
std::vector< double > fvis_distances_r_dome
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
static OpDetPhotonTable * Instance(bool LitePhotons=false)
void interpolate3(std::array< double, 3 > &inter, const std::vector< double > &xData, const std::vector< double > &yData1, const std::vector< double > &yData2, const std::vector< double > &yData3, double x, bool extrapolate) const
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:180
void LoadVisSemiAnalyticProperties(double &delta_angulo_vis, double &radius) const
std::vector< std::vector< std::vector< double > > > fcut_off_pars
A base class aware of world box coordinatesAn object describing a simple shape can inherit from this ...
Definition: BoxBoundedGeo.h:33
bool const fUseNhitsModel
Whether the semi-analytic model is being used for photon visibility.
A container for photon visibility mapping data.
static constexpr bool isApproximatelyEqual(TReal a, TReal b, TReal tolerance=std::numeric_limits< TReal >::epsilon())
std::vector< double > fvis_distances_x_dome
double Rectangle_SolidAngle(const double a, const double b, const double d) const
bool SetInSD
Whether the photon reaches the sensitive detector.
Definition: SimPhotons.h:84
void getVUVTimes(std::vector< double > &arrivalTimes, const double distance_in_cm, const size_t angle_bin)
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
double Omega_Dome_Model(const double distance, const double theta) const
Singleton to access a unified treatment of ionization and scintillation in LAr.
TVector3 toTVector3(Vector const &v)
Converts a vector into a TVector3.
double Disk_SolidAngle(const double d, const double h, const double b) const
double interpolate(const std::vector< double > &xData, const std::vector< double > &yData, double x, bool extrapolate, size_t i=0) const
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
void generateParam(const size_t index, const size_t angle_bin)
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
void ExtendToInclude(Coord_t x, Coord_t y, Coord_t z)
Extends the current box to also include the specified point.
constexpr T pi()
Returns the constant pi (up to 35 decimal digits of precision)
phot::PhotonVisibilityService const *const fPVS
Photon visibility service instance.
bool const fOpaqueCathode
Whether the cathodes are fully opaque; currently hard coded "no".
std::vector< std::vector< double > > VUV_max
MaybeLogger_< ELseverityLevel::ELsev_success, true > LogTrace
float Energy
Scintillation photon energy [GeV].
Definition: SimPhotons.h:78
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *)
static constexpr bool isDefinitelyLessThan(TReal a, TReal b, TReal tolerance=std::numeric_limits< TReal >::epsilon())
phot::MappedFunctions_t ParPropTimeTF1
Float_t e
Definition: plot.C:35
virtual bool ScintByParticleType() const =0
double Height() const
Definition: OpDetGeo.h:80
static constexpr bool isApproximatelyZero(TReal a, TReal tolerance=std::numeric_limits< TReal >::epsilon())
void AddEnergyDeposit(int n_photon, int n_elec, double scint_yield, double energy, float start_x, float start_y, float start_z, float end_x, float end_y, float end_z, double start_time, double end_time, int trackid, int pdgcode, int g4trackid, std::string const &vol="EMPTY")
TGeoShape const * Shape() const
Returns the geometry object as TGeoShape.
Definition: OpDetGeo.h:139
art framework interface to geometry description
G4double single_exp(const G4double t, const G4double tau2) const
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
bool UseLitePhotons() const
double LandauPlusExpoFinal(double *x, double *par)
double fast_acos(double x)
std::unique_ptr< G4PhysicsTable > theFastIntegralTable
int MotherTrackID
ID of the GEANT4 track causing the scintillation.
Definition: SimPhotons.h:81