All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
larg4::OpFastScintillation Class Reference

#include <OpFastScintillation.hh>

Inheritance diagram for larg4::OpFastScintillation:

Classes

struct  Dims
 
struct  OpticalDetector
 

Public Member Functions

 OpFastScintillation (const G4String &processName="Scintillation", G4ProcessType type=fElectromagnetic)
 
 ~OpFastScintillation ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *)
 
virtual G4VParticleChange * PostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &aTrack, const G4Step &aStep)
 
void SetTrackSecondariesFirst (const G4bool state)
 
void SetFiniteRiseTime (const G4bool state)
 
G4bool GetTrackSecondariesFirst () const
 
G4bool GetFiniteRiseTime () const
 
void SetScintillationYieldFactor (const G4double yieldfactor)
 
G4double GetScintillationYieldFactor () const
 
void SetScintillationExcitationRatio (const G4double excitationratio)
 
G4double GetScintillationExcitationRatio () const
 
G4PhysicsTable * GetFastIntegralTable () const
 
G4PhysicsTable * GetSlowIntegralTable () const
 
void AddSaturation (G4EmSaturation *sat)
 
void RemoveSaturation ()
 
G4EmSaturation * GetSaturation () const
 
void SetScintillationByParticleType (const G4bool)
 
G4bool GetScintillationByParticleType () const
 
void DumpPhysicsTable () const
 
void getVUVTimes (std::vector< double > &arrivalTimes, const double distance_in_cm, const size_t angle_bin)
 
void generateParam (const size_t index, const size_t angle_bin)
 
void getVISTimes (std::vector< double > &arrivalTimes, const TVector3 &ScintPoint, const TVector3 &OpDetPoint)
 
void detectedDirectHits (std::map< size_t, int > &DetectedNum, const double Num, geo::Point_t const &ScintPoint) const
 
void detectedReflecHits (std::map< size_t, int > &ReflDetectedNum, const double Num, geo::Point_t const &ScintPoint) const
 

Protected Member Functions

void BuildThePhysicsTable ()
 
bool RecordPhotonsProduced (const G4Step &aStep, double N)
 

Protected Attributes

std::unique_ptr< G4PhysicsTable > theSlowIntegralTable
 
std::unique_ptr< G4PhysicsTable > theFastIntegralTable
 
G4bool fTrackSecondariesFirst
 
G4bool fFiniteRiseTime
 
G4double YieldFactor
 
G4double ExcitationRatio
 
G4bool scintillationByParticleType
 

Private Member Functions

bool usesSemiAnalyticModel () const
 Returns whether the semi-analytic visibility parametrization is being used. More...
 
int VUVHits (const double Nphotons_created, geo::Point_t const &ScintPoint, OpticalDetector const &opDet) const
 
int VISHits (geo::Point_t const &ScintPoint, OpticalDetector const &opDet, const double cathode_hits_rec, const std::array< double, 3 > hotspot) const
 
G4double single_exp (const G4double t, const G4double tau2) const
 
G4double bi_exp (const G4double t, const G4double tau1, const G4double tau2) const
 
G4double scint_time (const G4Step &aStep, G4double ScintillationTime, G4double ScintillationRiseTime) const
 
void propagationTime (std::vector< double > &arrival_time_dist, G4ThreeVector x0, const size_t OpChannel, bool Reflected=false)
 
G4double sample_time (const G4double tau1, const G4double tau2) const
 
double reemission_energy () const
 
void average_position (G4Step const &aStep, double *xzyPos) const
 
double Rectangle_SolidAngle (const double a, const double b, const double d) const
 
double Rectangle_SolidAngle (Dims const &o, const std::array< double, 3 > v) const
 
double Disk_SolidAngle (const double d, const double h, const double b) const
 
double Omega_Dome_Model (const double distance, const double theta) const
 
G4double Gaisser_Hillas (const double x, const double *par) const
 
void ProcessStep (const G4Step &step)
 
bool isOpDetInSameTPC (geo::Point_t const &ScintPoint, geo::Point_t const &OpDetPoint) const
 
bool isScintInActiveVolume (geo::Point_t const &ScintPoint)
 
double interpolate (const std::vector< double > &xData, const std::vector< double > &yData, double x, bool extrapolate, size_t i=0) const
 
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
 

Static Private Member Functions

static std::vector
< geo::BoxBoundedGeo
extractActiveVolumes (geo::GeometryCore const &geom)
 

Private Attributes

std::map< double, double > tpbemission
 
std::unique_ptr
< CLHEP::RandGeneral > 
fTPBEm
 
G4EmSaturation * emSaturation
 
phot::MappedFunctions_t ParPropTimeTF1
 
phot::MappedT0s_t ReflT0s
 
double fstep_size
 
double fmin_d
 
double fmax_d
 
double fvuv_vgroup_mean
 
double fvuv_vgroup_max
 
double finflexion_point_distance
 
double fangle_bin_timing_vuv
 
std::vector< std::vector
< double > > 
fparameters [7]
 
std::vector< std::vector< TF1 > > VUV_timing
 
std::vector< std::vector
< double > > 
VUV_max
 
std::vector< std::vector
< double > > 
VUV_min
 
double fvis_vmean
 
double fangle_bin_timing_vis
 
std::vector< double > fdistances_refl
 
std::vector< double > fradial_distances_refl
 
std::vector< std::vector
< std::vector< double > > > 
fcut_off_pars
 
std::vector< std::vector
< std::vector< double > > > 
ftau_pars
 
double fdelta_angulo_vuv
 
bool fIsFlatPDCorr
 
std::vector< std::vector
< double > > 
fGHvuvpars_flat
 
std::vector< double > fborder_corr_angulo_flat
 
std::vector< std::vector
< double > > 
fborder_corr_flat
 
bool fIsDomePDCorr
 
std::vector< std::vector
< double > > 
fGHvuvpars_dome
 
std::vector< double > fborder_corr_angulo_dome
 
std::vector< std::vector
< double > > 
fborder_corr_dome
 
bool fStoreReflected
 
double fdelta_angulo_vis
 
std::vector< double > fvis_distances_x_flat
 
std::vector< double > fvis_distances_r_flat
 
std::vector< std::vector
< std::vector< double > > > 
fvispars_flat
 
std::vector< double > fvis_distances_x_dome
 
std::vector< double > fvis_distances_r_dome
 
std::vector< std::vector
< std::vector< double > > > 
fvispars_dome
 
double fplane_depth
 
double fcathode_zdimension
 
double fcathode_ydimension
 
TVector3 fcathode_centre
 
std::vector
< geo::BoxBoundedGeo > const 
fActiveVolumes
 
double fradius
 
Dims fcathode_plane
 
int fL_abs_vuv
 
std::vector< geo::Point_tfOpDetCenter
 
std::vector< int > fOpDetType
 
std::vector< double > fOpDetLength
 
std::vector< double > fOpDetHeight
 
bool const bPropagate
 Whether propagation of photons is enabled. More...
 
phot::PhotonVisibilityService
const *const 
fPVS
 Photon visibility service instance. More...
 
bool const fUseNhitsModel = false
 Whether the semi-analytic model is being used for photon visibility. More...
 
bool const fOnlyActiveVolume = false
 Whether photon propagation is performed only from active volumes. More...
 
bool const fOnlyOneCryostat = false
 
bool const fOpaqueCathode = false
 Whether the cathodes are fully opaque; currently hard coded "no". More...
 

Detailed Description

Definition at line 130 of file OpFastScintillation.hh.

Constructor & Destructor Documentation

larg4::OpFastScintillation::OpFastScintillation ( const G4String &  processName = "Scintillation",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 180 of file OpFastScintillation.cxx.

181  : G4VRestDiscreteProcess(processName, type)
182  , fActiveVolumes{extractActiveVolumes(*(lar::providerFrom<geo::Geometry>()))}
183  , bPropagate(!(art::ServiceHandle<sim::LArG4Parameters const>()->NoPhotonPropagation()))
184  , fPVS(bPropagate ? art::ServiceHandle<phot::PhotonVisibilityService const>().get() : nullptr)
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 {
234  throw art::Exception(art::errors::Configuration)
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(0, 0).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 = (fmax_d - fmin_d) / fstep_size; // for d < fmin_d, no parameterisaton, a delta function is used instead
291  size_t num_angles = std::round(90/fangle_bin_timing_vuv);
292  VUV_timing = std::vector(num_angles, std::vector(num_params, TF1()));
293 
294  // initialise vectors to contain range parameterisations sampled to in each case
295  // when using TF1->GetRandom(xmin,xmax), must be in same range otherwise sampling is regenerated, this is the slow part!
296  VUV_max = std::vector(num_angles, std::vector(num_params, 0.0));
297  VUV_min = std::vector(num_angles, std::vector(num_params, 0.0));
298 
299  // VIS time parameterisation
300  if (fPVS->StoreReflected()) {
301  // load parameters
305  ftau_pars,
306  fvis_vmean,
308  }
309  }
310  if (usesSemiAnalyticModel()) {
311  mf::LogVerbatim("OpFastScintillation")
312  << "OpFastScintillation: using semi-analytic model for number of hits";
313 
314  // LAr absorption length in cm
315  std::map<double, double> abs_length_spectrum =
316  lar::providerFrom<detinfo::LArPropertiesService>()->AbsLengthSpectrum();
317  std::vector<double> x_v, y_v;
318  for (auto elem : abs_length_spectrum) {
319  x_v.push_back(elem.first);
320  y_v.push_back(elem.second);
321  }
322  fL_abs_vuv =
323  std::round(interpolate(x_v, y_v, 9.7, false)); // 9.7 eV: peak of VUV emission spectrum
324 
325  // Load Gaisser-Hillas corrections for VUV semi-analytic hits
326  std::cout << "Loading the GH corrections" << std::endl;
330  fradius
331  );
332  if (!fIsFlatPDCorr && !fIsDomePDCorr) {
333  throw cet::exception("OpFastScintillation")
334  << "Both isFlatPDCorr and isDomePDCorr parameters are false, at least one type of parameterisation is required for the semi-analytic light simulation." << "\n";
335  }
336  if (fIsFlatPDCorr) {
340  );
341  }
342  if (fIsDomePDCorr) {
346  );
347  }
348  // cathode center coordinates required for corrections
349  fcathode_centre = geom.TPC(0, 0).GetCathodeCenter();
350  fcathode_centre[1] = fActiveVolumes[0].CenterY();
351  fcathode_centre[2] = fActiveVolumes[0].CenterZ(); // to get full cathode dimension rather than just single tpc
352 
353 
354  if (fPVS->StoreReflected()) {
355  fStoreReflected = true;
356  // Load corrections for VIS semi-anlytic hits
357  std::cout << "Loading visible light corrections" << std::endl;
359  fradius
360  );
361  if (fIsFlatPDCorr) {
365  );
366  }
367  if (fIsDomePDCorr) {
371  );
372  }
373 
374  // cathode dimensions
375  fcathode_ydimension = fActiveVolumes[0].SizeY();
376  fcathode_zdimension = fActiveVolumes[0].SizeZ();
377 
378  // set cathode plane struct for solid angle function
382  }
383  else
384  fStoreReflected = false;
385  }
386  }
387  tpbemission = lar::providerFrom<detinfo::LArPropertiesService>()->TpbEm();
388  std::vector<double> parent;
389  parent.reserve(tpbemission.size());
390  for (auto iter = tpbemission.begin(); iter != tpbemission.end(); ++iter) {
391  parent.push_back(iter->second);
392  }
393  fTPBEm = std::make_unique<CLHEP::RandGeneral>
394  (parent.data(), parent.size());
395  }
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
std::vector< std::vector< std::vector< double > > > fvispars_dome
OpDetGeo const & OpDetGeoFromOpDet(unsigned int OpDet) const
Returns the geo::OpDetGeo object for the given detector number.
std::vector< geo::Point_t > fOpDetCenter
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
std::vector< double > fOpDetHeight
static std::vector< geo::BoxBoundedGeo > extractActiveVolumes(geo::GeometryCore const &geom)
bool isBar() const
Returns whether the detector shape is a bar (TGeoBBox).
Definition: OpDetGeo.h:195
bool const fOnlyActiveVolume
Whether photon propagation is performed only from active volumes.
Point GetCathodeCenter() const
Definition: TPCGeo.h:298
std::map< double, double > tpbemission
std::vector< std::vector< std::vector< double > > > fvispars_flat
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
std::vector< double > fvis_distances_r_flat
void GetCenter(double *xyz, double localz=0.0) const
Definition: OpDetGeo.cxx:40
std::vector< double > fdistances_refl
std::vector< std::vector< TF1 > > VUV_timing
std::vector< std::vector< double > > fparameters[7]
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
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
std::vector< std::vector< double > > VUV_min
std::vector< std::vector< double > > fGHvuvpars_dome
auto enumerate(Iterables &&...iterables)
Range-for loop helper tracking the number of iteration.
Definition: enumerate.h:69
void LoadVUVSemiAnalyticProperties(bool &isFlatPDCorr, bool &isDomePDCorr, double &delta_angulo_vuv, double &radius) const
std::vector< double > fborder_corr_angulo_dome
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
bool const bPropagate
Whether propagation of photons is enabled.
std::vector< std::vector< double > > fborder_corr_dome
T abs(T value)
std::vector< std::vector< double > > fborder_corr_flat
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
BEGIN_PROLOG AbsLengthSpectrum
auto counter(T begin, T end)
Returns an object to iterate values from begin to end in a range-for loop.
Definition: counter.h:285
bool usesSemiAnalyticModel() const
Returns whether the semi-analytic visibility parametrization is being used.
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
std::vector< double > fOpDetLength
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:81
Description of geometry of one entire detector.
std::vector< double > fradial_distances_refl
std::vector< double > fvis_distances_r_dome
void LoadVisSemiAnalyticProperties(double &delta_angulo_vis, double &radius) const
std::vector< std::vector< std::vector< double > > > fcut_off_pars
bool const fUseNhitsModel
Whether the semi-analytic model is being used for photon visibility.
std::vector< double > fvis_distances_x_dome
double interpolate(const std::vector< double > &xData, const std::vector< double > &yData, double x, bool extrapolate, size_t i=0) const
phot::PhotonVisibilityService const *const fPVS
Photon visibility service instance.
std::vector< std::vector< double > > VUV_max
TPCGeo const & TPC(unsigned int const tpc=0, unsigned int const cstat=0) const
Returns the specified TPC.
virtual bool ScintByParticleType() const =0
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:184
double Height() const
Definition: OpDetGeo.h:83
TGeoShape const * Shape() const
Returns the geometry object as TGeoShape.
Definition: OpDetGeo.h:154
BEGIN_PROLOG could also be cout
larg4::OpFastScintillation::~OpFastScintillation ( )

