LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
LArVoxelReadout.cxx
Go to the documentation of this file.
1 
8 // C/C++ standard library
9 #include <cassert>
10 #include <cmath> // std::ceil()
11 #include <cstdio> // std::sscanf()
12 #include <map>
13 #include <string>
14 #include <utility> // std::move()
15 
16 // GEANT
17 #include "Geant4/G4Step.hh"
18 #include "Geant4/G4StepPoint.hh"
19 #include "Geant4/G4ThreeVector.hh"
20 #include "Geant4/G4VPhysicalVolume.hh"
21 
22 // framework libraries
23 #include "cetlib_except/exception.h"
25 
26 // LArSoft code
28 #include "larcoreobj/SimpleTypesAndConstants/RawTypes.h" // raw::ChannelID_t
40 
41 // CLHEP
42 #include "CLHEP/Random/RandGauss.h"
43 
44 namespace larg4 {
45 
46  //---------------------------------------------------------------------------------------
47  // Constructor. Note that we force the name of this sensitive
48  // detector to be the value expected by LArVoxelListAction.
49  LArVoxelReadout::LArVoxelReadout(std::string const& name) : G4VSensitiveDetector(name)
50  {
51  // Initialize values for the electron-cluster calculation.
53 
54  // the standard name contains cryostat and TPC;
55  // if we don't find it, we will detect the TPC at each Geant hit
56  unsigned int cryostat, tpc;
57  if (std::sscanf(name.c_str(), "%*19s%1u%*4s%u", &cryostat, &tpc) == 2)
58  SetSingleTPC(cryostat, tpc);
59  else
61 
62  } // LArVoxelReadout::LArVoxelReadout()
63 
64  //---------------------------------------------------------------------------------------
65  // Constructor. Note that we force the name of this sensitive
66  // detector to be the value expected by LArVoxelListAction.
67  LArVoxelReadout::LArVoxelReadout(std::string const& name, unsigned int cryostat, unsigned int tpc)
68  : LArVoxelReadout(name)
69  {
70  SetSingleTPC(cryostat, tpc);
71  }
72 
73  //---------------------------------------------------------------------------------------
74  void LArVoxelReadout::Setup(Setup_t const& setupData)
75  {
77  SetRandomEngines(setupData.propGen);
78  }
79 
80  //---------------------------------------------------------------------------------------
81  void LArVoxelReadout::SetSingleTPC(unsigned int cryostat, unsigned int tpc)
82  {
83  bSingleTPC = true;
84  fCstat = cryostat;
85  fTPC = tpc;
86  MF_LOG_DEBUG("LArVoxelReadout") << GetName() << "covers C=" << fCstat << " T=" << fTPC;
87  }
88 
90  {
91  bSingleTPC = false;
92  fCstat = 0;
93  fTPC = 0;
94  MF_LOG_DEBUG("LArVoxelReadout") << GetName() << " autodetects TPC";
95  }
96 
97  //---------------------------------------------------------------------------------------
98  // Called at the start of each event.
100  {
101  assert(fClockData != nullptr &&
102  "You must use set the clock data pointer at the beginning "
103  "of each module's event-level call. This might be done through the "
104  "LArVoxelReadoutGeometry::Sentry class.");
105  assert(fDetProp != nullptr &&
106  "You must use set the detector-properties data pointer at the beginning "
107  "of each module's event-level call. This might be done through the "
108  "LArVoxelReadoutGeometry::Sentry class.");
109 
111  for (int i = 0; i < 3; ++i)
112  fDriftVelocity[i] =
114 
121 
122  MF_LOG_DEBUG("LArVoxelReadout")
123  << " e lifetime: " << fElectronLifetime << "\n Temperature: " << fDetProp->Temperature()
124  << "\n Drift velocity: " << fDriftVelocity[0] << " " << fDriftVelocity[1] << " "
125  << fDriftVelocity[2];
126 
128 
129  fNSteps = 0;
130  }
131 
132  //---------------------------------------------------------------------------------------
133  // Called at the end of each event.
135  {
136  MF_LOG_DEBUG("LArVoxelReadout") << "Total number of steps was " << fNSteps << std::endl;
137  }
138 
139  //---------------------------------------------------------------------------------------
141 
142  //--------------------------------------------------------------------------------------
144  {
146  size_t cryo = 0;
147  for (auto& cryoData : fChannelMaps) { // each, a vector of maps
148  cryoData.resize(fGeoHandle->NTPC(geo::CryostatID(cryo++)));
149  for (auto& channelsMap : cryoData)
150  channelsMap.clear(); // each, a map
151  } // for cryostats
152  } // LArVoxelReadout::ClearSimChannels()
153 
155  {
156  if (bSingleTPC) return GetSimChannelMap(fCstat, fTPC);
157  throw cet::exception("LArVoxelReadout") << "TPC not specified";
158  }
159 
161  {
162  if (bSingleTPC) return GetSimChannelMap(fCstat, fTPC);
163  throw cet::exception("LArVoxelReadout") << "TPC not specified";
164  }
165 
167  unsigned short tpc) const
168  {
169  return fChannelMaps.at(cryo).at(tpc);
170  }
171 
173  unsigned short tpc)
174  {
175  return fChannelMaps.at(cryo).at(tpc);
176  }
177 
178  std::vector<sim::SimChannel> LArVoxelReadout::GetSimChannels() const
179  {
180  if (bSingleTPC) return GetSimChannels(fCstat, fTPC);
181  throw cet::exception("LArVoxelReadout") << "TPC not specified";
182  }
183 
184  std::vector<sim::SimChannel> LArVoxelReadout::GetSimChannels(unsigned short cryo,
185  unsigned short tpc) const
186  {
187  std::vector<sim::SimChannel> channels;
188  const ChannelMap_t& chmap = fChannelMaps.at(cryo).at(tpc);
189  channels.reserve(chmap.size());
190  for (const auto& chpair : chmap)
191  channels.push_back(chpair.second);
192  return channels;
193  }
194 
195  //---------------------------------------------------------------------------------------
196  // Called for each step.
197  G4bool LArVoxelReadout::ProcessHits(G4Step* step, G4TouchableHistory* pHistory)
198  {
199  // All work done for the "parallel world" "box of voxels" in
200  // LArVoxelReadoutGeometry makes this a fairly simple routine.
201  // First, the usual check for non-zero energy:
202 
203  // Only process the hit if the step is inside the active volume and
204  // it has deposited energy. The hit being inside the active volume
205  // is virtually sure to happen because the LArVoxelReadoutGeometry
206  // that this class makes use of only has voxels for inside the TPC.
207 
208  // The step can be no bigger than the size of the voxel,
209  // because of the geometry set up in LArVoxelGeometry and the
210  // transportation set up in PhysicsList. Find the mid-point
211  // of the step.
212 
213  if (step->GetTotalEnergyDeposit() > 0) {
214 
215  // Make sure we have the IonizationAndScintillation singleton
216  // reset to this step
218  fNSteps++;
219  if (!fDontDriftThem) {
220 
221  G4ThreeVector midPoint =
222  0.5 * (step->GetPreStepPoint()->GetPosition() + step->GetPostStepPoint()->GetPosition());
223  double g4time = step->GetPreStepPoint()->GetGlobalTime();
224 
225  // Find the Geant4 track ID for the particle responsible for depositing the
226  // energy. if we are only storing primary EM shower particles, and this energy
227  // is from a secondary etc EM shower particle, the ID returned is the primary
228  const int trackID = ParticleListAction::GetCurrentTrackID();
229  // For all particles but shower daughters, origTrackID is the same as trackID.
230  // For shower daughters it contains their original trackID instead of the
231  // shower primary's trackID.
232  const int origTrackID = ParticleListAction::GetCurrentOrigTrackID();
233 
234  // Find out which TPC we are in.
235  // If this readout object covers just one, we already know it.
236  // Otherwise, we have to ask Geant where we are.
237  unsigned short int cryostat = 0, tpc = 0;
238  if (bSingleTPC) {
239  cryostat = fCstat;
240  tpc = fTPC;
241  }
242  else {
243  // detect the TPC we are in
244  const G4VTouchable* pTouchable = step->GetPreStepPoint()->GetTouchable();
245  if (!pTouchable) {
246  throw cet::exception("LArG4") << "Untouchable step in LArVoxelReadout::ProcessHits()";
247  }
248 
249  // one of the ancestors of the touched volume is supposed to be
250  // actually a G4PVPlacementInTPC that knows which TPC it covers;
251  // currently, it's the 4th in the ladder:
252  // [0] voxel [1] voxel tower [2] voxel plane [3] the full box;
253  G4int depth = 0;
254  while (depth < pTouchable->GetHistoryDepth()) {
255  const G4PVPlacementInTPC* pPVinTPC =
256  dynamic_cast<const G4PVPlacementInTPC*>(pTouchable->GetVolume(depth++));
257  if (!pPVinTPC) continue;
258  cryostat = pPVinTPC->ID.Cryostat;
259  tpc = pPVinTPC->ID.TPC;
260  if (Has(fSkipWireSignalInTPCs, tpc)) { return true; }
261  break;
262  } // while
263  if (depth < pTouchable->GetHistoryDepth()) {
264  // this is a fundamental error where the step does not happen in
265  // any TPC; this should not happen in the readout geometry!
266  throw cet::exception("LArG4") << "No TPC ID found in LArVoxelReadout::ProcessHits()";
267  } // if
268  MF_LOG_DEBUG("LArVoxelReadoutHit") << " hit in C=" << cryostat << " T=" << tpc;
269  } // if more than one TPC
270 
271  // Note that if there is no particle ID for this energy deposit, the
272  // trackID will be sim::NoParticleId.
273 
275  *fClockData, midPoint, g4time, trackID, cryostat, tpc, origTrackID);
276  } // end we are drifting
277  } // end there is non-zero energy deposition
278 
279  return true;
280  }
281 
282  //----------------------------------------------------------------------------
283  void LArVoxelReadout::SetRandomEngines(CLHEP::HepRandomEngine* pPropGen)
284  {
285  assert(pPropGen); // random engine must be present
286  fPropGen = pPropGen;
287  }
288 
289  //----------------------------------------------------------------------------
291  geo::PlaneGeo const& plane) const
292  {
293  //
294  // translate the landing position on the two frame coordinates
295  // ("width" and "depth")
296  //
297  auto const landingPos = plane.PointWidthDepthProjection(pos);
298 
299  //
300  // compute the distance of the landing position on the two frame
301  // coordinates ("width" and "depth");
302  // keep the point within 10 micrometers (0.001 cm) from the border
303  //
304  auto const offPlane = plane.DeltaFromActivePlane(landingPos, 0.001);
305 
306  //
307  // if both the distances are below the margin, move the point to
308  // the border
309  //
310 
311  // nothing to recover: landing is inside
312  if ((offPlane.X() == 0.0) && (offPlane.Y() == 0.0)) return pos;
313 
314  // won't recover: too far
315  if ((std::abs(offPlane.X()) > fOffPlaneMargin) || (std::abs(offPlane.Y()) > fOffPlaneMargin))
316  return pos;
317 
318  // we didn't fully decompose because it might be unnecessary;
319  // now we need the full thing
320  auto const distance = plane.DistanceFromPlane(pos);
321 
322  return plane.ComposePoint(distance, landingPos + offPlane);
323  } // LArVoxelReadout::RecoverOffPlaneDeposit()
324 
325  //----------------------------------------------------------------------------
326  // energy is passed in with units of MeV, dx has units of cm
328  G4ThreeVector stepMidPoint,
329  const double simTime,
330  int trackID,
331  unsigned short int cryostat,
332  unsigned short int tpc,
333  int origTrackID)
334  {
335  auto const tpcClock = clockData.TPCClock();
336 
337  // this must be always true, unless caller has been sloppy
338  assert(fPropGen); // No propagation random generator provided?!
339 
340  CLHEP::RandGauss PropRand(*fPropGen);
341 
342  // This routine gets called frequently, once per every particle
343  // traveling through every voxel. Use whatever tricks we can to
344  // increase its execution speed.
345 
346  static double LifetimeCorr_const = -1000. * fElectronLifetime;
347  static double LDiff_const = std::sqrt(2. * fLongitudinalDiffusion);
348  static double TDiff_const = std::sqrt(2. * fTransverseDiffusion);
349  static double RecipDriftVel[3] = {
350  1. / fDriftVelocity[0], 1. / fDriftVelocity[1], 1. / fDriftVelocity[2]};
351 
352  struct Deposit_t {
353  double energy = 0.;
354  double electrons = 0.;
355 
356  void add(double more_energy, double more_electrons)
357  {
358  energy += more_energy;
359  electrons += more_electrons;
360  }
361  }; // Deposit_t
362 
363  // Map of electrons to store - catalogued by map[channel][tdc]
364  std::map<raw::ChannelID_t, std::map<unsigned int, Deposit_t>> DepositsToStore;
365 
366  geo::Point_t xyz1;
367 
368  double const xyz[3] = {
369  stepMidPoint.x() / CLHEP::cm, stepMidPoint.y() / CLHEP::cm, stepMidPoint.z() / CLHEP::cm};
370 
371  // Already know which TPC we're in because we have been told
372 
373  try {
374  const geo::TPCGeo& tpcg = fGeoHandle->TPC(geo::TPCID{cryostat, tpc});
375 
376  // X drift distance - the drift direction can be either in
377  // the positive or negative direction, so use std::abs
378 
380  auto const& plane_center = tpcg.Plane(0).GetCenter();
381  double XDrift = std::abs(stepMidPoint.x() / CLHEP::cm - plane_center.X());
382  //std::cout<<tpcg.DriftDirection()<<std::endl;
383  if (tpcg.DriftDirection() == geo::kNegX)
384  XDrift = stepMidPoint.x() / CLHEP::cm - plane_center.X();
385  else if (tpcg.DriftDirection() == geo::kPosX)
386  XDrift = plane_center.X() - stepMidPoint.x() / CLHEP::cm;
387 
388  if (XDrift < 0.) return;
389 
390  // Get SCE {x,y,z} offsets for particular location in TPC
391  geo::Vector_t posOffsets;
392  auto const* SCE = lar::providerFrom<spacecharge::SpaceChargeService>();
393  if (SCE->EnableSimSpatialSCE() == true) {
394  posOffsets = SCE->GetPosOffsets({xyz[0], xyz[1], xyz[2]});
395  if (larsim::Utils::SCE::out_of_bounds(posOffsets)) { return; }
396  }
397  posOffsets.SetX(-posOffsets.X());
398 
399  // Drift time (nano-sec)
400  double TDrift;
401  XDrift += posOffsets.X();
402 
403  // Space charge distortion could push the energy deposit beyond the wire
404  // plane (see issue #15131). Given that we don't have any subtlety in the
405  // simulation of this region, bringing the deposit exactly on the plane
406  // should be enough for the time being.
407  if (XDrift < 0.) XDrift = 0.;
408 
409  TDrift = XDrift * RecipDriftVel[0];
410  if (tpcg.Nplanes() == 2) { // special case for ArgoNeuT (plane 0 is the second wire plane)
411  TDrift = ((XDrift - tpcg.PlanePitch(0, 1)) * RecipDriftVel[0] +
412  tpcg.PlanePitch(0, 1) * RecipDriftVel[1]);
413  }
414 
415  const double lifetimecorrection = TMath::Exp(TDrift / LifetimeCorr_const);
416  const int nIonizedElectrons =
419 
420  // if we have no electrons (too small energy or too large recombination)
421  // we are done already here
422  if (nIonizedElectrons <= 0) {
423  MF_LOG_DEBUG("LArVoxelReadout")
424  << "No electrons drifted to readout, " << energy << " MeV lost.";
425  return;
426  }
427  // includes the effect of lifetime
428  const double nElectrons = nIonizedElectrons * lifetimecorrection;
429 
430  // Longitudinal & transverse diffusion sigma (cm)
431  double SqrtT = std::sqrt(TDrift);
432  double LDiffSig = SqrtT * LDiff_const;
433  double TDiffSig = SqrtT * TDiff_const;
434  double electronclsize = fElectronClusterSize;
435 
436  int nClus = (int)std::ceil(nElectrons / electronclsize);
437  if (nClus < fMinNumberOfElCluster) {
438  electronclsize = nElectrons / fMinNumberOfElCluster;
439  if (electronclsize < 1.0) { electronclsize = 1.0; }
440  nClus = (int)std::ceil(nElectrons / electronclsize);
441  }
442 
443  // Compute arrays of values as quickly as possible.
444  std::vector<double> XDiff(nClus);
445  std::vector<double> YDiff(nClus);
446  std::vector<double> ZDiff(nClus);
447  std::vector<double> nElDiff(nClus, electronclsize);
448  std::vector<double> nEnDiff(nClus);
449 
450  // fix the number of electrons in the last cluster, that has smaller size
451  nElDiff.back() = nElectrons - (nClus - 1) * electronclsize;
452 
453  for (size_t xx = 0; xx < nElDiff.size(); ++xx) {
454  if (nElectrons > 0)
455  nEnDiff[xx] = energy / nElectrons * nElDiff[xx];
456  else
457  nEnDiff[xx] = 0.;
458  }
459 
460  double const averageYtransversePos = (stepMidPoint.y() / CLHEP::cm) + posOffsets.Y();
461  double const averageZtransversePos = (stepMidPoint.z() / CLHEP::cm) + posOffsets.Z();
462 
463  // Smear drift times by x position and drift time
464  if (LDiffSig > 0.0)
465  PropRand.fireArray(nClus, &XDiff[0], 0., LDiffSig);
466  else
467  XDiff.assign(nClus, 0.0);
468 
469  if (TDiffSig > 0.0) {
470  // Smear the Y,Z position by the transverse diffusion
471  PropRand.fireArray(nClus, &YDiff[0], averageYtransversePos, TDiffSig);
472  PropRand.fireArray(nClus, &ZDiff[0], averageZtransversePos, TDiffSig);
473  }
474  else {
475  YDiff.assign(nClus, averageYtransversePos);
476  ZDiff.assign(nClus, averageZtransversePos);
477  }
478 
479  // make a collection of electrons for each plane
480  for (size_t p = 0; p < tpcg.Nplanes(); ++p) {
481 
482  geo::PlaneGeo const& plane = tpcg.Plane(p);
483 
484  double Plane0Pitch = tpcg.Plane0Pitch(p);
485 
486  // "-" sign is because Plane0Pitch output is positive. Andrzej
487  xyz1.SetX(tpcg.Plane(0).GetCenter().X() - Plane0Pitch);
488 
489  // Drift nClus electron clusters to the induction plane
490  for (int k = 0; k < nClus; ++k) {
491  // Correct drift time for longitudinal diffusion and plane
492  double TDiff = TDrift + XDiff[k] * RecipDriftVel[0];
493  // Take into account different Efields between planes
494  // Also take into account special case for ArgoNeuT where Nplanes = 2.
495  for (size_t ip = 0; ip < p; ++ip) {
496  TDiff +=
497  tpcg.PlanePitch(ip, ip + 1) * RecipDriftVel[tpcg.Nplanes() == 3 ? ip + 1 : ip + 2];
498  }
499  xyz1.SetY(YDiff[k]);
500  xyz1.SetZ(ZDiff[k]);
501 
503 
504  // grab the nearest channel to the xyz1 position
505  try {
506  if (fOffPlaneMargin != 0) {
507  // get the effective position where to consider the charge landed;
508  //
509  // Some optimisations are possible; in particular, this method
510  // could be extended to inform us if the point was too far.
511  // Currently, if that is the case the code will proceed, find the
512  // point is off plane, emit a warning and skip the deposition.
513  //
514  xyz1 = RecoverOffPlaneDeposit(xyz1, plane);
515  } // if charge lands off plane
516  uint32_t channel = fGeoHandle->NearestChannel(xyz1, plane.ID());
517 
520  // Add potential decay/capture/etc delay effect, simTime.
521  unsigned int tdc = tpcClock.Ticks(clockData.G4ToElecTime(TDiff + simTime));
522 
523  // Add electrons produced by each cluster to the map
524  DepositsToStore[channel][tdc].add(nEnDiff[k], nElDiff[k]);
525  }
526  catch (cet::exception& e) {
527  MF_LOG_DEBUG("LArVoxelReadout")
528  << "unable to drift electrons from point (" << xyz[0] << "," << xyz[1] << ","
529  << xyz[2] << ") with exception " << e;
530  }
531  } // end loop over clusters
532  } // end loop over planes
533 
534  // Now store them in SimChannels
535  ChannelMap_t& ChannelDataMap = fChannelMaps[cryostat][tpc];
536 
537  // browse deposited data on each channel: (channel; deposit data in time)
538  for (auto const& deposit_per_channel : DepositsToStore) {
539 
540  raw::ChannelID_t channel = deposit_per_channel.first;
541 
542  // find whether we already have this channel
543  auto iChannelData = ChannelDataMap.find(channel);
544 
545  // channelData is the SimChannel these deposits are going to be added to
546  // If there is such a channel already, use it (first beanch).
547  // If it's a new channel, the inner assignment creates a new SimChannel
548  // in the map, and we save its reference in channelData
549  sim::SimChannel& channelData = (iChannelData == ChannelDataMap.end()) ?
550  (ChannelDataMap[channel] = sim::SimChannel(channel)) :
551  iChannelData->second;
552 
553  // go through all deposits, one for each TDC: (TDC, deposit data)
554  for (auto const& deposit_per_tdc : deposit_per_channel.second) {
555  channelData.AddIonizationElectrons(trackID,
556  deposit_per_tdc.first,
557  deposit_per_tdc.second.electrons,
558  xyz,
559  deposit_per_tdc.second.energy,
560  origTrackID);
561 
562  } // for deposit on TDCs
563  } // for deposit on channels
564 
565  } // end try intended to catch points where TPC can't be found
566  catch (cet::exception& e) {
567  MF_LOG_DEBUG("LArVoxelReadout") << "step cannot be found in a TPC\n" << e;
568  }
569 
570  return;
571  }
572 
573  //---------------------------------------------------------------------------------------
574  // Never used but still have to be defined for G4
577 
578 } // namespace larg4
CLHEP::HepRandomEngine * propGen
Random engine for charge propagation.
unsigned int fTPC
which TPC this LArVoxelReadout corresponds to
unsigned int NTPC(CryostatID const &cryoid=cryostat_zero) const
Returns the total number of TPCs in the specified cryostat.
Definition: GeometryCore.h:686
Point_t ComposePoint(WireDecomposedVector_t const &decomp) const
Returns the 3D point from composition of projection and distance.
Definition: PlaneGeo.h:924
double PlanePitch(unsigned int p1=0, unsigned int p2=1) const
Returns the center of the TPC volume in world coordinates [cm].
Definition: TPCGeo.cxx:361
Double_t xx
Definition: macro.C:12
Utilities related to art service access.
Collection of what it takes to set a LArVoxelReadout up.
Energy deposited on a readout channel by simulated tracks.
Definition: SimChannel.h:136
A G4PVPlacement with an additional identificator.
art::ServiceHandle< sim::LArG4Parameters const > fLgpHandle
Handle to the LArG4 parameters service.
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Vector_t
Type for representation of momenta in 3D space.
Definition: geo_vectors.h:160
Point_t const & GetCenter() const
Returns the centre of the wire plane in world coordinates [cm].
Definition: PlaneGeo.h:435
bool out_of_bounds(geo::Vector_t const &offset)
unsigned int Nplanes() const
Number of planes in this tpc.
Definition: TPCGeo.h:137
art::ServiceHandle< geo::Geometry const > fGeoHandle
Handle to the Geometry service.
Geometry information for a single TPC.
Definition: TPCGeo.h:36
double Temperature() const
In kelvin.
const ChannelMap_t & GetSimChannelMap() const
Returns the accumulated channel -> SimChannel map for the single TPC.
constexpr auto abs(T v)
Returns the absolute value of the argument.
LArVoxelReadout(std::string const &name)
Constructor. Can detect which TPC to cover by the name.
Geant4 interface.
for(Int_t i=0;i< nentries;i++)
Definition: comparison.C:30
void Setup(Setup_t const &setupData)
Reads all the configuration elements from setupData
WidthDepthProjection_t PointWidthDepthProjection(geo::Point_t const &point) const
Returns the projection of the specified point on the plane.
Definition: PlaneGeo.h:1000
void SetRandomEngines(CLHEP::HepRandomEngine *pPropGen)
Sets the random generators to be used.
detinfo::DetectorClocksData const * fClockData
Drift towards negative X values.
Definition: geo_types.h:167
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
Definition: GeometryCore.h:430
pure virtual base interface for detector clocks
double Efield(unsigned int planegap=0) const
kV/cm
TPCGeo const & TPC(TPCID const &tpcid=tpc_zero) const
Returns the specified TPC.
Definition: GeometryCore.h:722
bool bSingleTPC
true if this readout is associated with a single TPC
ID_t ID
Physical Volume identificator.
geo::Point_t RecoverOffPlaneDeposit(geo::Point_t const &pos, geo::PlaneGeo const &plane) const
Returns the point on the specified plane closest to position.
double TransverseDiffusion() const
void SetSingleTPC(unsigned int cryostat, unsigned int tpc)
Associates this readout to one specific TPC.
std::map< unsigned int, sim::SimChannel > ChannelMap_t
Type of map channel -> sim::SimChannel.
Use Geant4&#39;s user "hooks" to maintain a list of particles generated by Geant4.
bool NoElectronPropagation() const
std::vector< unsigned short int > fSkipWireSignalInTPCs
void DriftIonizationElectrons(detinfo::DetectorClocksData const &clockData, G4ThreeVector stepMidPoint, const double simTime, int trackID, unsigned short int cryostat, unsigned short int tpc, int origTrackID)
const std::vector< unsigned short int > & SkipWireSignalInTPCs() const
Utility function for testing if Space Charge offsets are out of bounds.
double ElectronClusterSize() const
Geometry information for a single wire plane.The plane is represented in the geometry by a solid whic...
Definition: PlaneGeo.h:78
double energy
Definition: plottest35.C:25
static IonizationAndScintillation * Instance()
ElecClock const & TPCClock() const noexcept
Borrow a const TPC clock with time set to Trigger time [us].
std::vector< sim::SimChannel > GetSimChannels() const
Creates a list with the accumulated information for the single TPC.
detinfo::DetectorPropertiesData const * fDetProp
double DriftVelocity(double efield=0., double temperature=0.) const
cm/us
bool Has(std::vector< unsigned short int > v, unsigned short int tpc) const
The data type to uniquely identify a TPC.
Definition: geo_types.h:381
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
Definition of data types for geometry description.
double Plane0Pitch(unsigned int p) const
Returns the center of the TPC volume in world coordinates [cm].
Definition: TPCGeo.cxx:314
DriftDirection_t DriftDirection() const
Returns an enumerator value describing the drift direction.
Definition: TPCGeo.h:121
std::vector< std::vector< ChannelMap_t > > fChannelMaps
Maps of cryostat, tpc to channel data.
virtual void EndOfEvent(G4HCofThisEvent *)
double fOffPlaneMargin
Charge deposited within this many [cm] from the plane is lead onto it.
void SetDiscoverTPC()
Sets this readout to discover the TPC of each processed hit.
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:180
A Geant4 sensitive detector that accumulates voxel information.
double offPlaneMargin
Margin for charge recovery (see LArVoxelReadout).
WidthDepthProjection_t DeltaFromActivePlane(WidthDepthProjection_t const &proj, double wMargin, double dMargin) const
Returns a projection vector that, added to the argument, gives a projection inside (or at the border ...
Definition: PlaneGeo.cxx:515
double G4ToElecTime(double const g4_time) const
Drift towards positive X values.
Definition: geo_types.h:166
Contains all timing reference information for the detector.
int MinNumberOfElCluster() const
double DistanceFromPlane(geo::Point_t const &point) const
Returns the distance of the specified point from the wire plane.
Definition: PlaneGeo.h:558
unsigned int fCstat
and in which cryostat (if bSingleTPC is true)
void SetOffPlaneChargeRecoveryMargin(double margin)
Sets the margin for recovery of charge drifted off-plane.
Singleton to access a unified treatment of ionization and scintillation in LAr.
geo::PlaneID const & ID() const
Returns the identifier of this plane.
Definition: PlaneGeo.h:184
CLHEP::HepRandomEngine * fPropGen
random engine for charge propagation
#define MF_LOG_DEBUG(id)
virtual void Initialize(G4HCofThisEvent *)
raw::ChannelID_t NearestChannel(Point_t const &worldLoc, PlaneID const &planeid) const
Returns the ID of the channel nearest to the specified position.
Transports energy depositions from GEANT4 to TPC channels.
PlaneGeo const & Plane(geo::View_t view) const
Return the plane in the tpc with View_t view.
Definition: TPCGeo.cxx:252
void AddIonizationElectrons(TrackID_t trackID, TDC_t tdc, double numberElectrons, double const *xyz, double energy, TrackID_t origTrackID=util::kBogusI)
Add ionization electrons and energy to this channel.
Definition: SimChannel.cxx:49
double LongitudinalDiffusion() const
unsigned int ChannelID_t
Type representing the ID of a readout channel.
Definition: RawTypes.h:28
Float_t e
Definition: plot.C:35
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
bool DisableWireplanes() const
The data type to uniquely identify a cryostat.
Definition: geo_types.h:192