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