Definition at line 400 of file OpFastScintillation.cxx.

401  {
402  if (theFastIntegralTable) theFastIntegralTable->clearAndDestroy();
403  if (theSlowIntegralTable) theSlowIntegralTable->clearAndDestroy();
404  }
std::unique_ptr< G4PhysicsTable > theSlowIntegralTable
std::unique_ptr< G4PhysicsTable > theFastIntegralTable

Member Function Documentation

void larg4::OpFastScintillation::AddSaturation ( G4EmSaturation *  sat)
inline

Definition at line 216 of file OpFastScintillation.hh.

217  {
218  emSaturation = sat;
219  }
G4VParticleChange * larg4::OpFastScintillation::AtRestDoIt ( const G4Track &  aTrack,
const G4Step &  aStep 
)
virtual

Definition at line 414 of file OpFastScintillation.cxx.

417  {
418  return OpFastScintillation::PostStepDoIt(aTrack, aStep);
419  }
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
void larg4::OpFastScintillation::average_position ( G4Step const &  aStep,
double *  xzyPos 
) const
private

Definition at line 1103 of file OpFastScintillation.cxx.

1104  {
1105  CLHEP::Hep3Vector prePoint = (aStep.GetPreStepPoint())->GetPosition();
1106  CLHEP::Hep3Vector postPoint = (aStep.GetPostStepPoint())->GetPosition();
1107  xyzPos[0] = (((prePoint.getX() + postPoint.getX()) / 2.0) / CLHEP::cm);
1108  xyzPos[1] = (((prePoint.getY() + postPoint.getY()) / 2.0) / CLHEP::cm);
1109  xyzPos[2] = (((prePoint.getZ() + postPoint.getZ()) / 2.0) / CLHEP::cm);
1110  }
G4double larg4::OpFastScintillation::bi_exp ( const G4double  t,
const G4double  tau1,
const G4double  tau2 
) const
private

Definition at line 1850 of file OpFastScintillation.cxx.

1851  { // TODO: what's up with this? ... / tau2 / tau2 ...
1852  return std::exp(-1.0 * t / tau2) * (1 - std::exp(-1.0 * t / tau1)) / tau2 / tau2 *
1853  (tau1 + tau2);
1854  }
void larg4::OpFastScintillation::BuildThePhysicsTable ( )
protected

Definition at line 871 of file OpFastScintillation.cxx.

872  {
874 
875  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
876  G4int numOfMaterials = G4Material::GetNumberOfMaterials();
877 
878  // create new physics table
880  theFastIntegralTable = std::make_unique<G4PhysicsTable>(numOfMaterials);
882  theSlowIntegralTable = std::make_unique<G4PhysicsTable>(numOfMaterials);
883 
884  // loop for materials
885  for (G4int i = 0; i < numOfMaterials; i++) {
886  G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector = new G4PhysicsOrderedFreeVector();
887  G4PhysicsOrderedFreeVector* bPhysicsOrderedFreeVector = new G4PhysicsOrderedFreeVector();
888 
889  // Retrieve vector of scintillation wavelength intensity for
890  // the material from the material's optical properties table.
891  G4Material* aMaterial = (*theMaterialTable)[i];
892 
893  G4MaterialPropertiesTable* aMaterialPropertiesTable = aMaterial->GetMaterialPropertiesTable();
894 
895  if (aMaterialPropertiesTable) {
896 
897  G4MaterialPropertyVector* theFastLightVector =
898  aMaterialPropertiesTable->GetProperty("FASTCOMPONENT");
899 
900  if (theFastLightVector) {
901  // Retrieve the first intensity point in vector
902  // of (photon energy, intensity) pairs
903  G4double currentIN = (*theFastLightVector)[0];
904  if (currentIN >= 0.0) {
905  // Create first (photon energy, Scintillation
906  // Integral pair
907  G4double currentPM = theFastLightVector->Energy(0);
908  G4double currentCII = 0.0;
909 
910  aPhysicsOrderedFreeVector->InsertValues(currentPM, currentCII);
911 
912  // Set previous values to current ones prior to loop
913  G4double prevPM = currentPM;
914  G4double prevCII = currentCII;
915  G4double prevIN = currentIN;
916 
917  // loop over all (photon energy, intensity)
918  // pairs stored for this material
919  for (size_t i = 1; i < theFastLightVector->GetVectorLength(); i++) {
920  currentPM = theFastLightVector->Energy(i);
921  currentIN = (*theFastLightVector)[i];
922 
923  currentCII = 0.5 * (prevIN + currentIN);
924 
925  currentCII = prevCII + (currentPM - prevPM) * currentCII;
926 
927  aPhysicsOrderedFreeVector->InsertValues(currentPM, currentCII);
928 
929  prevPM = currentPM;
930  prevCII = currentCII;
931  prevIN = currentIN;
932  }
933  }
934  }
935 
936  G4MaterialPropertyVector* theSlowLightVector =
937  aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT");
938  if (theSlowLightVector) {
939  // Retrieve the first intensity point in vector
940  // of (photon energy, intensity) pairs
941  G4double currentIN = (*theSlowLightVector)[0];
942  if (currentIN >= 0.0) {
943  // Create first (photon energy, Scintillation
944  // Integral pair
945  G4double currentPM = theSlowLightVector->Energy(0);
946  G4double currentCII = 0.0;
947 
948  bPhysicsOrderedFreeVector->InsertValues(currentPM, currentCII);
949 
950  // Set previous values to current ones prior to loop
951  G4double prevPM = currentPM;
952  G4double prevCII = currentCII;
953  G4double prevIN = currentIN;
954 
955  // loop over all (photon energy, intensity)
956  // pairs stored for this material
957  for (size_t i = 1; i < theSlowLightVector->GetVectorLength(); i++) {
958  currentPM = theSlowLightVector->Energy(i);
959  currentIN = (*theSlowLightVector)[i];
960 
961  currentCII = 0.5 * (prevIN + currentIN);
962 
963  currentCII = prevCII + (currentPM - prevPM) * currentCII;
964 
965  bPhysicsOrderedFreeVector->InsertValues(currentPM, currentCII);
966 
967  prevPM = currentPM;
968  prevCII = currentCII;
969  prevIN = currentIN;
970  }
971  }
972  }
973  }
974  // The scintillation integral(s) for a given material
975  // will be inserted in the table(s) according to the
976  // position of the material in the material table.
977  theFastIntegralTable->insertAt(i, aPhysicsOrderedFreeVector);
978  theSlowIntegralTable->insertAt(i, bPhysicsOrderedFreeVector);
979  }
980  }
std::unique_ptr< G4PhysicsTable > theSlowIntegralTable
std::unique_ptr< G4PhysicsTable > theFastIntegralTable
void larg4::OpFastScintillation::detectedDirectHits ( std::map< size_t, int > &  DetectedNum,
const double  Num,
geo::Point_t const &  ScintPoint 
) const

Definition at line 1500 of file OpFastScintillation.cxx.

1503  {
1504  for (size_t const OpDet : util::counter(fPVS->NOpChannels())) {
1505  if (!isOpDetInSameTPC(ScintPoint, fOpDetCenter.at(OpDet))) continue;
1506 
1507  // set detector struct for solid angle function
1508  const OpFastScintillation::OpticalDetector op{
1509  fOpDetHeight.at(OpDet), fOpDetLength.at(OpDet),
1510  fOpDetCenter.at(OpDet), fOpDetType.at(OpDet)};
1511  const int DetThis = VUVHits(Num, ScintPoint, op);
1512  if (DetThis > 0) {
1513  DetectedNum[OpDet] = DetThis;
1514  // mf::LogInfo("OpFastScintillation") << "FastScint: " <<
1515  // // it->second<<" " << Num << " " << DetThisPMT;
1516  //det_photon_ctr += DetThisPMT; // CASE-DEBUG DO NOT REMOVE THIS COMMENT
1517  }
1518  }
1519  }
std::vector< geo::Point_t > fOpDetCenter
std::vector< double > fOpDetHeight
bool isOpDetInSameTPC(geo::Point_t const &ScintPoint, geo::Point_t const &OpDetPoint) const
int VUVHits(const double Nphotons_created, geo::Point_t const &ScintPoint, OpticalDetector const &opDet) 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:285
std::vector< double > fOpDetLength
phot::PhotonVisibilityService const *const fPVS
Photon visibility service instance.
void larg4::OpFastScintillation::detectedReflecHits ( std::map< size_t, int > &  ReflDetectedNum,
const double  Num,
geo::Point_t const &  ScintPoint 
) const

Definition at line 1522 of file OpFastScintillation.cxx.

1525  {
1526  // 1). calculate total number of hits of VUV photons on
1527  // reflective foils via solid angle + Gaisser-Hillas
1528  // corrections:
1529 
1530  // set plane_depth for correct TPC:
1531  double plane_depth;
1532  if (ScintPoint.X() < 0.) { plane_depth = -fplane_depth; }
1533  else {
1534  plane_depth = fplane_depth;
1535  }
1536 
1537  // get scintpoint coords relative to centre of cathode plane
1538  std::array<double, 3> ScintPoint_relative = {std::abs(ScintPoint.X() - plane_depth),
1539  std::abs(ScintPoint.Y() - fcathode_centre[1]),
1540  std::abs(ScintPoint.Z() - fcathode_centre[2])};
1541  // calculate solid angle of cathode from the scintillation point
1542  double solid_angle_cathode = Rectangle_SolidAngle(fcathode_plane, ScintPoint_relative);
1543 
1544  // calculate distance and angle between ScintPoint and hotspot
1545  // vast majority of hits in hotspot region directly infront of scintpoint,
1546  // therefore consider attenuation for this distance and on axis GH instead of for the centre coordinate
1547  double distance_cathode = std::abs(plane_depth - ScintPoint.X());
1548  // calculate hits on cathode plane via geometric acceptance
1549  double cathode_hits_geo = std::exp(-1. * distance_cathode / fL_abs_vuv) *
1550  (solid_angle_cathode / (4. * CLHEP::pi)) * Num;
1551  // determine Gaisser-Hillas correction including border effects
1552  // use flat correction
1553  double r = std::sqrt(std::pow(ScintPoint.Y() - fcathode_centre[1], 2) + std::pow(ScintPoint.Z() - fcathode_centre[2], 2));
1554  double pars_ini[4] = {0, 0, 0, 0};
1555  double s1 = 0; double s2 = 0; double s3 = 0;
1556  if(fIsFlatPDCorr) {
1557  pars_ini[0] = fGHvuvpars_flat[0][0];
1558  pars_ini[1] = fGHvuvpars_flat[1][0];
1559  pars_ini[2] = fGHvuvpars_flat[2][0];
1560  pars_ini[3] = fGHvuvpars_flat[3][0];
1564  }
1565  else std::cout << "Error: flat optical detector VUV correction required for reflected semi-analytic hits." << std::endl;
1566 
1567  // add border correction
1568  pars_ini[0] = pars_ini[0] + s1 * r;
1569  pars_ini[1] = pars_ini[1] + s2 * r;
1570  pars_ini[2] = pars_ini[2] + s3 * r;
1571  pars_ini[3] = pars_ini[3];
1572 
1573 
1574  // calculate corrected number of hits
1575  double GH_correction = Gaisser_Hillas(distance_cathode, pars_ini);
1576  const double cathode_hits_rec = GH_correction * cathode_hits_geo;
1577 
1578  // detemine hits on each PD
1579  const std::array<double, 3> hotspot = {plane_depth, ScintPoint.Y(), ScintPoint.Z()};
1580  for (size_t const OpDet : util::counter(fPVS->NOpChannels())) {
1581  if (!isOpDetInSameTPC(ScintPoint, fOpDetCenter.at(OpDet))) continue;
1582 
1583  // set detector struct for solid angle function
1584  const OpFastScintillation::OpticalDetector op{
1585  fOpDetHeight.at(OpDet), fOpDetLength.at(OpDet),
1586  fOpDetCenter.at(OpDet), fOpDetType.at(OpDet)};
1587 
1588  int const ReflDetThis =
1589  VISHits(ScintPoint, op, cathode_hits_rec, hotspot);
1590  if (ReflDetThis > 0) { ReflDetectedNum[OpDet] = ReflDetThis; }
1591  }
1592  }
std::vector< geo::Point_t > fOpDetCenter
G4double Gaisser_Hillas(const double x, const double *par) const
std::vector< double > fOpDetHeight
int VISHits(geo::Point_t const &ScintPoint, OpticalDetector const &opDet, const double cathode_hits_rec, const std::array< double, 3 > hotspot) const
bool isOpDetInSameTPC(geo::Point_t const &ScintPoint, geo::Point_t const &OpDetPoint) const
std::vector< std::vector< double > > fGHvuvpars_flat
T abs(T value)
std::vector< std::vector< double > > fborder_corr_flat
auto counter(T begin, T end)
Returns an object to iterate values from begin to end in a range-for loop.
Definition: counter.h:285
std::vector< double > fborder_corr_angulo_flat
std::vector< double > fOpDetLength
pdgs pi
Definition: selectors.fcl:22
double Rectangle_SolidAngle(const double a, const double b, const double d) const
double interpolate(const std::vector< double > &xData, const std::vector< double > &yData, double x, bool extrapolate, size_t i=0) const
phot::PhotonVisibilityService const *const fPVS
Photon visibility service instance.
esac echo uname r
BEGIN_PROLOG could also be cout
double larg4::OpFastScintillation::Disk_SolidAngle ( const double  d,
const double  h,
const double  b 
) const
private

