LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
LArVoxelReadoutGeometry.cxx
Go to the documentation of this file.
1 
9 
10 // C/C++ libraries
11 #include <cmath>
12 #include <map>
13 #include <memory> // std::unique_ptr()
14 #include <sstream> // std::ostringstream
15 #include <typeinfo>
16 
17 // Framework includes
19 #include "cetlib_except/demangle.h"
21 
22 // LArSoft includes
27 
28 // G4 includes
29 #include "Geant4/G4Box.hh"
30 #include "Geant4/G4LogicalVolume.hh"
31 #include "Geant4/G4PVReplica.hh"
32 #include "Geant4/G4Point3D.hh"
33 #include "Geant4/G4RotationMatrix.hh"
34 #include "Geant4/G4SDManager.hh"
35 #include "Geant4/G4ThreeVector.hh"
36 #include "Geant4/G4Tubs.hh"
37 #include "Geant4/G4VSolid.hh"
38 #include "Geant4/G4VisAttributes.hh"
39 
40 constexpr bool DisableVoxelCaching = false;
41 
42 namespace {
43  template <typename T>
44  std::string demangle_cxx_symbol(const T& obj)
45  {
46  return cet::demangle_symbol(typeid(obj).name());
47  }
48 
49  template <class STREAM>
50  void DumpPhysicalVolume(STREAM& out, const G4VPhysicalVolume& PV, std::string indentstr = "")
51  {
52  const G4ThreeVector& pos = PV.GetTranslation();
53  const G4LogicalVolume* LV = PV.GetLogicalVolume();
54 
55  out << indentstr << PV.GetName() << " [" << demangle_cxx_symbol(PV) << "]"
56  << " at (" << pos.x() << ", " << pos.y() << ", " << pos.z() << ")";
57 
58  if (!LV) {
59  out << " with no logical volume (!?)\n";
60  return;
61  }
62 
63  out << ", a " << LV->GetName();
64 
65  const G4VSolid* Solid = LV->GetSolid();
66  if (Solid) {
67  out << " shaped as " << Solid->GetName();
68  if (auto pBox = dynamic_cast<const G4Box*>(Solid)) {
69  out << ", a (" << (2. * pBox->GetXHalfLength()) << " x " << (2. * pBox->GetYHalfLength())
70  << " x " << (2. * pBox->GetZHalfLength()) << ") cm box";
71  }
72  else {
73  out << ", a " << demangle_cxx_symbol(*Solid);
74  }
75  }
76  else {
77  out << " with no shape (!?!)";
78  }
79 
80  G4int nDaughters = LV->GetNoDaughters();
81  if (nDaughters > 0) {
82  out << " with " << nDaughters << " subvolumes:\n";
83  for (G4int i = 0; i < nDaughters; ++i) {
84  DumpPhysicalVolume(out, *LV->GetDaughter(i), indentstr + " ");
85  }
86  }
87  else
88  out << '\n';
89  } // DumpPhysicalVolume()
90 
91 }
92 
93 namespace larg4 {
94 
95  // Constructor and destructor.
96  LArVoxelReadoutGeometry::LArVoxelReadoutGeometry(const G4String& name, Setup_t const& setupData)
97  : G4VUserParallelWorld(name), fReadoutSetupData(setupData.readoutSetup)
98  {
100  fStepLimit = std::make_unique<G4UserLimits>(ios->StepSizeLimit());
101  }
102 
105  {
106  // With a "parallel geometry", Geant4 has already created a clone
107  // of the world physical and logical volumes. We want to place
108  // the LAr TPC, and only the LAr TPC, within this cloned world.
109 
110  // Get the parallel world physical volume.
111  G4VPhysicalVolume* parallelPhysical = GetWorld();
112 
113  // Now we want to place a parallel LAr TPC volume within this
114  // parallel world volume. We only want to duplicate the LAr TPC
115  // volume; any other volumes in the "official" geometry are going
116  // to be ignored. Our parallel world will consist only of the
117  // world volume and the LAr TPC volume.
118 
119  G4Transform3D worldTransform;
120 
121  // first get the volDetEnclosure
122  std::string daughterName("volDetEnclosure");
123  G4Transform3D detEnclosureTransform;
124  G4VPhysicalVolume* detEnclosureVolume = FindNestedVolume(g4b::DetectorConstruction::GetWorld(),
125  worldTransform,
126  detEnclosureTransform,
127  daughterName,
128  0);
129 
130  // we keep track of all the voxel boxes with different geometries we create
131  struct VoxelSpecs_t {
132  double w, h, d;
133  unsigned int nw, nh, nd;
134 
136  bool operator<(const VoxelSpecs_t& vs) const
137  {
138  if (w < vs.w) return true;
139  if (w > vs.w) return false;
140  if (h < vs.h) return true;
141  if (h > vs.h) return false;
142  if (d < vs.d) return true;
143  if (d > vs.d) return false;
144  if (nw < vs.nw) return true;
145  if (nw > vs.nw) return false;
146  if (nh < vs.nh) return true;
147  if (nh > vs.nh) return false;
148  if (nd < vs.nd) return true;
149  return false;
150  } // operator<
151  }; // VoxelSpecs_t
152 
153  // TODO we don't need to cache the cell any more: a simple map will suffice
154  struct VoxelVolumes_t {
155  G4LogicalVolume* pBox;
156  G4LogicalVolume* pCell;
157 
158  VoxelVolumes_t(G4LogicalVolume* box = nullptr, G4LogicalVolume* cell = nullptr)
159  : pBox(box), pCell(cell)
160  {}
161  }; // VoxelVolumes_t
162 
163  // Define the sensitive detector for the voxel readout. This class
164  // routines will be called every time a particle deposits energy in
165  // a voxel that overlaps the LAr TPC.
166  flarVoxelReadout = new LArVoxelReadout("LArVoxelSD");
168  if ((fGeo->Ncryostats() == 1) && (fGeo->Cryostat().NTPC() == 1))
169  flarVoxelReadout->SetSingleTPC(0, 0); // just one TPC in the detector...
170 
171  // Tell Geant4's sensitive-detector manager that the voxel SD
172  // class exists.
173  G4SDManager* sdManager = G4SDManager::GetSDMpointer();
174  sdManager->AddNewDetector(flarVoxelReadout);
175 
176  // hope hashing doubles is reliable...
177  using VoxelCache_t = std::map<VoxelSpecs_t, VoxelVolumes_t>;
178  VoxelCache_t VoxelCache;
179 
180  for (unsigned int c = 0; c < fGeo->Ncryostats(); ++c) {
181 
182  // next get the cryostat
183  daughterName = "volCryostat";
184  G4Transform3D cryostatTransform;
185  G4VPhysicalVolume* cryostatVolume = FindNestedVolume(
186  detEnclosureVolume, detEnclosureTransform, cryostatTransform, daughterName, c);
187 
188  for (unsigned int t = 0; t < fGeo->Cryostat(geo::CryostatID(c)).NTPC(); ++t) {
189 
190  // now for the TPC
191  daughterName = "volTPC";
192  G4Transform3D tpcTransform;
193  G4VPhysicalVolume* tpcVolume =
194  FindNestedVolume(cryostatVolume, cryostatTransform, tpcTransform, daughterName, t);
195 
196  daughterName = "volTPCActive";
197  G4Transform3D transform;
198  G4VPhysicalVolume* larTPCPhysical =
199  FindNestedVolume(tpcVolume, tpcTransform, transform, daughterName, t);
200 
201  // Get the LAr TPC volume, and its shape.
202  G4LogicalVolume* larTPCLogical = larTPCPhysical->GetLogicalVolume();
203  larTPCLogical->SetUserLimits(fStepLimit.get());
204 
205  G4VSolid* larTPCShape = larTPCLogical->GetSolid();
206 
207  // We're not going to exactly duplicate the LAr TPC in our
208  // parallel world. We're going to construct a box of voxels.
209  // What should the size and position of that box be?
210 
211  // To get our first hints, we need the overall dimensions of a
212  // "bounding box" that contains the shape. For now, we'll allow
213  // two possible shapes: a box (by the far the most likely) and a
214  // cylinder (for bizarre future detectors that I know nothing
215  // about).
216 
217  G4double larTPCHalfXLength = 0;
218  G4double larTPCHalfYLength = 0;
219  G4double larTPCHalfZLength = 0;
220  G4Box* tpcBox = dynamic_cast<G4Box*>(larTPCShape);
221  if (tpcBox != 0) {
222  larTPCHalfXLength = tpcBox->GetXHalfLength();
223  larTPCHalfYLength = tpcBox->GetYHalfLength();
224  larTPCHalfZLength = tpcBox->GetZHalfLength();
225  }
226  else {
227  // It's not a box. Try a cylinder.
228  G4Tubs* tube = dynamic_cast<G4Tubs*>(larTPCShape);
229  if (tube != 0) {
230  larTPCHalfXLength = tube->GetOuterRadius();
231  larTPCHalfYLength = tube->GetOuterRadius();
232  larTPCHalfZLength = tube->GetZHalfLength();
233  }
234  else {
235  throw cet::exception("LArVoxelReadoutGeometry")
236  << "Unknown shape in readout geometry"
237  << "The LAr TPC volume is not a box or a tube. "
238  << "This routine can't convert any other shapes.\n";
239  }
240  }
241 
242  MF_LOG_DEBUG("LArVoxelReadoutGeometry") << ": larTPCHalfXLength=" << larTPCHalfXLength
243  << ": larTPCHalfYLength=" << larTPCHalfYLength
244  << ": larTPCHalfZLength=" << larTPCHalfZLength;
245 
246  // Get some constants from the LAr voxel information object.
247  // Remember, ROOT uses cm.
249  G4double voxelSizeX = lvc->VoxelSizeX() * CLHEP::cm;
250  G4double voxelSizeY = lvc->VoxelSizeY() * CLHEP::cm;
251  G4double voxelSizeZ = lvc->VoxelSizeZ() * CLHEP::cm;
252  G4double voxelOffsetX = lvc->VoxelOffsetX() * CLHEP::cm;
253  G4double voxelOffsetY = lvc->VoxelOffsetY() * CLHEP::cm;
254  G4double voxelOffsetZ = lvc->VoxelOffsetZ() * CLHEP::cm;
255 
256  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
257  << ": voxelSizeX=" << voxelSizeX << ", voxelSizeY=" << voxelSizeY
258  << ", voxelSizeZ=" << voxelSizeZ;
259 
260  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
261  << ": voxelOffsetX=" << voxelOffsetX << ", voxelOffsetY=" << voxelOffsetY
262  << ", voxelOffsetZ=" << voxelOffsetZ;
263 
264  // We want our voxelization region to be an integer multiple of
265  // the voxel sizes in all directions; if we didn't do this, we
266  // might get into trouble when we start playing with replicas.
267  // Compute the the dimensions of our voxelization to be about the
268  // size of the LAr TPC region, adjusted to be an integer number of
269  // voxels in all directions.
270 
271  G4double numberXvoxels = 2. * larTPCHalfXLength / voxelSizeX;
272  G4double numberYvoxels = 2. * larTPCHalfYLength / voxelSizeY;
273  G4double numberZvoxels = 2. * larTPCHalfZLength / voxelSizeZ;
274  numberXvoxels = trunc(numberXvoxels) + 1.;
275  numberYvoxels = trunc(numberYvoxels) + 1.;
276  numberZvoxels = trunc(numberZvoxels) + 1.;
277 
278  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
279  << "Active volume of cryo #" << c << " TPC #" << t << " will be split in "
280  << numberXvoxels << " x " << numberYvoxels << " x " << numberZvoxels << " = "
281  << (numberXvoxels * numberYvoxels * numberZvoxels) << " voxels of size " << voxelSizeX
282  << " x " << voxelSizeY << " x " << voxelSizeZ << " cm";
283 
284  VoxelSpecs_t VoxelSpecs{/* w */ 2. * larTPCHalfXLength,
285  /* h */ 2. * larTPCHalfYLength,
286  /* d */ 2. * larTPCHalfZLength,
287  /* nw */ (unsigned int)numberXvoxels,
288  /* nh */ (unsigned int)numberYvoxels,
289  /* nd */ (unsigned int)numberZvoxels};
290 
291  VoxelCache_t::iterator iVoxelVol =
292  DisableVoxelCaching ? VoxelCache.end() : VoxelCache.find(VoxelSpecs);
293  if (iVoxelVol == VoxelCache.end()) {
294  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
295  << "Creating a new voxel volume " << VoxelSpecs.w << " x " << VoxelSpecs.h << " x "
296  << VoxelSpecs.d << " cm, " << VoxelSpecs.nw << " x " << VoxelSpecs.nh << " x "
297  << VoxelSpecs.nd << " voxels";
298 
299  G4double voxelBoxHalfX = numberXvoxels * voxelSizeX / 2.;
300  G4double voxelBoxHalfY = numberYvoxels * voxelSizeY / 2.;
301  G4double voxelBoxHalfZ = numberZvoxels * voxelSizeZ / 2.;
302 
303  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
304  << ": voxelBoxHalfX=" << voxelBoxHalfX << ", voxelBoxHalfY=" << voxelBoxHalfY
305  << ", voxelBoxHalfZ=" << voxelBoxHalfZ;
306 
307  // Now we have a box that will include an integer number of voxels
308  // in each direction. Note that the material is irrelevant for a
309  // "parallel world."
310  auto voxelBox = new G4Box("VoxelBox", voxelBoxHalfX, voxelBoxHalfY, voxelBoxHalfZ);
311  auto voxelBoxLogical = new G4LogicalVolume(voxelBox, 0, "VoxelizationLogicalVolume");
312 
313  // If we generate an event display within Geant4, we won't want to
314  // see this box.
315  auto invisible = new G4VisAttributes();
316  invisible->SetVisibility(false);
317  voxelBoxLogical->SetVisAttributes(invisible);
318 
319  // Now we've fill our "box of voxels" with the voxels themselves.
320  // We'll do this by sub-dividing the volume in x, then y, then z.
321 
322  // Create an "x-slice".
323  G4Box* xSlice = new G4Box("xSlice", voxelSizeX / 2., voxelBoxHalfY, voxelBoxHalfZ);
324  G4LogicalVolume* xSliceLogical = new G4LogicalVolume(xSlice, 0, "xLArVoxelSlice");
325  xSliceLogical->SetVisAttributes(invisible);
326 
327  // Use replication to slice up the "box of voxels" along the x-axis.
328  new G4PVReplica("VoxelSlicesInX",
329  xSliceLogical,
330  voxelBoxLogical,
331  kXAxis,
332  G4int(numberXvoxels),
333  voxelSizeX);
334 
335  // Now do the same thing, dividing that x-slice along the y-axis.
336  auto ySlice = new G4Box("ySlice", voxelSizeX / 2., voxelSizeY / 2., voxelBoxHalfZ);
337  auto ySliceLogical = new G4LogicalVolume(ySlice, 0, "yLArVoxelSlice");
338  ySliceLogical->SetVisAttributes(invisible);
339  new G4PVReplica("VoxelSlicesInY",
340  ySliceLogical,
341  xSliceLogical,
342  kYAxis,
343  G4int(numberYvoxels),
344  voxelSizeY);
345 
346  // Now divide the y-slice along the z-axis, giving us our actual
347  // voxels.
348  auto zSlice = new G4Box("zSlice", voxelSizeX / 2., voxelSizeY / 2., voxelSizeZ / 2.);
349  auto voxelLogical = new G4LogicalVolume(zSlice, 0, "LArVoxel");
350  voxelLogical->SetVisAttributes(invisible);
351  new G4PVReplica(
352  "LArVoxel", voxelLogical, ySliceLogical, kZAxis, G4int(numberZvoxels), voxelSizeZ);
353 
354  iVoxelVol = VoxelCache.emplace(VoxelSpecs, VoxelVolumes_t(voxelBoxLogical, voxelLogical))
355  .first; // iterator to the inserted element
356 
357  // Set the sensitive detector of this LAr TPC (logical) volume
358  // to be the voxel readout.
359  voxelLogical->SetSensitiveDetector(flarVoxelReadout);
360 
361  } // if not cached yet
362  G4LogicalVolume* voxelBoxLogical = iVoxelVol->second.pBox;
363 
364  // We have a "box of voxels" that's the right dimensions, but we
365  // have to know exactly where to put it. The user has the option
366  // to offset the voxel co-ordinate system. We want to place our
367  // box so the edges of our voxels align with that co-ordinate
368  // system. In effect, we want to offset our "box of voxels" by
369  // the user's offsets, modulo the size of the voxel in each
370  // direction.
371 
372  G4double offsetInVoxelsX = voxelOffsetX / voxelSizeX;
373  G4double offsetInVoxelsY = voxelOffsetY / voxelSizeY;
374  G4double offsetInVoxelsZ = voxelOffsetZ / voxelSizeZ;
375  G4double fractionOffsetX = offsetInVoxelsX - trunc(offsetInVoxelsX);
376  G4double fractionOffsetY = offsetInVoxelsY - trunc(offsetInVoxelsY);
377  G4double fractionOffsetZ = offsetInVoxelsZ - trunc(offsetInVoxelsZ);
378  G4double offsetX = fractionOffsetX * voxelSizeX;
379  G4double offsetY = fractionOffsetY * voxelSizeY;
380  G4double offsetZ = fractionOffsetZ * voxelSizeZ;
381 
382  // Now we know how much to offset the "box of voxels". Include
383  // that in the transformation of the co-ordinates from world
384  // volume to LAr TPC volume.
385  transform = G4Translate3D(offsetX, offsetY, offsetZ) * transform;
386 
387  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
388  << ": offsetX=" << offsetX << ", offsetY=" << offsetY << ", offsetZ=" << offsetZ;
389 
390  // suffix representing this TPC
391  std::ostringstream nums;
392  nums << "_Cryostat" << c << "_TPC" << t;
393 
394  // Place the box of voxels, with the accumulated transformations
395  // computed above.
396  new G4PVPlacementInTPC(transform,
397  "VoxelizationPhysicalVolume" + nums.str(),
398  voxelBoxLogical,
399  parallelPhysical,
400  false, // Only one volume
401  0, // Copy number
402  false, // No surface check
403  {(unsigned short int)c, (unsigned short int)t} // TPC ID
404  );
405 
406  } // end loop over tpcs
407  } // end loop over cryostats
408 
409  if (mf::isDebugEnabled()) {
410  MF_LOG_DEBUG("LArVoxelReadoutGeometryDump") << "Dump of voxelized volume";
411  {
412  mf::LogDebug log("LArVoxelReadoutGeometryDump");
413  DumpPhysicalVolume(log, *parallelPhysical);
414  }
415  MF_LOG_DEBUG("LArVoxelReadoutGeometryDump") << "End of dump of voxelized volume";
416  }
417  }
418 
419  //---------------------------------------------------------------
421  // We know the ordering of the volumes in the Geometry,
422  // see https://cdcvs.fnal.gov/redmine/projects/larsoftsvn/wiki/Geometry
423  // Make use of that knowledge to efficiently get the desired volumes and
424  // their total transforms.
425  // the daughterTransform is the total transform to the world coordinate system
426  G4VPhysicalVolume* LArVoxelReadoutGeometry::FindNestedVolume(G4VPhysicalVolume* mother,
427  G4Transform3D& motherTransform,
428  G4Transform3D& daughterTransform,
429  std::string const& daughterName,
430  unsigned int expectedNum)
431  {
432  G4LogicalVolume* logicalVolume = mother->GetLogicalVolume();
433  G4int numberDaughters = logicalVolume->GetNoDaughters();
434  for (G4int i = 0; i != numberDaughters; ++i) {
435  G4VPhysicalVolume* d = logicalVolume->GetDaughter(i);
436 
437  if (d->GetName().contains(daughterName)) {
438 
439  // check that this cryostat is the requested one using fCryostat
440  G4ThreeVector translation = d->GetObjectTranslation();
441  G4RotationMatrix rotation = d->GetObjectRotationValue();
442  G4Transform3D transform(rotation, translation);
443  daughterTransform = motherTransform * transform;
444 
445  // take the origin of the volume and transform it to
446  // world coordinated
447  G4Point3D local(0., 0., 0.);
448  G4Point3D world = daughterTransform * local;
449 
450  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
451  << "current " << daughterName << " origin is at (" << world.x() / CLHEP::cm << ","
452  << world.y() / CLHEP::cm << "," << world.z() / CLHEP::cm << ")";
453 
454  // we don't bother with the cryostat number when calling Geometry::PositionToTPC
455  // because we know we have already started off with the correct cryostat volume
456  // G4 uses mm, we want cm
457  geo::Point_t const worldPos{
458  world.x() / CLHEP::cm, world.y() / CLHEP::cm, world.z() / CLHEP::cm};
459  unsigned int daughterNum = 0;
460  if (daughterName.compare("volCryostat") == 0)
461  daughterNum = fGeo->PositionToCryostatID(worldPos).Cryostat;
462  else if (daughterName.compare("volTPC") == 0)
463  daughterNum = fGeo->PositionToTPCID(worldPos).TPC;
464  else if (daughterName.compare("volTPCActive") == 0 ||
465  daughterName.compare("volDetEnclosure") == 0) {
466  // for either of these volumes, we know there is only 1 in the mother volume
467  MF_LOG_DEBUG("LArVoxelReadoutGeometry") << "found the desired " << daughterName;
468  return d;
469  }
470 
471  // if we found the desired volume, stop looking
472  if (daughterNum == expectedNum) {
473  MF_LOG_DEBUG("LArVoxelReadoutGeometry") << "found the desired " << daughterName;
474  return d;
475  }
476  } // end if the volume has the right name
477  } // end loop over volumes
478 
479  throw cet::exception("LArVoxelReadoutGeometry")
480  << "could not find the desired " << daughterName << " to make LArVoxelReadoutGeometry\n";
481  }
482 
483 } // namespace larg4
Build Geant4 geometry from GDML.
intermediate_table::iterator iterator
Collection of all it takes to set up this object.
double VoxelSizeX() const
Access to voxel dimensions and offsets.
larg4::LArVoxelReadout::Setup_t fReadoutSetupData
Encapsulates calculation of LArVoxelID and LArVoxel parameters.
G4VPhysicalVolume * FindNestedVolume(G4VPhysicalVolume *mother, G4Transform3D &motherTransform, G4Transform3D &daughterTransform, std::string const &daughterName, unsigned int expectedNum)
constexpr bool operator<(CryostatID const &a, CryostatID const &b)
Order cryostats with increasing ID.
Definition: geo_types.h:692
static G4VPhysicalVolume * GetWorld()
Geant4 interface.
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:211
void Setup(Setup_t const &setupData)
Reads all the configuration elements from setupData
larg4::LArVoxelReadout * flarVoxelReadout
Data for LArVoxelReadout setup.
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
Definition: GeometryCore.h:430
CryostatGeo const & Cryostat(CryostatID const &cryoid=cryostat_zero) const
Returns the specified cryostat.
Define the "parallel" geometry that&#39;s seen by the LAr Voxels.
Float_t voxelSizeY
Definition: plot.C:39
void SetSingleTPC(unsigned int cryostat, unsigned int tpc)
Associates this readout to one specific TPC.
std::unique_ptr< G4UserLimits > fStepLimit
constexpr bool DisableVoxelCaching
Float_t voxelSizeZ
Definition: plot.C:39
static IonizationAndScintillation * Instance()
Float_t d
Definition: plot.C:235
unsigned int NTPC() const
Number of TPCs in this cryostat.
Definition: CryostatGeo.h:175
G4PVPlacementWithID< TPCID_t > G4PVPlacementInTPC
A physical volume with a TPC ID.
bool isDebugEnabled()
double VoxelOffsetX() const
art::ServiceHandle< geo::Geometry const > fGeo
Handle to the geometry.
Float_t voxelSizeX
Definition: plot.C:39
double VoxelOffsetZ() const
CryostatID PositionToCryostatID(Point_t const &point) const
Returns the ID of the cryostat at specified location.
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.
Singleton to access a unified treatment of ionization and scintillation in LAr.
#define MF_LOG_DEBUG(id)
Transports energy depositions from GEANT4 to TPC channels.
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:399
Float_t w
Definition: plot.C:20
TPCID PositionToTPCID(Point_t const &point) const
Returns the ID of the TPC at specified location.
double VoxelOffsetY() const
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
The data type to uniquely identify a cryostat.
Definition: geo_types.h:192