Definition at line 2014 of file OpFastScintillation.cxx.

2015  {
2016  if (b <= 0. || d < 0. || h <= 0.) return 0.;
2017  const double leg2 = (b + d) * (b + d);
2018  const double aa = std::sqrt(h * h / (h * h + leg2));
2019  if (isApproximatelyZero(d)) { return 2. * CLHEP::pi * (1. - aa); }
2020  double bb = 2. * std::sqrt(b * d / (h * h + leg2));
2021  double cc = 4. * b * d / leg2;
2022 
2023  if (isDefinitelyGreaterThan(d, b)) {
2024  try {
2025  return 2. * aa *
2026  (std::sqrt(1. - cc) * boost::math::ellint_3(bb, cc, noLDoublePromote()) -
2027  boost::math::ellint_1(bb, noLDoublePromote()));
2028  }
2029  catch (std::domain_error& e) {
2030  if (isApproximatelyEqual(d, b, 1e-9)) {
2031  mf::LogWarning("OpFastScintillation")
2032  << "Elliptic Integral in Disk_SolidAngle() given parameters outside domain."
2033  << "\nbb: " << bb << "\ncc: " << cc << "\nException message: " << e.what()
2034  << "\nRelax condition and carry on.";
2035  return CLHEP::pi - 2. * aa * boost::math::ellint_1(bb, noLDoublePromote());
2036  }
2037  else {
2038  mf::LogError("OpFastScintillation")
2039  << "Elliptic Integral inside Disk_SolidAngle() given parameters outside domain.\n"
2040  << "\nbb: " << bb << "\ncc: " << cc << "Exception message: " << e.what();
2041  return 0.;
2042  }
2043  }
2044  }
2045  if (isDefinitelyLessThan(d, b)) {
2046  try {
2047  return 2. * CLHEP::pi -
2048  2. * aa *
2049  (boost::math::ellint_1(bb, noLDoublePromote()) +
2050  std::sqrt(1. - cc) * boost::math::ellint_3(bb, cc, noLDoublePromote()));
2051  }
2052  catch (std::domain_error& e) {
2053  if (isApproximatelyEqual(d, b, 1e-9)) {
2054  mf::LogWarning("OpFastScintillation")
2055  << "Elliptic Integral in Disk_SolidAngle() given parameters outside domain."
2056  << "\nbb: " << bb << "\ncc: " << cc << "\nException message: " << e.what()
2057  << "\nRelax condition and carry on.";
2058  return CLHEP::pi - 2. * aa * boost::math::ellint_1(bb, noLDoublePromote());
2059  }
2060  else {
2061  mf::LogError("OpFastScintillation")
2062  << "Elliptic Integral inside Disk_SolidAngle() given parameters outside domain.\n"
2063  << "\nbb: " << bb << "\ncc: " << cc << "Exception message: " << e.what();
2064  return 0.;
2065  }
2066  }
2067  }
2068  if (isApproximatelyEqual(d, b)) {
2069  return CLHEP::pi - 2. * aa * boost::math::ellint_1(bb, noLDoublePromote());
2070  }
2071  return 0.;
2072  }
while getopts h
boost::math::policies::policy< boost::math::policies::promote_double< false > > noLDoublePromote
pdgs pi
Definition: selectors.fcl:22
do i e
void larg4::OpFastScintillation::DumpPhysicsTable ( ) const
inline

Definition at line 541 of file OpFastScintillation.hh.

542  {
543  if (theFastIntegralTable) {
544  G4int PhysicsTableSize = theFastIntegralTable->entries();
545  G4PhysicsOrderedFreeVector* v;
546  for (G4int i = 0; i < PhysicsTableSize; i++) {
547  v = (G4PhysicsOrderedFreeVector*)(*theFastIntegralTable)[i];
548  v->DumpValues();
549  }
550  }
551  if (theSlowIntegralTable) {
552  G4int PhysicsTableSize = theSlowIntegralTable->entries();
553  G4PhysicsOrderedFreeVector* v;
554  for (G4int i = 0; i < PhysicsTableSize; i++) {
555  v = (G4PhysicsOrderedFreeVector*)(*theSlowIntegralTable)[i];
556  v->DumpValues();
557  }
558  }
559  }
std::unique_ptr< G4PhysicsTable > theSlowIntegralTable
std::unique_ptr< G4PhysicsTable > theFastIntegralTable
std::vector< geo::BoxBoundedGeo > larg4::OpFastScintillation::extractActiveVolumes ( geo::GeometryCore const &  geom)
staticprivate

Definition at line 2173 of file OpFastScintillation.cxx.

2174  {
2175  std::vector<geo::BoxBoundedGeo> activeVolumes;
2176  activeVolumes.reserve(geom.Ncryostats());
2177 
2178  for (geo::CryostatGeo const& cryo : geom.IterateCryostats()) {
2179 
2180  // can't use it default-constructed since it would always include origin
2181  geo::BoxBoundedGeo box{cryo.TPC(0).ActiveBoundingBox()};
2182 
2183  for (geo::TPCGeo const& TPC : cryo.IterateTPCs())
2184  box.ExtendToInclude(TPC.ActiveBoundingBox());
2185 
2186  activeVolumes.push_back(std::move(box));
2187 
2188  } // for cryostats
2189 
2190  return activeVolumes;
2191  } // OpFastScintillation::extractActiveVolumes()
Geometry information for a single TPC.
Definition: TPCGeo.h:38
Geometry information for a single cryostat.
Definition: CryostatGeo.h:43
BEGIN_PROLOG TPC
A base class aware of world box coordinatesAn object describing a simple shape can inherit from this ...
Definition: BoxBoundedGeo.h:33
void ExtendToInclude(Coord_t x, Coord_t y, Coord_t z)
Extends the current box to also include the specified point.
G4double larg4::OpFastScintillation::Gaisser_Hillas ( const double  x,
const double *  par 
) const
private

Definition at line 1857 of file OpFastScintillation.cxx.

1858  {
1859  double X_mu_0 = par[3];
1860  double Normalization = par[0];
1861  double Diff = par[1] - X_mu_0;
1862  double Term = std::pow((x - X_mu_0) / Diff, Diff / par[2]);
1863  double Exponential = std::exp((par[1] - x) / par[2]);
1864  return (Normalization * Term * Exponential);
1865  }
process_name opflash particleana ie x
void larg4::OpFastScintillation::generateParam ( const size_t  index,
const size_t  angle_bin 
)

Definition at line 1253 of file OpFastScintillation.cxx.

1254  {
1255  // get distance
1256  double distance_in_cm = (index * fstep_size) + fmin_d;
1257 
1258  // time range
1259  const double signal_t_range = 5000.; // TODO: unhardcode
1260 
1261  // parameterisation TF1
1262  TF1 fVUVTiming;
1263 
1264  // For very short distances the time correction is just a shift
1265  double t_direct_mean = distance_in_cm / fvuv_vgroup_mean;
1266  double t_direct_min = distance_in_cm / fvuv_vgroup_max;
1267 
1268  // Defining the model function(s) describing the photon transportation timing vs distance
1269  // Getting the landau parameters from the time parametrization
1270  std::array<double, 3> pars_landau;
1271  interpolate3(pars_landau,
1272  fparameters[0][0],
1273  fparameters[2][angle_bin],
1274  fparameters[3][angle_bin],
1275  fparameters[1][angle_bin],
1276  distance_in_cm,
1277  true);
1278  // Deciding which time model to use (depends on the distance)
1279  // defining useful times for the VUV arrival time shapes
1280  if (distance_in_cm >= finflexion_point_distance) {
1281  double pars_far[4] = {t_direct_min, pars_landau[0], pars_landau[1], pars_landau[2]};
1282  // Set model: Landau
1283  fVUVTiming = TF1("fVUVTiming", model_far, 0, signal_t_range, 4);
1284  fVUVTiming.SetParameters(pars_far);
1285  }
1286  else {
1287  // Set model: Landau + Exponential
1288  fVUVTiming = TF1("fVUVTiming", model_close, 0, signal_t_range, 7);
1289  // Exponential parameters
1290  double pars_expo[2];
1291  // Getting the exponential parameters from the time parametrization
1292  pars_expo[1] = interpolate(fparameters[4][0], fparameters[5][angle_bin], distance_in_cm, true);
1293  pars_expo[0] = interpolate(fparameters[4][0], fparameters[6][angle_bin], distance_in_cm, true);
1294  pars_expo[0] *= pars_landau[2];
1295  pars_expo[0] = std::log(pars_expo[0]);
1296  // this is to find the intersection point between the two functions:
1297  TF1 fint = TF1("fint", finter_d, pars_landau[0], 4 * t_direct_mean, 5);
1298  double parsInt[5] = {
1299  pars_landau[0], pars_landau[1], pars_landau[2], pars_expo[0], pars_expo[1]};
1300  fint.SetParameters(parsInt);
1301  double t_int = fint.GetMinimumX();
1302  double minVal = fint.Eval(t_int);
1303  // the functions must intersect - output warning if they don't
1304  if (minVal > 0.015) {
1305  std::cout << "WARNING: Parametrization of VUV light discontinuous for distance = "
1306  << distance_in_cm << std::endl;
1307  std::cout << "WARNING: This shouldn't be happening " << std::endl;
1308  }
1309  double parsfinal[7] = {t_int,
1310  pars_landau[0],
1311  pars_landau[1],
1312  pars_landau[2],
1313  pars_expo[0],
1314  pars_expo[1],
1315  t_direct_min};
1316  fVUVTiming.SetParameters(parsfinal);
1317  }
1318 
1319  // set the number of points used to sample parameterisation
1320  // for shorter distances, peak is sharper so more sensitive sampling required
1321  int fsampling; // TODO: unhardcode
1322  if (distance_in_cm < 50) { fsampling = 10000; }
1323  else if (distance_in_cm < 100) {
1324  fsampling = 5000;
1325  }
1326  else {
1327  fsampling = 1000;
1328  }
1329  fVUVTiming.SetNpx(fsampling);
1330 
1331  // calculate max and min distance relevant to sample parameterisation
1332  // max
1333  const size_t nq_max = 1;
1334  double xq_max[nq_max];
1335  double yq_max[nq_max];
1336  xq_max[0] = 0.995; // include 99.5%
1337  fVUVTiming.GetQuantiles(nq_max, yq_max, xq_max);
1338  double max = yq_max[0];
1339  // min
1340  double min = t_direct_min;
1341 
1342  // store TF1 and min/max, this allows identical TF1 to be used every time sampling
1343  // the first call of GetRandom generates the timing sampling and stores it in the TF1 object, this is the slow part
1344  // all subsequent calls check if it has been generated previously and are ~100+ times quicker
1345  VUV_timing[angle_bin][index] = fVUVTiming;
1346  VUV_max[angle_bin][index] = max;
1347  VUV_min[angle_bin][index] = min;
1348  }
double finter_d(double *x, double *par)
double model_far(double *x, double *par)
std::vector< std::vector< TF1 > > VUV_timing
std::vector< std::vector< double > > fparameters[7]
double model_close(double *x, double *par)
std::vector< std::vector< double > > VUV_min
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
double interpolate(const std::vector< double > &xData, const std::vector< double > &yData, double x, bool extrapolate, size_t i=0) const
std::vector< std::vector< double > > VUV_max
BEGIN_PROLOG could also be cout
G4PhysicsTable * larg4::OpFastScintillation::GetFastIntegralTable ( ) const
inline

Definition at line 535 of file OpFastScintillation.hh.

536  {
537  return theFastIntegralTable.get();
538  }
std::unique_ptr< G4PhysicsTable > theFastIntegralTable
G4bool larg4::OpFastScintillation::GetFiniteRiseTime ( ) const
inline

Definition at line 499 of file OpFastScintillation.hh.

500  {
501  return fFiniteRiseTime;
502  }
G4double larg4::OpFastScintillation::GetMeanFreePath ( const G4Track &  aTrack,
G4double  ,
G4ForceCondition *  condition 
)

Definition at line 998 of file OpFastScintillation.cxx.

999  {
1000  *condition = StronglyForced;
1001  return DBL_MAX;
1002  }
G4double larg4::OpFastScintillation::GetMeanLifeTime ( const G4Track &  aTrack,
G4ForceCondition *  condition 
)

Definition at line 1005 of file OpFastScintillation.cxx.

1006  {
1007  *condition = Forced;
1008  return DBL_MAX;
1009  }
G4EmSaturation* larg4::OpFastScintillation::GetSaturation ( ) const
inline

Definition at line 230 of file OpFastScintillation.hh.

231  {
232  return emSaturation;
233  }
G4bool larg4::OpFastScintillation::GetScintillationByParticleType ( ) const
inline

Definition at line 241 of file OpFastScintillation.hh.

242  {
244  }
G4double larg4::OpFastScintillation::GetScintillationExcitationRatio ( ) const
inline

Definition at line 523 of file OpFastScintillation.hh.

524  {
525  return ExcitationRatio;
526  }
G4double larg4::OpFastScintillation::GetScintillationYieldFactor ( ) const
inline

Definition at line 511 of file OpFastScintillation.hh.

512  {
513  return YieldFactor;
514  }
G4PhysicsTable * larg4::OpFastScintillation::GetSlowIntegralTable ( ) const
inline

Definition at line 529 of file OpFastScintillation.hh.

530  {
531  return theSlowIntegralTable.get();
532  }
std::unique_ptr< G4PhysicsTable > theSlowIntegralTable
G4bool larg4::OpFastScintillation::GetTrackSecondariesFirst ( ) const
inline

Definition at line 493 of file OpFastScintillation.hh.

494  {
495  return fTrackSecondariesFirst;
496  }
void larg4::OpFastScintillation::getVISTimes ( std::vector< double > &  arrivalTimes,
const TVector3 &  ScintPoint,
const TVector3 &  OpDetPoint 
)

Definition at line 1375 of file OpFastScintillation.cxx.

1378  {
1379  // *************************************************************************************************
1380  // Calculation of earliest arrival times and corresponding unsmeared distribution
1381  // *************************************************************************************************
1382 
1383  // set plane_depth for correct TPC:
1384  double plane_depth;
1385  if (ScintPoint[0] < 0) { plane_depth = -fplane_depth; }
1386  else {
1387  plane_depth = fplane_depth;
1388  }
1389 
1390  // calculate point of reflection for shortest path
1391  TVector3 bounce_point(plane_depth,ScintPoint[1],ScintPoint[2]);
1392 
1393  // calculate distance travelled by VUV light and by vis light
1394  double VUVdist = (bounce_point - ScintPoint).Mag();
1395  double Visdist = (OpDetPoint - bounce_point).Mag();
1396 
1397  // calculate times taken by VUV part of path
1398  int angle_bin_vuv = 0; // on-axis by definition
1399  getVUVTimes(arrivalTimes, VUVdist, angle_bin_vuv);
1400 
1401  // add visible direct path transport time
1402  for (size_t i = 0; i < arrivalTimes.size(); ++i) {
1403  arrivalTimes[i] += Visdist / fvis_vmean;
1404  }
1405 
1406  // *************************************************************************************************
1407  // Smearing of arrival time distribution
1408  // *************************************************************************************************
1409  // calculate fastest time possible
1410  // vis part
1411  double vis_time = Visdist / fvis_vmean;
1412  // vuv part
1413  double vuv_time;
1414  if (VUVdist < fmin_d) {
1415  vuv_time = VUVdist / fvuv_vgroup_max;
1416  }
1417  else {
1418  // find index of required parameterisation
1419  const size_t index = std::round((VUVdist - fmin_d) / fstep_size);
1420  // find shortest time
1421  vuv_time = VUV_min[angle_bin_vuv][index];
1422  }
1423  // sum
1424  double fastest_time = vis_time + vuv_time;
1425 
1426  // calculate angle theta between bound_point and optical detector
1427  double cosine_theta = std::abs(OpDetPoint[0] - bounce_point[0]) / Visdist;
1428  double theta = fast_acos(cosine_theta) * 180. / CLHEP::pi;
1429 
1430  // determine smearing parameters using interpolation of generated points:
1431  // 1). tau = exponential smearing factor, varies with distance and angle
1432  // 2). cutoff = largest smeared time allowed, preventing excessively large
1433  // times caused by exponential distance to cathode
1434  double distance_cathode_plane = std::abs(plane_depth - ScintPoint[0]);
1435  // angular bin
1436  size_t theta_bin = theta / fangle_bin_timing_vis;
1437  // radial distance from centre of TPC (y,z plane)
1438  double r = std::sqrt(std::pow(ScintPoint[1] - fcathode_centre[1], 2) + std::pow(ScintPoint[2] - fcathode_centre[2], 2));
1439 
1440  // cut-off and tau
1441  // cut-off
1442  // interpolate in d_c for each r bin
1443  std::vector<double> interp_vals(fcut_off_pars[theta_bin].size(), 0.0);
1444  for (size_t i = 0; i < fcut_off_pars[theta_bin].size(); i++){
1445  interp_vals[i] = interpolate(fdistances_refl, fcut_off_pars[theta_bin][i], distance_cathode_plane, true);
1446  }
1447  // interpolate in r
1448  double cutoff = interpolate(fradial_distances_refl, interp_vals, r, true);
1449 
1450  // tau
1451  // interpolate in x for each r bin
1452  std::vector<double> interp_vals_tau(ftau_pars[theta_bin].size(), 0.0);
1453  for (size_t i = 0; i < ftau_pars[theta_bin].size(); i++){
1454  interp_vals_tau[i] = interpolate(fdistances_refl, ftau_pars[theta_bin][i], distance_cathode_plane, true);
1455  }
1456  // interpolate in r
1457  double tau = interpolate(fradial_distances_refl, interp_vals_tau, r, true);
1458 
1459  if (tau < 0) { tau = 0; } // failsafe if tau extrapolate goes wrong
1460 
1461  // apply smearing:
1462  for (size_t i = 0; i < arrivalTimes.size(); ++i) {
1463  double arrival_time_smeared;
1464  // if time is already greater than cutoff, do not apply smearing
1465  if (arrivalTimes[i] >= cutoff) { continue; }
1466  // otherwise smear
1467  else {
1468  unsigned int counter = 0;
1469  // loop until time generated is within cutoff limit
1470  // most are within single attempt, very few take more than two
1471  do {
1472  // don't attempt smearings too many times
1473  if (counter >= 10) { // TODO: unhardcode
1474  arrival_time_smeared = arrivalTimes[i]; // don't smear
1475  break;
1476  }
1477  else {
1478  // generate random number in appropriate range
1479  double x = gRandom->Uniform(0.5, 1.0); // TODO: unhardcode
1480  // apply the exponential smearing
1481  arrival_time_smeared =
1482  arrivalTimes[i] + (arrivalTimes[i] - fastest_time) * (std::pow(x, -tau) - 1);
1483  }
1484  counter++;
1485  } while (arrival_time_smeared > cutoff);
1486  }
1487  arrivalTimes[i] = arrival_time_smeared;
1488  }
1489  }
process_name opflash particleana ie x
std::size_t size(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:561
std::vector< double > fdistances_refl
std::vector< std::vector< double > > VUV_min
T abs(T value)
auto counter(T begin, T end)
Returns an object to iterate values from begin to end in a range-for loop.
Definition: counter.h:285
pdgs pi
Definition: selectors.fcl:22
std::vector< std::vector< std::vector< double > > > ftau_pars
std::vector< double > fradial_distances_refl
std::vector< std::vector< std::vector< double > > > fcut_off_pars
void getVUVTimes(std::vector< double > &arrivalTimes, const double distance_in_cm, const size_t angle_bin)
double interpolate(const std::vector< double > &xData, const std::vector< double > &yData, double x, bool extrapolate, size_t i=0) const
esac echo uname r
double fast_acos(double x)
void larg4::OpFastScintillation::getVUVTimes ( std::vector< double > &  arrivalTimes,
const double  distance_in_cm,
const size_t  angle_bin 
)

Definition at line 1352 of file OpFastScintillation.cxx.

1353  {
1354  if (distance < fmin_d) {
1355  // times are fixed shift i.e. direct path only
1356  double t_prop_correction = distance / fvuv_vgroup_mean;
1357  for (size_t i = 0; i < arrivalTimes.size(); ++i) {
1358  arrivalTimes[i] = t_prop_correction;
1359  }
1360  }
1361  else { // distance >= fmin_d
1362  // determine nearest parameterisation in discretisation
1363  int index = std::round((distance - fmin_d) / fstep_size);
1364  // check whether required parameterisation has been generated, generating if not
1365  if (VUV_timing[angle_bin][index].GetNdim() == 0) { generateParam(index, angle_bin); }
1366  // randomly sample parameterisation for each photon
1367  for (size_t i = 0; i < arrivalTimes.size(); ++i) {
1368  arrivalTimes[i] = VUV_timing[angle_bin][index].GetRandom(VUV_min[angle_bin][index], VUV_max[angle_bin][index]);
1369  }
1370  }
1371  }
std::vector< std::vector< TF1 > > VUV_timing
std::vector< std::vector< double > > VUV_min
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
void generateParam(const size_t index, const size_t angle_bin)
std::vector< std::vector< double > > VUV_max
double larg4::OpFastScintillation::interpolate ( const std::vector< double > &  xData,
const std::vector< double > &  yData,
double  x,
bool  extrapolate,
size_t  i = 0 
) const
private

Definition at line 1933 of file OpFastScintillation.cxx.

1938  {
1939  if (i == 0) {
1940  size_t size = xData.size();
1941  if (x >= xData[size - 2]) { // special case: beyond right end
1942  i = size - 2;
1943  }
1944  else {
1945  while (x > xData[i + 1])
1946  i++;
1947  }
1948  }
1949  double xL = xData[i];
1950  double xR = xData[i + 1];
1951  double yL = yData[i];
1952  double yR = yData[i + 1]; // points on either side (unless beyond ends)
1953  if (!extrapolate) { // if beyond ends of array and not extrapolating
1954  if (x < xL) return yL;
1955  if (x > xR) return yL;
1956  }
1957  const double dydx = (yR - yL) / (xR - xL); // gradient
1958  return yL + dydx * (x - xL); // linear interpolation
1959  }
process_name opflash particleana ie x
std::size_t size(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:561
void larg4::OpFastScintillation::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
private

Definition at line 1962 of file OpFastScintillation.cxx.

1969  {
1970  size_t size = xData.size();
1971  size_t i = 0; // find left end of interval for interpolation
1972  if (x >= xData[size - 2]) { // special case: beyond right end
1973  i = size - 2;
1974  }
1975  else {
1976  while (x > xData[i + 1])
1977  i++;
1978  }
1979  double xL = xData[i];
1980  double xR = xData[i + 1]; // points on either side (unless beyond ends)
1981  double yL1 = yData1[i];
1982  double yR1 = yData1[i + 1];
1983  double yL2 = yData2[i];
1984  double yR2 = yData2[i + 1];
1985  double yL3 = yData3[i];
1986  double yR3 = yData3[i + 1];
1987 
1988  if (!extrapolate) { // if beyond ends of array and not extrapolating
1989  if (x < xL) {
1990  inter[0] = yL1;
1991  inter[1] = yL2;
1992  inter[2] = yL3;
1993  return;
1994  }
1995  if (x > xR) {
1996  inter[0] = yL1;
1997  inter[1] = yL2;
1998  inter[2] = yL3;
1999  return;
2000  }
2001  }
2002  const double m = (x - xL) / (xR - xL);
2003  inter[0] = m * (yR1 - yL1) + yL1;
2004  inter[1] = m * (yR2 - yL2) + yL2;
2005  inter[2] = m * (yR3 - yL3) + yL3;
2006  }
process_name opflash particleana ie x
std::size_t size(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:561
tuple m
now if test mode generate materials, CRT shell, world, gdml header else just generate CRT shell for u...
G4bool larg4::OpFastScintillation::IsApplicable ( const G4ParticleDefinition &  aParticleType)
inlinevirtual

Definition at line 472 of file OpFastScintillation.hh.

473  {
474  if (aParticleType.GetParticleName() == "opticalphoton") return false;
475  if (aParticleType.IsShortLived()) return false;
476 
477  return true;
478  }
bool larg4::OpFastScintillation::isOpDetInSameTPC ( geo::Point_t const &  ScintPoint,
geo::Point_t const &  OpDetPoint 
) const
private

Definition at line 1823 of file OpFastScintillation.cxx.

1825  {
1826  // check optical channel is in same TPC as scintillation light, if not return 0 hits
1827  // temporary method working for SBND, uBooNE, DUNE 1x2x6; to be replaced to work in full DUNE geometry
1828  // check x coordinate has same sign or is close to zero, otherwise return 0 hits
1829  if (((ScintPoint.X() < 0.) != (OpDetPoint.X() < 0.)) &&
1830  std::abs(OpDetPoint.X()) > 10.) { // TODO: unhardcode
1831  return false;
1832  }
1833  return true;
1834  }
T abs(T value)
bool larg4::OpFastScintillation::isScintInActiveVolume ( geo::Point_t const &  ScintPoint)
private

Definition at line 1837 of file OpFastScintillation.cxx.

1838  {
1839  //semi-analytic approach only works in the active volume
1840  return fActiveVolumes[0].ContainsPosition(ScintPoint);
1841  }
std::vector< geo::BoxBoundedGeo > const fActiveVolumes
double larg4::OpFastScintillation::Omega_Dome_Model ( const double  distance,
const double  theta 
) const
private

Definition at line 2141 of file OpFastScintillation.cxx.

2141  {
2142  // this function calculates the solid angle of a semi-sphere of radius b,
2143  // as a correction to the analytic formula of the on-axix solid angle,
2144  // as we move off-axis an angle theta. We have used 9-angular bins
2145  // with delta_theta width.
2146 
2147  // par0 = Radius correction close
2148  // par1 = Radius correction far
2149  // par2 = breaking distance betwween "close" and "far"
2150 
2151  double par0[9] = {0., 0., 0., 0., 0., 0.597542, 1.00872, 1.46993, 2.04221};
2152  double par1[9] = {0, 0, 0.19569, 0.300449, 0.555598, 0.854939, 1.39166, 2.19141, 2.57732};
2153  const double delta_theta = 10.;
2154  int j = int(theta/delta_theta);
2155  // PMT radius
2156  const double b = fradius; // cm
2157  // distance form which the model parameters break (empirical value)
2158  const double d_break = 5*b; //par2
2159 
2160  if(distance >= d_break) {
2161  double R_apparent_far = b - par1[j];
2162  return (2*CLHEP::pi * (1 - std::sqrt(1 - std::pow(R_apparent_far/distance,2))));
2163 
2164  }
2165  else {
2166  double R_apparent_close = b - par0[j];
2167  return (2*CLHEP::pi * (1 - std::sqrt(1 - std::pow(R_apparent_close/distance,2))));
2168  }
2169 }
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
pdgs pi
Definition: selectors.fcl:22
G4VParticleChange * larg4::OpFastScintillation::PostStepDoIt ( const G4Track &  aTrack,
const G4Step &  aStep 
)
virtual

Definition at line 425 of file OpFastScintillation.cxx.

430  {
431  aParticleChange.Initialize(aTrack);
432  // Check that we are in a material with a properties table, if not
433  // just return
434  const G4Material* aMaterial = aTrack.GetMaterial();
435  G4MaterialPropertiesTable* aMaterialPropertiesTable = aMaterial->GetMaterialPropertiesTable();
436  if (!aMaterialPropertiesTable) return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
437 
438  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
439  G4ThreeVector x0 = pPreStepPoint->GetPosition();
440  G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
441 
442  ///////////////////////////////////////////////////////////////////////////////////
443  // This is the old G4 way - but we do things differently - Ben J, Oct Nov 2012.
444  ///////////////////////////////////////////////////////////////////////////////////
445  //
446  // if (MeanNumberOfPhotons > 10.)
447  // {
448  // G4double sigma = ResolutionScale * std::sqrt(MeanNumberOfPhotons);
449  // NumPhotons = G4int(G4RandGauss::shoot(MeanNumberOfPhotons,sigma)+0.5);
450  // }
451  // else
452  // {
453  // NumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
454  // }
455  //
456  //
457  //
458  // if (NumPhotons <= 0)
459  // {
460  // // return unchanged particle and no secondaries
461  //
462  // aParticleChange.SetNumberOfSecondaries(0);
463  //
464  // return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
465  // }
466  //
467  //
468  // aParticleChange.SetNumberOfSecondaries(NumPhotons);
469  //
470  //
471  // if (fTrackSecondariesFirst) {
472  // if (aTrack.GetTrackStatus() == fAlive )
473  // aParticleChange.ProposeTrackStatus(fSuspend);
474  // }
475  //
476  //
477  //
478  //
479  ////////////////////////////////////////////////////////////////////////////////////
480  //
481 
482  ////////////////////////////////////////////////////////////////////////////////////
483  // The fast sim way - Ben J, Nov 2012
484  ////////////////////////////////////////////////////////////////////////////////////
485  //
486  //
487  // We don't want to produce any trackable G4 secondaries
488  aParticleChange.SetNumberOfSecondaries(0);
489 
490  // Retrieve the Scintillation Integral for this material
491  // new G4PhysicsOrderedFreeVector allocated to hold CII's
492 
493  // Some explanation for later improvements to scint yield code:
494  //
495  // What does G4 do here?
496  // It produces light in 2 steps, fast (scnt=1) then slow (scnt=2)
497  //
498  // The ratio of slow photons to fast photons is related by the yieldratio
499  // parameter. G4's poisson fluctuating scheme is a bit different to ours
500  // - we should check that they are equivalent.
501  //
502  // G4 poisson fluctuates the number of initial photons then divides them
503  // with a constant factor between fast + slow, whereas we poisson
504  // fluctuate separateyly the fast and slow detection numbers.
505  //
506 
507  // get the number of photons produced from the IonizationAndScintillation
508  // singleton
510  double MeanNumberOfPhotons =
512  // double stepEnergy = larg4::IonizationAndScintillation::Instance()->VisibleEnergyDeposit()/CLHEP::MeV;
513  RecordPhotonsProduced(aStep, MeanNumberOfPhotons); //, stepEnergy);
514  if (verboseLevel > 0) {
515  G4cout << "\n Exiting from OpFastScintillation::DoIt -- NumberOfSecondaries = "
516  << aParticleChange.GetNumberOfSecondaries() << G4endl;
517  }
518  return G4VRestDiscreteProcess::PostStepDoIt(aTrack, aStep);
519  }
bool RecordPhotonsProduced(const G4Step &aStep, double N)
static IonizationAndScintillation * Instance()
void larg4::OpFastScintillation::ProcessStep ( const G4Step &  step)
private

Definition at line 522 of file OpFastScintillation.cxx.

523  {
524  if (step.GetTotalEnergyDeposit() <= 0) return;
525 
527  -1,
528  -1,
529  1.0, //scintillation yield
530  (double)(step.GetTotalEnergyDeposit() / CLHEP::MeV), //energy in MeV
531  (float)(step.GetPreStepPoint()->GetPosition().x() / CLHEP::cm),
532  (float)(step.GetPreStepPoint()->GetPosition().y() / CLHEP::cm),
533  (float)(step.GetPreStepPoint()->GetPosition().z() / CLHEP::cm),
534  (float)(step.GetPostStepPoint()->GetPosition().x() / CLHEP::cm),
535  (float)(step.GetPostStepPoint()->GetPosition().y() / CLHEP::cm),
536  (float)(step.GetPostStepPoint()->GetPosition().z() / CLHEP::cm),
537  (double)(step.GetPreStepPoint()->GetGlobalTime()),
538  (double)(step.GetPostStepPoint()->GetGlobalTime()),
539  //step.GetTrack()->GetTrackID(),
541  step.GetTrack()->GetParticleDefinition()->GetPDGEncoding(),
543  step.GetPreStepPoint()->GetPhysicalVolume()->GetName());
544  }
util::quantities::megaelectronvolt MeV
static OpDetPhotonTable * Instance(bool LitePhotons=false)
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")
void larg4::OpFastScintillation::propagationTime ( std::vector< double > &  arrival_time_dist,
G4ThreeVector  x0,
const size_t  OpChannel,
bool  Reflected = false 
)
private

Definition at line 1029 of file OpFastScintillation.cxx.

1033  {
1034  assert(fPVS);
1036  throw cet::exception("OpFastScintillation")
1037  << "Cannot have both propagation time models simultaneously.";
1038  }
1039  else if (fPVS->IncludeParPropTime() &&
1040  !(ParPropTimeTF1 && (ParPropTimeTF1[OpChannel].GetNdim() == 1))) {
1041  //Warning: TF1::GetNdim()==1 will tell us if the TF1 is really defined or it is the default one.
1042  //This will fix a segfault when using timing and interpolation.
1043  G4cout << "WARNING: Requested parameterized timing, but no function found. Not applying "
1044  "propagation time."
1045  << G4endl;
1046  }
1047  else if (fPVS->IncludeParPropTime()) {
1048  if (Reflected)
1049  throw cet::exception("OpFastScintillation")
1050  << "No parameterized propagation time for reflected light";
1051  for (size_t i = 0; i < arrival_time_dist.size(); ++i) {
1052  arrival_time_dist[i] = ParPropTimeTF1[OpChannel].GetRandom();
1053  }
1054  }
1055  else if (fPVS->IncludePropTime()) {
1056  // Get VUV photons arrival time distribution from the parametrization
1057  geo::Point_t const& opDetCenter = fOpDetCenter.at(OpChannel);
1058  if (!Reflected) {
1059  const G4ThreeVector OpDetPoint(
1060  opDetCenter.X() * CLHEP::cm, opDetCenter.Y() * CLHEP::cm, opDetCenter.Z() * CLHEP::cm);
1061  double distance_in_cm = (x0 - OpDetPoint).mag() / CLHEP::cm; // this must be in CENTIMETERS!
1062  double cosine = std::abs(x0[0]*CLHEP::cm - OpDetPoint[0]*CLHEP::cm) / distance_in_cm;
1063  double theta = fast_acos(cosine)*180./CLHEP::pi;
1064  int angle_bin = theta/fangle_bin_timing_vuv;
1065  getVUVTimes(arrival_time_dist, distance_in_cm, angle_bin); // in ns
1066  }
1067  else {
1068  TVector3 const ScintPoint(x0[0] / CLHEP::cm, x0[1] / CLHEP::cm, x0[2] / CLHEP::cm); // in cm
1069  getVISTimes(arrival_time_dist, ScintPoint, geo::vect::toTVector3(opDetCenter)); // in ns
1070  }
1071  }
1072  }
std::vector< geo::Point_t > fOpDetCenter
void getVISTimes(std::vector< double > &arrivalTimes, const TVector3 &ScintPoint, const TVector3 &OpDetPoint)
T abs(T value)
pdgs pi
Definition: selectors.fcl:22
void getVUVTimes(std::vector< double > &arrivalTimes, const double distance_in_cm, const size_t angle_bin)
TVector3 toTVector3(Vector const &v)
Converts a vector into a TVector3.
phot::PhotonVisibilityService const *const fPVS
Photon visibility service instance.
phot::MappedFunctions_t ParPropTimeTF1
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:184
double fast_acos(double x)
bool larg4::OpFastScintillation::RecordPhotonsProduced ( const G4Step &  aStep,
double  N 
)
protected

Definition at line 546 of file OpFastScintillation.cxx.

548  {
549  // make sure that whatever happens afterwards, the energy deposition is stored
550  art::ServiceHandle<sim::LArG4Parameters const> lgp;
551  if (lgp->FillSimEnergyDeposits()) ProcessStep(aStep);
552 
553  // Get the pointer to the fast scintillation table
554  OpDetPhotonTable* fst = OpDetPhotonTable::Instance();
555 
556  const G4Track* aTrack = aStep.GetTrack();
557 
558  G4StepPoint const* pPreStepPoint = aStep.GetPreStepPoint();
559  // unused G4StepPoint const* pPostStepPoint = aStep.GetPostStepPoint();
560 
561  const G4DynamicParticle* aParticle = aTrack->GetDynamicParticle();
562  const G4Material* aMaterial = aTrack->GetMaterial();
563 
564  G4int materialIndex = aMaterial->GetIndex();
565  G4int tracknumber = aStep.GetTrack()->GetTrackID();
566 
567  G4ThreeVector x0 = pPreStepPoint->GetPosition();
568  G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
569  //G4double t0 = pPreStepPoint->GetGlobalTime() - fGlobalTimeOffset;
570  G4double t0 = pPreStepPoint->GetGlobalTime();
571 
572  G4MaterialPropertiesTable* aMaterialPropertiesTable = aMaterial->GetMaterialPropertiesTable();
573 
574  G4MaterialPropertyVector* Fast_Intensity =
575  aMaterialPropertiesTable->GetProperty("FASTCOMPONENT");
576  G4MaterialPropertyVector* Slow_Intensity =
577  aMaterialPropertiesTable->GetProperty("SLOWCOMPONENT");
578 
579  if (!Fast_Intensity && !Slow_Intensity) return 1;
580 
581  if (!bPropagate) return 0;
582 
583  // Get the visibility vector for this point
584  assert(fPVS);
585 
586  G4int nscnt = 1;
587  if (Fast_Intensity && Slow_Intensity) nscnt = 2;
588 
589  double Num = 0;
590  double YieldRatio = 0;
591 
593  // The scintillation response is a function of the energy
594  // deposited by particle types.
595 
596  // Get the definition of the current particle
597  G4ParticleDefinition* pDef = aParticle->GetDefinition();
598 
599  // Obtain the G4MaterialPropertyVectory containing the
600  // scintillation light yield as a function of the deposited
601  // energy for the current particle type
602 
603  // Protons
604  if (pDef == G4Proton::ProtonDefinition()) {
605  YieldRatio = aMaterialPropertiesTable->GetConstProperty("PROTONYIELDRATIO");
606  }
607  // Muons
608  else if (pDef == G4MuonPlus::MuonPlusDefinition() ||
609  pDef == G4MuonMinus::MuonMinusDefinition()) {
610  YieldRatio = aMaterialPropertiesTable->GetConstProperty("MUONYIELDRATIO");
611  }
612  // Pions
613  else if (pDef == G4PionPlus::PionPlusDefinition() ||
614  pDef == G4PionMinus::PionMinusDefinition()) {
615  YieldRatio = aMaterialPropertiesTable->GetConstProperty("PIONYIELDRATIO");
616  }
617  // Kaons
618  else if (pDef == G4KaonPlus::KaonPlusDefinition() ||
619  pDef == G4KaonMinus::KaonMinusDefinition()) {
620  YieldRatio = aMaterialPropertiesTable->GetConstProperty("KAONYIELDRATIO");
621  }
622  // Alphas
623  else if (pDef == G4Alpha::AlphaDefinition()) {
624  YieldRatio = aMaterialPropertiesTable->GetConstProperty("ALPHAYIELDRATIO");
625  }
626  // Electrons (must also account for shell-binding energy
627  // attributed to gamma from standard PhotoElectricEffect)
628  else if (pDef == G4Electron::ElectronDefinition() || pDef == G4Gamma::GammaDefinition()) {
629  YieldRatio = aMaterialPropertiesTable->GetConstProperty("ELECTRONYIELDRATIO");
630  }
631  // Default for particles not enumerated/listed above
632  else {
633  YieldRatio = aMaterialPropertiesTable->GetConstProperty("ELECTRONYIELDRATIO");
634  }
635  // If the user has not specified yields for (p,d,t,a,carbon)
636  // then these unspecified particles will default to the
637  // electron's scintillation yield
638  if (YieldRatio == 0) {
639  YieldRatio = aMaterialPropertiesTable->GetConstProperty("ELECTRONYIELDRATIO");
640  }
641  }
642 
643  geo::Point_t const ScintPoint = {x0[0] / CLHEP::cm, x0[1] / CLHEP::cm, x0[2] / CLHEP::cm};
644  if (fOnlyActiveVolume && !isScintInActiveVolume(ScintPoint)) return 0;
645  const phot::MappedCounts_t& Visibilities = fPVS->GetAllVisibilities(ScintPoint);
646 
647  phot::MappedCounts_t ReflVisibilities;
648 
649  // Store timing information in the object for use in propagationTime method
650  if (fPVS->StoreReflected()) {
651  ReflVisibilities = fPVS->GetAllVisibilities(ScintPoint, true);
652  if (fPVS->StoreReflT0()) ReflT0s = fPVS->GetReflT0s(ScintPoint);
653  }
654  if (fPVS->IncludeParPropTime()) { ParPropTimeTF1 = fPVS->GetTimingTF1(ScintPoint); }
655 
656  /*
657  // For Kazu to debug # photons generated using csv file, by default should be commented out
658  // but don't remove as it's useful. Separated portion of code relevant to this block
659  // is labeled as "CASE-DEBUG DO NOT REMOVE THIS COMMENT"
660  // (2017 May 2nd ... "kazuhiro at nevis dot columbia dot edu")
661  //
662  static bool first_time=false;
663  if(!first_time) {
664  std::cout<<"LOGMEid,pdg,mass,ke,te,de,x,y,z,t,dr,mean_npe,gen_npe,det_npe"<<std::endl;
665  first_time=true;
666  }
667 
668  std::cout<<"LOGME"
669  <<aStep.GetTrack()->GetTrackID()<<","
670  <<aParticle->GetDefinition()->GetPDGEncoding()<<","
671  <<aParticle->GetDefinition()->GetPDGMass()/CLHEP::MeV<<","
672  <<pPreStepPoint->GetKineticEnergy()<<","
673  <<pPreStepPoint->GetTotalEnergy()<<","
674  <<aStep.GetTotalEnergyDeposit()<<","
675  <<ScintPoint<<","<<t0<<","
676  <<aStep.GetDeltaPosition().mag()<<","
677  <<MeanNumberOfPhotons<<","<<std::flush;
678 
679  double gen_photon_ctr=0;
680  double det_photon_ctr=0;
681  */
682  for (G4int scnt = 1; scnt <= nscnt; scnt++) {
683  G4double ScintillationTime = 0. * CLHEP::ns;
684  G4double ScintillationRiseTime = 0. * CLHEP::ns;
685  G4PhysicsOrderedFreeVector* ScintillationIntegral = NULL;
686  if (scnt == 1) {
687  if (nscnt == 1) {
688  if (Fast_Intensity) {
689  ScintillationTime = aMaterialPropertiesTable->GetConstProperty("FASTTIMECONSTANT");
690  if (fFiniteRiseTime) {
691  ScintillationRiseTime =
692  aMaterialPropertiesTable->GetConstProperty("FASTSCINTILLATIONRISETIME");
693  }
694  ScintillationIntegral =
695  (G4PhysicsOrderedFreeVector*)((*theFastIntegralTable)(materialIndex));
696  }
697  if (Slow_Intensity) {
698  ScintillationTime = aMaterialPropertiesTable->GetConstProperty("SLOWTIMECONSTANT");
699  if (fFiniteRiseTime) {
700  ScintillationRiseTime =
701  aMaterialPropertiesTable->GetConstProperty("SLOWSCINTILLATIONRISETIME");
702  }
703  ScintillationIntegral =
704  (G4PhysicsOrderedFreeVector*)((*theSlowIntegralTable)(materialIndex));
705  }
706  } //endif nscnt=1
707  else {
708  if (YieldRatio == 0)
709  YieldRatio = aMaterialPropertiesTable->GetConstProperty("YIELDRATIO");
710 
711  if (ExcitationRatio == 1.0) { Num = std::min(YieldRatio, 1.0) * MeanNumberOfPhotons; }
712  else {
713  Num = std::min(ExcitationRatio, 1.0) * MeanNumberOfPhotons;
714  }
715  ScintillationTime = aMaterialPropertiesTable->GetConstProperty("FASTTIMECONSTANT");
716  if (fFiniteRiseTime) {
717  ScintillationRiseTime =
718  aMaterialPropertiesTable->GetConstProperty("FASTSCINTILLATIONRISETIME");
719  }
720  ScintillationIntegral =
721  (G4PhysicsOrderedFreeVector*)((*theFastIntegralTable)(materialIndex));
722  } //endif nscnt!=1
723  } //end scnt=1
724 
725  else {
726  Num = MeanNumberOfPhotons - Num;
727  ScintillationTime = aMaterialPropertiesTable->GetConstProperty("SLOWTIMECONSTANT");
728  if (fFiniteRiseTime) {
729  ScintillationRiseTime =
730  aMaterialPropertiesTable->GetConstProperty("SLOWSCINTILLATIONRISETIME");
731  }
732  ScintillationIntegral =
733  (G4PhysicsOrderedFreeVector*)((*theSlowIntegralTable)(materialIndex));
734  }
735 
736  if (!ScintillationIntegral) continue;
737  //gen_photon_ctr += Num; // CASE-DEBUG DO NOT REMOVE THIS COMMENT
738  // Max Scintillation Integral
739  // G4double CIImax = ScintillationIntegral->GetMaxValue();
740  //std::cout << "++++++++++++" << Num << "++++++++++" << std::endl;
741 
742  // here we go: now if visibilities are invalid, we are in trouble
743  //if (!Visibilities && (fPVS->NOpChannels() > 0)) {
744  // throw cet::exception("OpFastScintillator")
745  // << "Photon library does not cover point " << ScintPoint << " cm.\n";
746  //}
747 
748  if (!Visibilities && !usesSemiAnalyticModel()) continue;
749 
750  // detected photons from direct light
751  std::map<size_t, int> DetectedNum;
752  if (Visibilities && !usesSemiAnalyticModel()) {
753  for (size_t const OpDet : util::counter(fPVS->NOpChannels())) {
754  if (fOpaqueCathode && !isOpDetInSameTPC(ScintPoint, fOpDetCenter.at(OpDet))) continue;
755  int const DetThis = std::round(G4Poisson(Visibilities[OpDet] * Num));
756  if (DetThis > 0) DetectedNum[OpDet] = DetThis;
757  }
758  }
759  else {
760  detectedDirectHits(DetectedNum, Num, ScintPoint);
761  }
762 
763  // detected photons from reflected light
764  std::map<size_t, int> ReflDetectedNum;
765  if (fPVS->StoreReflected()) {
766  if (!usesSemiAnalyticModel()) {
767  for (size_t const OpDet : util::counter(fPVS->NOpChannels())) {
768  if (fOpaqueCathode && !isOpDetInSameTPC(ScintPoint, fOpDetCenter.at(OpDet))) continue;
769  int const ReflDetThis = std::round(G4Poisson(ReflVisibilities[OpDet] * Num));
770  if (ReflDetThis > 0) ReflDetectedNum[OpDet] = ReflDetThis;
771  }
772  }
773  else {
774  detectedReflecHits(ReflDetectedNum, Num, ScintPoint);
775  }
776  }
777 
778  std::vector<double> arrival_time_dist;
779  // Now we run through each PMT figuring out num of detected photons
780  for (size_t Reflected = 0; Reflected <= 1; ++Reflected) {
781  // Only do the reflected loop if we have reflected visibilities
782  if (Reflected && !fPVS->StoreReflected()) continue;
783 
784  std::map<size_t, int>::const_iterator itstart;
785  std::map<size_t, int>::const_iterator itend;
786  if (Reflected) {
787  itstart = ReflDetectedNum.begin();
788  itend = ReflDetectedNum.end();
789  }
790  else {
791  itstart = DetectedNum.begin();
792  itend = DetectedNum.end();
793  }
794  for (auto itdetphot = itstart; itdetphot != itend; ++itdetphot) {
795  const size_t OpChannel = itdetphot->first;
796  const int NPhotons = itdetphot->second;
797 
798  // Set up the OpDetBTR information
799  sim::OpDetBacktrackerRecord tmpOpDetBTRecord(OpChannel);
800  int thisG4TrackID = ParticleListAction::GetCurrentTrackID();
801  double xyzPos[3];
802  average_position(aStep, xyzPos);
803  double Edeposited = 0;
805  //We use this when it is the only sensical information. It may be of limited use to end users.
806  Edeposited = aStep.GetTotalEnergyDeposit();
807  }
808  else if (emSaturation) {
809  //If Birk Coefficient used, log VisibleEnergies.
810  Edeposited =
812  }
813  else {
814  //We use this when it is the only sensical information. It may be of limited use to end users.
815  Edeposited = aStep.GetTotalEnergyDeposit();
816  }
817 
818  // Get the transport time distribution
819  arrival_time_dist.resize(NPhotons);
820  propagationTime(arrival_time_dist, x0, OpChannel, Reflected);
821 
822  //We need to split the energy up by the number of photons so that we never try to write a 0 energy.
823  Edeposited = Edeposited / double(NPhotons);
824 
825  // Loop through the photons
826  for (G4int i = 0; i < NPhotons; ++i) {
827  //std::cout<<"VUV time correction: "<<arrival_time_dist[i]<<std::endl;
828  G4double Time = t0 + scint_time(aStep, ScintillationTime, ScintillationRiseTime) +
829  arrival_time_dist[i] * CLHEP::ns;
830 
831  // Always store the BTR
832  tmpOpDetBTRecord.AddScintillationPhotons(thisG4TrackID, Time, 1, xyzPos, Edeposited);
833 
834  // Store as lite photon or as OnePhoton
835  if (lgp->UseLitePhotons()) {
836  fst->AddLitePhoton(OpChannel, static_cast<int>(Time), 1, Reflected);
837  }
838  else {
839  // The sim photon in this case stores its production point and time
840  TVector3 PhotonPosition(x0[0], x0[1], x0[2]);
841 
842  float PhotonEnergy = 0;
843  if (Reflected)
844  PhotonEnergy = reemission_energy() * CLHEP::eV;
845  else
846  PhotonEnergy = 9.7 * CLHEP::eV; // 9.7 eV peak of VUV emission spectrum
847 
848  // Make a photon object for the collection
849  sim::OnePhoton PhotToAdd;
850  PhotToAdd.InitialPosition = PhotonPosition;
851  PhotToAdd.Energy = PhotonEnergy;
852  PhotToAdd.Time = Time;
853  PhotToAdd.SetInSD = false;
854  PhotToAdd.MotherTrackID = tracknumber;
855 
856  fst->AddPhoton(OpChannel, std::move(PhotToAdd), Reflected);
857  }
858  }
859  fst->AddOpDetBacktrackerRecord(tmpOpDetBTRecord, Reflected);
860  }
861  }
862  }
863  //std::cout<<gen_photon_ctr<<","<<det_photon_ctr<<std::endl; // CASE-DEBUG DO NOT REMOVE THIS COMMENT
864  return 0;
865  }
void detectedDirectHits(std::map< size_t, int > &DetectedNum, const double Num, geo::Point_t const &ScintPoint) const
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
bool isScintInActiveVolume(geo::Point_t const &ScintPoint)
std::vector< geo::Point_t > fOpDetCenter
void propagationTime(std::vector< double > &arrival_time_dist, G4ThreeVector x0, const size_t OpChannel, bool Reflected=false)
bool const fOnlyActiveVolume
Whether photon propagation is performed only from active volumes.
All information of a photon entering the sensitive optical detector volume.
Definition: SimPhotons.h:64
bool isOpDetInSameTPC(geo::Point_t const &ScintPoint, geo::Point_t const &OpDetPoint) const
void ProcessStep(const G4Step &step)
std::unique_ptr< G4PhysicsTable > theSlowIntegralTable
G4double scint_time(const G4Step &aStep, G4double ScintillationTime, G4double ScintillationRiseTime) const
util::quantities::megaelectronvolt MeV
Energy deposited on a readout Optical Detector by simulated tracks.
geo::Point_t InitialPosition
Scintillation position in world coordinates [cm].
Definition: SimPhotons.h:67
bool const bPropagate
Whether propagation of photons is enabled.
auto counter(T begin, T end)
Returns an object to iterate values from begin to end in a range-for loop.
Definition: counter.h:285
bool usesSemiAnalyticModel() const
Returns whether the semi-analytic visibility parametrization is being used.
static IonizationAndScintillation * Instance()
MappedCounts_t GetAllVisibilities(Point const &p, bool wantReflected=false) const
static OpDetPhotonTable * Instance(bool LitePhotons=false)
A container for photon visibility mapping data.
bool SetInSD
Whether the photon reaches the sensitive detector.
Definition: SimPhotons.h:88
phot::PhotonVisibilityService const *const fPVS
Photon visibility service instance.
bool const fOpaqueCathode
Whether the cathodes are fully opaque; currently hard coded &quot;no&quot;.
float Energy
Scintillation photon energy [GeV].
Definition: SimPhotons.h:82
phot::MappedFunctions_t ParPropTimeTF1
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:184
std::unique_ptr< G4PhysicsTable > theFastIntegralTable
int MotherTrackID
ID of the GEANT4 track causing the scintillation.
Definition: SimPhotons.h:85
double larg4::OpFastScintillation::Rectangle_SolidAngle ( const double  a,
const double  b,
const double  d 
) const
private

Definition at line 2076 of file OpFastScintillation.cxx.

2077  {
2078  double aa = a / (2. * d);
2079  double bb = b / (2. * d);
2080  double aux = (1. + aa * aa + bb * bb) / ((1. + aa * aa) * (1. + bb * bb));
2081  // return 4 * std::acos(std::sqrt(aux));
2082  return 4. * fast_acos(std::sqrt(aux));
2083  }
process_name gaushit a
double fast_acos(double x)
double larg4::OpFastScintillation::Rectangle_SolidAngle ( Dims const &  o,
const std::array< double, 3 >  v 
) const
private

Definition at line 2087 of file OpFastScintillation.cxx.

2088  {
2089  // v is the position of the track segment with respect to
2090  // the center position of the arapuca window
2091 
2092  // arapuca plane fixed in x direction
2093  if (isApproximatelyZero(v[1]) && isApproximatelyZero(v[2])) {
2094  return Rectangle_SolidAngle(o.h, o.w, v[0]);
2095  }
2096  if (isDefinitelyGreaterThan(v[1], o.h * .5) && isDefinitelyGreaterThan(v[2], o.w * .5)) {
2097  double A = v[1] - o.h * .5;
2098  double B = v[2] - o.w * .5;
2099  double to_return = (Rectangle_SolidAngle(2. * (A + o.h), 2. * (B + o.w), v[0]) -
2100  Rectangle_SolidAngle(2. * A, 2. * (B + o.w), v[0]) -
2101  Rectangle_SolidAngle(2. * (A + o.h), 2. * B, v[0]) +
2102  Rectangle_SolidAngle(2. * A, 2. * B, v[0])) *
2103  .25;
2104  return to_return;
2105  }
2106  if ((v[1] <= o.h * .5) && (v[2] <= o.w * .5)) {
2107  double A = -v[1] + o.h * .5;
2108  double B = -v[2] + o.w * .5;
2109  double to_return = (Rectangle_SolidAngle(2. * (o.h - A), 2. * (o.w - B), v[0]) +
2110  Rectangle_SolidAngle(2. * A, 2. * (o.w - B), v[0]) +
2111  Rectangle_SolidAngle(2. * (o.h - A), 2. * B, v[0]) +
2112  Rectangle_SolidAngle(2. * A, 2. * B, v[0])) *
2113  .25;
2114  return to_return;
2115  }
2116  if (isDefinitelyGreaterThan(v[1], o.h * .5) && (v[2] <= o.w * .5)) {
2117  double A = v[1] - o.h * .5;
2118  double B = -v[2] + o.w * .5;
2119  double to_return = (Rectangle_SolidAngle(2. * (A + o.h), 2. * (o.w - B), v[0]) -
2120  Rectangle_SolidAngle(2. * A, 2. * (o.w - B), v[0]) +
2121  Rectangle_SolidAngle(2. * (A + o.h), 2. * B, v[0]) -
2122  Rectangle_SolidAngle(2. * A, 2. * B, v[0])) *
2123  .25;
2124  return to_return;
2125  }
2126  if ((v[1] <= o.h * .5) && isDefinitelyGreaterThan(v[2], o.w * .5)) {
2127  double A = -v[1] + o.h * .5;
2128  double B = v[2] - o.w * .5;
2129  double to_return = (Rectangle_SolidAngle(2. * (o.h - A), 2. * (B + o.w), v[0]) -
2130  Rectangle_SolidAngle(2. * (o.h - A), 2. * B, v[0]) +
2131  Rectangle_SolidAngle(2. * A, 2. * (B + o.w), v[0]) -
2132  Rectangle_SolidAngle(2. * A, 2. * B, v[0])) *
2133  .25;
2134  return to_return;
2135  }
2136  // error message if none of these cases, i.e. something has gone wrong!
2137  // std::cout << "Warning: invalid solid angle call." << std::endl;
2138  return 0.;
2139  }
double Rectangle_SolidAngle(const double a, const double b, const double d) const
float A
Definition: dedx.py:137
double larg4::OpFastScintillation::reemission_energy ( ) const
private

Definition at line 1096 of file OpFastScintillation.cxx.

1097  {
1098  return fTPBEm->fire() * ((*(--tpbemission.end())).first - (*tpbemission.begin()).first) +
1099  (*tpbemission.begin()).first;
1100  }
std::map< double, double > tpbemission
std::unique_ptr< CLHEP::RandGeneral > fTPBEm
void larg4::OpFastScintillation::RemoveSaturation ( )
inline

Definition at line 223 of file OpFastScintillation.hh.

224  {
225  emSaturation = NULL;
226  }
G4double larg4::OpFastScintillation::sample_time ( const G4double  tau1,
const G4double  tau2 
) const
private

Definition at line 1075 of file OpFastScintillation.cxx.

1076  {
1077  // tau1: rise time and tau2: decay time
1078  while (1) {
1079  // two random numbers
1080  G4double ran1 = G4UniformRand();
1081  G4double ran2 = G4UniformRand();
1082  //
1083  // exponential distribution as envelope function: very efficient
1084  //
1085  G4double d = (tau1 + tau2) / tau2;
1086  // make sure the envelope function is
1087  // always larger than the bi-exponential
1088  G4double t = -1.0 * tau2 * std::log(1 - ran1);
1089  G4double g = d * single_exp(t, tau2);
1090  if (ran2 <= bi_exp(t, tau1, tau2) / g) return t;
1091  }
1092  return -1.0;
1093  }
G4double bi_exp(const G4double t, const G4double tau1, const G4double tau2) const
BEGIN_PROLOG g
G4double single_exp(const G4double t, const G4double tau2) const
G4double larg4::OpFastScintillation::scint_time ( const G4Step &  aStep,
G4double  ScintillationTime,
G4double  ScintillationRiseTime 
) const
private

Definition at line 1012 of file OpFastScintillation.cxx.

1015  {
1016  G4StepPoint const* pPreStepPoint = aStep.GetPreStepPoint();
1017  G4StepPoint const* pPostStepPoint = aStep.GetPostStepPoint();
1018  G4double avgVelocity = (pPreStepPoint->GetVelocity() + pPostStepPoint->GetVelocity()) / 2.;
1019  G4double deltaTime = aStep.GetStepLength() / avgVelocity;
1020  if (ScintillationRiseTime == 0.0) {
1021  deltaTime = deltaTime - ScintillationTime * std::log(G4UniformRand());
1022  }
1023  else {
1024  deltaTime = deltaTime + sample_time(ScintillationRiseTime, ScintillationTime);
1025  }
1026  return deltaTime;
1027  }
G4double sample_time(const G4double tau1, const G4double tau2) const
void larg4::OpFastScintillation::SetFiniteRiseTime ( const G4bool  state)
inline

Definition at line 487 of file OpFastScintillation.hh.

void larg4::OpFastScintillation::SetScintillationByParticleType ( const G4bool  scintType)

Definition at line 985 of file OpFastScintillation.cxx.

986  {
987  if (emSaturation) {
988  G4Exception("OpFastScintillation::SetScintillationByParticleType",
989  "Scint02",
990  JustWarning,
991  "Redefinition: Birks Saturation is replaced by ScintillationByParticleType!");
993  }
994  scintillationByParticleType = scintType;
995  }
void larg4::OpFastScintillation::SetScintillationExcitationRatio ( const G4double  excitationratio)
inline

Definition at line 517 of file OpFastScintillation.hh.

518  {
519  ExcitationRatio = excitationratio;
520  }
void larg4::OpFastScintillation::SetScintillationYieldFactor ( const G4double  yieldfactor)
inline

Definition at line 505 of file OpFastScintillation.hh.

506  {
507  YieldFactor = yieldfactor;
508  }
void larg4::OpFastScintillation::SetTrackSecondariesFirst ( const G4bool  state)
inline

Definition at line 481 of file OpFastScintillation.hh.

G4double larg4::OpFastScintillation::single_exp ( const G4double  t,
const G4double  tau2 
) const
private

Definition at line 1844 of file OpFastScintillation.cxx.

1845  {
1846  return std::exp(-1.0 * t / tau2) / tau2;
1847  }
bool larg4::OpFastScintillation::usesSemiAnalyticModel ( ) const
private

Returns whether the semi-analytic visibility parametrization is being used.

Definition at line 1493 of file OpFastScintillation.cxx.

1494  {
1495  return fUseNhitsModel;
1496  } // OpFastScintillation::usesSemiAnalyticModel()
bool const fUseNhitsModel
Whether the semi-analytic model is being used for photon visibility.
int larg4::OpFastScintillation::VISHits ( geo::Point_t const &  ScintPoint,
OpticalDetector const &  opDet,
const double  cathode_hits_rec,
const std::array< double, 3 >  hotspot 
) const
private

Definition at line 1692 of file OpFastScintillation.cxx.

1696  {
1697  // 1). calculate total number of hits of VUV photons on reflective
1698  // foils via solid angle + Gaisser-Hillas corrections.
1699  // Done outside as it doesn't depend on OpDetPoint
1700 
1701  // set plane_depth for correct TPC:
1702  double plane_depth;
1703  if (ScintPoint_v.X() < 0.) { plane_depth = -fplane_depth; }
1704  else {
1705  plane_depth = fplane_depth;
1706  }
1707 
1708  // 2). calculate number of these hits which reach the optical
1709  // detector from the hotspot via solid angle
1710 
1711  // the interface has been converted into geo::Point_t, the implementation not yet
1712  std::array<double, 3U> ScintPoint;
1713  std::array<double, 3U> OpDetPoint;
1714  geo::vect::fillCoords(ScintPoint, ScintPoint_v);
1715  geo::vect::fillCoords(OpDetPoint, opDet.OpDetPoint);
1716 
1717  // calculate distances and angles for application of corrections
1718  // distance from hotspot to optical detector
1719  double distance_vis = dist(&hotspot[0], &OpDetPoint[0], 3);
1720  // angle between hotspot and optical detector
1721  double cosine_vis = std::abs(hotspot[0] - OpDetPoint[0]) / distance_vis;
1722  // double theta_vis = std::acos(cosine_vis) * 180. / CLHEP::pi;
1723  double theta_vis = fast_acos(cosine_vis) * 180. / CLHEP::pi;
1724 
1725  // calculate solid angle of optical detector
1726  double solid_angle_detector = 0;
1727  // rectangular aperture
1728  if (opDet.type == 0) {
1729  // get hotspot coordinates relative to opDet
1730  std::array<double, 3> emission_relative = {std::abs(hotspot[0] - OpDetPoint[0]),
1731  std::abs(hotspot[1] - OpDetPoint[1]),
1732  std::abs(hotspot[2] - OpDetPoint[2])};
1733  // calculate solid angle
1734  solid_angle_detector = Rectangle_SolidAngle(Dims{opDet.h, opDet.w}, emission_relative);
1735  }
1736  // dome aperture
1737  else if (opDet.type == 1) {
1738  solid_angle_detector = Omega_Dome_Model(distance_vis, theta_vis);
1739  }
1740  // disk aperture
1741  else if (opDet.type == 2) {
1742  // offset in z-y plane
1743  double d = dist(&hotspot[1], &OpDetPoint[1], 2);
1744  // drift distance (in x)
1745  double h = std::abs(hotspot[0] - OpDetPoint[0]);
1746  // calculate solid angle
1747  solid_angle_detector = Disk_SolidAngle(d, h, fradius);
1748  }
1749  else {
1750  std::cout << "Error: Invalid optical detector type. 0 = rectangular, 1 = dome, 2 = disk" << std::endl;
1751  }
1752 
1753  // calculate number of hits via geometeric acceptance
1754  double hits_geo = (solid_angle_detector / (2. * CLHEP::pi)) * cathode_hits_rec; // 2*pi due to presence of reflective foils
1755 
1756  // determine correction factor, depending on PD type
1757  const size_t k = (theta_vis / fdelta_angulo_vis); // off-set angle bin
1758  double r = dist(ScintPoint, fcathode_centre, 2, 1); // radial distance from centre of detector (Y-Z)
1759  double d_c = std::abs(ScintPoint[0] - plane_depth); // distance to cathode
1760  double border_correction = 0;
1761  // flat PDs
1762  if ((opDet.type == 0 || opDet.type == 2) && fIsFlatPDCorr){
1763  // interpolate in d_c for each r bin
1764  const size_t nbins_r = fvispars_flat[k].size();
1765  std::vector<double> interp_vals(nbins_r, 0.0);
1766  {
1767  size_t idx = 0;
1768  size_t size = fvis_distances_x_flat.size();
1769  if (d_c >= fvis_distances_x_flat[size - 2])
1770  idx = size - 2;
1771  else {
1772  while (d_c > fvis_distances_x_flat[idx + 1])
1773  idx++;
1774  }
1775  for (size_t i = 0; i < nbins_r; ++i) {
1776  interp_vals[i] = interpolate(fvis_distances_x_flat,
1777  fvispars_flat[k][i],
1778  d_c,
1779  false,
1780  idx);
1781  }
1782  }
1783  // interpolate in r
1784  border_correction = interpolate(fvis_distances_r_flat, interp_vals, r, false);
1785  }
1786  // dome PDs
1787  else if (opDet.type == 1 && fIsDomePDCorr) {
1788  // interpolate in d_c for each r bin
1789  const size_t nbins_r = fvispars_dome[k].size();
1790  std::vector<double> interp_vals(nbins_r, 0.0);
1791  {
1792  size_t idx = 0;
1793  size_t size = fvis_distances_x_dome.size();
1794  if (d_c >= fvis_distances_x_dome[size - 2])
1795  idx = size - 2;
1796  else {
1797  while (d_c > fvis_distances_x_dome[idx + 1])
1798  idx++;
1799  }
1800  for (size_t i = 0; i < nbins_r; ++i) {
1801  interp_vals[i] = interpolate(fvis_distances_x_dome,
1802  fvispars_dome[k][i],
1803  d_c,
1804  false,
1805  idx);
1806  }
1807  }
1808  // interpolate in r
1809  border_correction = interpolate(fvis_distances_r_dome, interp_vals, r, false);
1810  }
1811  else {
1812  std::cout << "Error: Invalid optical detector type. 0 = rectangular, 1 = dome, 2 = disk. Or corrections for chosen optical detector type missing." << std::endl;
1813  }
1814 
1815  // apply correction factor
1816  double hits_rec = border_correction * hits_geo / cosine_vis;
1817  double hits_vis = std::round(G4Poisson(hits_rec));
1818 
1819  return hits_vis;
1820  }
std::vector< std::vector< std::vector< double > > > fvispars_dome
std::vector< std::vector< std::vector< double > > > fvispars_flat
std::vector< double > fvis_distances_r_flat
std::size_t size(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:561
unsigned int fillCoords(Coords &dest, Vector const &src)
Fills a coordinate array with the coordinates of a vector.
while getopts h
T abs(T value)
std::vector< double > fvis_distances_x_flat
pdgs pi
Definition: selectors.fcl:22
std::vector< double > fvis_distances_r_dome
std::vector< double > fvis_distances_x_dome
double Rectangle_SolidAngle(const double a, const double b, const double d) const
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
double Omega_Dome_Model(const double distance, const double theta) const
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
pdgs k
Definition: selectors.fcl:22
esac echo uname r
BEGIN_PROLOG could also be cout
double fast_acos(double x)
int larg4::OpFastScintillation::VUVHits ( const double  Nphotons_created,
geo::Point_t const &  ScintPoint,
OpticalDetector const &  opDet 
) const
private

Definition at line 1596 of file OpFastScintillation.cxx.

1599  {
1600  // the interface has been converted into geo::Point_t, the implementation not yet
1601  std::array<double, 3U> ScintPoint;
1602  std::array<double, 3U> OpDetPoint;
1603  geo::vect::fillCoords(ScintPoint, ScintPoint_v);
1604  geo::vect::fillCoords(OpDetPoint, opDet.OpDetPoint);
1605 
1606  // distance and angle between ScintPoint and OpDetPoint
1607  double distance = dist(&ScintPoint[0], &OpDetPoint[0], 3);
1608  double cosine = std::abs(ScintPoint[0] - OpDetPoint[0]) / distance;
1609  // double theta = std::acos(cosine) * 180. / CLHEP::pi;
1610  double theta = fast_acos(cosine) * 180. / CLHEP::pi;
1611 
1612  // calculate solid angle:
1613  double solid_angle = 0;
1614  // Arapucas/Bars (rectangle)
1615  if (opDet.type == 0) {
1616  // get scintillation point coordinates relative to arapuca window centre
1617  std::array<double, 3> ScintPoint_rel = {std::abs(ScintPoint[0] - OpDetPoint[0]),
1618  std::abs(ScintPoint[1] - OpDetPoint[1]),
1619  std::abs(ScintPoint[2] - OpDetPoint[2])};
1620  // calculate solid angle
1621  solid_angle = Rectangle_SolidAngle(Dims{opDet.h, opDet.w}, ScintPoint_rel);
1622  }
1623  // PMTs (dome)
1624  else if (opDet.type == 1) {
1625  solid_angle = Omega_Dome_Model(distance, theta);
1626  }
1627  // PMTs (disk approximation)
1628  else if (opDet.type == 2) {
1629  // offset in z-y plane
1630  double d = dist(&ScintPoint[1], &OpDetPoint[1], 2);
1631  // drift distance (in x)
1632  double h = std::abs(ScintPoint[0] - OpDetPoint[0]);
1633  // Solid angle of a disk
1634  solid_angle = Disk_SolidAngle(d, h, fradius);
1635  }
1636  else {
1637  std::cout << "Error: Invalid optical detector type. 0 = rectangular, 1 = dome, 2 = disk" << std::endl;
1638  }
1639 
1640  // calculate number of photons hits by geometric acceptance: accounting for solid angle and LAr absorbtion length
1641  double hits_geo = std::exp(-1. * distance / fL_abs_vuv) * (solid_angle / (4 * CLHEP::pi)) * Nphotons_created;
1642 
1643  // apply Gaisser-Hillas correction for Rayleigh scattering distance
1644  // and angular dependence offset angle bin
1645  const size_t j = (theta / fdelta_angulo_vuv);
1646 
1647  // determine GH parameters, accounting for border effects
1648  // radial distance from centre of detector (Y-Z)
1649  double r = dist(ScintPoint, fcathode_centre, 2, 1);
1650  double pars_ini[4] = {0, 0, 0, 0};
1651  double s1 = 0; double s2 = 0; double s3 = 0;
1652  // flat PDs
1653  if ((opDet.type == 0 || opDet.type == 2) && fIsFlatPDCorr){
1654  pars_ini[0] = fGHvuvpars_flat[0][j];
1655  pars_ini[1] = fGHvuvpars_flat[1][j];
1656  pars_ini[2] = fGHvuvpars_flat[2][j];
1657  pars_ini[3] = fGHvuvpars_flat[3][j];
1658  s1 = interpolate( fborder_corr_angulo_flat, fborder_corr_flat[0], theta, true);
1659  s2 = interpolate( fborder_corr_angulo_flat, fborder_corr_flat[1], theta, true);
1660  s3 = interpolate( fborder_corr_angulo_flat, fborder_corr_flat[2], theta, true);
1661  }
1662  // dome PDs
1663  else if (opDet.type == 1 && fIsDomePDCorr) {
1664  pars_ini[0] = fGHvuvpars_dome[0][j];
1665  pars_ini[1] = fGHvuvpars_dome[1][j];
1666  pars_ini[2] = fGHvuvpars_dome[2][j];
1667  pars_ini[3] = fGHvuvpars_dome[3][j];
1668  s1 = interpolate( fborder_corr_angulo_dome, fborder_corr_dome[0], theta, true);
1669  s2 = interpolate( fborder_corr_angulo_dome, fborder_corr_dome[1], theta, true);
1670  s3 = interpolate( fborder_corr_angulo_dome, fborder_corr_dome[2], theta, true);
1671  }
1672  else std::cout << "Error: Invalid optical detector type. 0 = rectangular, 1 = dome, 2 = disk. Or corrections for chosen optical detector type missing." << std::endl;
1673 
1674  // add border correction
1675  pars_ini[0] = pars_ini[0] + s1 * r;
1676  pars_ini[1] = pars_ini[1] + s2 * r;
1677  pars_ini[2] = pars_ini[2] + s3 * r;
1678  pars_ini[3] = pars_ini[3];
1679 
1680  // calculate correction
1681  double GH_correction = Gaisser_Hillas(distance, pars_ini);
1682 
1683  // calculate number photons
1684  double hits_rec = GH_correction * hits_geo / cosine;
1685  int hits_vuv = std::round(G4Poisson(hits_rec));
1686 
1687  return hits_vuv;
1688  }
G4double Gaisser_Hillas(const double x, const double *par) const
std::vector< std::vector< double > > fGHvuvpars_flat
std::vector< std::vector< double > > fGHvuvpars_dome
unsigned int fillCoords(Coords &dest, Vector const &src)
Fills a coordinate array with the coordinates of a vector.
std::vector< double > fborder_corr_angulo_dome
while getopts h
std::vector< std::vector< double > > fborder_corr_dome
T abs(T value)
std::vector< std::vector< double > > fborder_corr_flat
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
std::vector< double > fborder_corr_angulo_flat
pdgs pi
Definition: selectors.fcl:22
double Rectangle_SolidAngle(const double a, const double b, const double d) const
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
double Omega_Dome_Model(const double distance, const double theta) const
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
esac echo uname r
BEGIN_PROLOG could also be cout
double fast_acos(double x)

Member Data Documentation

bool const larg4::OpFastScintillation::bPropagate
private

Whether propagation of photons is enabled.

Definition at line 425 of file OpFastScintillation.hh.

G4EmSaturation* larg4::OpFastScintillation::emSaturation
private

Definition at line 338 of file OpFastScintillation.hh.

G4double larg4::OpFastScintillation::ExcitationRatio
protected

Definition at line 290 of file OpFastScintillation.hh.

std::vector<geo::BoxBoundedGeo> const larg4::OpFastScintillation::fActiveVolumes
private

Definition at line 409 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fangle_bin_timing_vis
private

Definition at line 360 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fangle_bin_timing_vuv
private

Definition at line 351 of file OpFastScintillation.hh.

std::vector<double> larg4::OpFastScintillation::fborder_corr_angulo_dome
private

Definition at line 390 of file OpFastScintillation.hh.

std::vector<double> larg4::OpFastScintillation::fborder_corr_angulo_flat
private

Definition at line 385 of file OpFastScintillation.hh.

std::vector<std::vector<double> > larg4::OpFastScintillation::fborder_corr_dome
private

Definition at line 391 of file OpFastScintillation.hh.

std::vector<std::vector<double> > larg4::OpFastScintillation::fborder_corr_flat
private

Definition at line 386 of file OpFastScintillation.hh.

TVector3 larg4::OpFastScintillation::fcathode_centre
private

Definition at line 408 of file OpFastScintillation.hh.

Dims larg4::OpFastScintillation::fcathode_plane
private

Definition at line 415 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fcathode_ydimension
private

Definition at line 407 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fcathode_zdimension
private

Definition at line 407 of file OpFastScintillation.hh.

std::vector<std::vector<std::vector<double> > > larg4::OpFastScintillation::fcut_off_pars
private

Definition at line 363 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fdelta_angulo_vis
private

Definition at line 396 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fdelta_angulo_vuv
private

Definition at line 381 of file OpFastScintillation.hh.

std::vector<double> larg4::OpFastScintillation::fdistances_refl
private

Definition at line 361 of file OpFastScintillation.hh.

G4bool larg4::OpFastScintillation::fFiniteRiseTime
protected

Definition at line 286 of file OpFastScintillation.hh.

std::vector<std::vector<double> > larg4::OpFastScintillation::fGHvuvpars_dome
private

Definition at line 389 of file OpFastScintillation.hh.

std::vector<std::vector<double> > larg4::OpFastScintillation::fGHvuvpars_flat
private

Definition at line 384 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::finflexion_point_distance
private

Definition at line 351 of file OpFastScintillation.hh.

bool larg4::OpFastScintillation::fIsDomePDCorr
private

Definition at line 388 of file OpFastScintillation.hh.

bool larg4::OpFastScintillation::fIsFlatPDCorr
private

Definition at line 383 of file OpFastScintillation.hh.

int larg4::OpFastScintillation::fL_abs_vuv
private

Definition at line 416 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fmax_d
private

Definition at line 351 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fmin_d
private

Definition at line 351 of file OpFastScintillation.hh.

bool const larg4::OpFastScintillation::fOnlyActiveVolume = false
private

Whether photon propagation is performed only from active volumes.

Definition at line 433 of file OpFastScintillation.hh.

bool const larg4::OpFastScintillation::fOnlyOneCryostat = false
private

Allows running even if light on cryostats C:1 and higher is not supported. Currently hard coded "no"

Definition at line 436 of file OpFastScintillation.hh.

bool const larg4::OpFastScintillation::fOpaqueCathode = false
private

Whether the cathodes are fully opaque; currently hard coded "no".

Definition at line 438 of file OpFastScintillation.hh.

std::vector<geo::Point_t> larg4::OpFastScintillation::fOpDetCenter
private

Definition at line 417 of file OpFastScintillation.hh.

std::vector<double> larg4::OpFastScintillation::fOpDetHeight
private

Definition at line 420 of file OpFastScintillation.hh.

std::vector<double> larg4::OpFastScintillation::fOpDetLength
private

Definition at line 419 of file OpFastScintillation.hh.

std::vector<int> larg4::OpFastScintillation::fOpDetType
private

Definition at line 418 of file OpFastScintillation.hh.

std::vector<std::vector<double> > larg4::OpFastScintillation::fparameters[7]
private

Definition at line 352 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fplane_depth
private

Definition at line 407 of file OpFastScintillation.hh.

phot::PhotonVisibilityService const* const larg4::OpFastScintillation::fPVS
private

Photon visibility service instance.

Definition at line 428 of file OpFastScintillation.hh.

std::vector<double> larg4::OpFastScintillation::fradial_distances_refl
private

Definition at line 362 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fradius
private

Definition at line 414 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fstep_size
private

Definition at line 351 of file OpFastScintillation.hh.

bool larg4::OpFastScintillation::fStoreReflected
private

Definition at line 394 of file OpFastScintillation.hh.

std::vector<std::vector<std::vector<double> > > larg4::OpFastScintillation::ftau_pars
private

Definition at line 364 of file OpFastScintillation.hh.

std::unique_ptr<CLHEP::RandGeneral> larg4::OpFastScintillation::fTPBEm
private

Definition at line 334 of file OpFastScintillation.hh.

G4bool larg4::OpFastScintillation::fTrackSecondariesFirst
protected

Definition at line 285 of file OpFastScintillation.hh.

bool const larg4::OpFastScintillation::fUseNhitsModel = false
private

Whether the semi-analytic model is being used for photon visibility.

Definition at line 431 of file OpFastScintillation.hh.

std::vector<double> larg4::OpFastScintillation::fvis_distances_r_dome
private

Definition at line 403 of file OpFastScintillation.hh.

std::vector<double> larg4::OpFastScintillation::fvis_distances_r_flat
private

Definition at line 399 of file OpFastScintillation.hh.

std::vector<double> larg4::OpFastScintillation::fvis_distances_x_dome
private

Definition at line 402 of file OpFastScintillation.hh.

std::vector<double> larg4::OpFastScintillation::fvis_distances_x_flat
private

Definition at line 398 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fvis_vmean
private

Definition at line 360 of file OpFastScintillation.hh.

std::vector<std::vector<std::vector<double> > > larg4::OpFastScintillation::fvispars_dome
private

Definition at line 404 of file OpFastScintillation.hh.

std::vector<std::vector<std::vector<double> > > larg4::OpFastScintillation::fvispars_flat
private

Definition at line 400 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fvuv_vgroup_max
private

Definition at line 351 of file OpFastScintillation.hh.

double larg4::OpFastScintillation::fvuv_vgroup_mean
private

Definition at line 351 of file OpFastScintillation.hh.

phot::MappedFunctions_t larg4::OpFastScintillation::ParPropTimeTF1
private

Definition at line 340 of file OpFastScintillation.hh.

phot::MappedT0s_t larg4::OpFastScintillation::ReflT0s
private

Definition at line 341 of file OpFastScintillation.hh.

G4bool larg4::OpFastScintillation::scintillationByParticleType
protected

Definition at line 292 of file OpFastScintillation.hh.

std::unique_ptr<G4PhysicsTable> larg4::OpFastScintillation::theFastIntegralTable
protected

Definition at line 283 of file OpFastScintillation.hh.

std::unique_ptr<G4PhysicsTable> larg4::OpFastScintillation::theSlowIntegralTable
protected

Definition at line 282 of file OpFastScintillation.hh.

std::map<double, double> larg4::OpFastScintillation::tpbemission
private

Definition at line 333 of file OpFastScintillation.hh.

std::vector<std::vector<double> > larg4::OpFastScintillation::VUV_max
private

Definition at line 356 of file OpFastScintillation.hh.

std::vector<std::vector<double> > larg4::OpFastScintillation::VUV_min
private

Definition at line 357 of file OpFastScintillation.hh.

std::vector<std::vector<TF1> > larg4::OpFastScintillation::VUV_timing
private

Definition at line 354 of file OpFastScintillation.hh.

G4double larg4::OpFastScintillation::YieldFactor
protected

Definition at line 288 of file OpFastScintillation.hh.


The documentation for this class was generated from the following files: