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

#include <LArVoxelReadoutGeometry.h>

Inheritance diagram for larg4::LArVoxelReadoutGeometry:

Classes

class  Sentry
 
struct  Setup_t
 Collection of all it takes to set up this object. More...
 

Public Member Functions

 LArVoxelReadoutGeometry (const G4String name, Setup_t const &setupData)
 Constructor: sets up all its LArVoxelReadout instances. More...
 
void Construct () override
 

Private Member Functions

G4VPhysicalVolume * FindNestedVolume (G4VPhysicalVolume *mother, G4Transform3D &motherTransform, G4Transform3D &daughterTransform, std::string &daughterName, unsigned int expectedNum)
 

Private Attributes

art::ServiceHandle
< geo::Geometry const > 
fGeo
 Handle to the geometry. More...
 
std::unique_ptr< G4UserLimits > fStepLimit
 
larg4::LArVoxelReadoutflarVoxelReadout {nullptr}
 Data for LArVoxelReadout setup. More...
 
larg4::LArVoxelReadout::Setup_t fReadoutSetupData
 

Detailed Description

Definition at line 41 of file LArVoxelReadoutGeometry.h.

Constructor & Destructor Documentation

larg4::LArVoxelReadoutGeometry::LArVoxelReadoutGeometry ( const G4String  name,
Setup_t const &  setupData 
)

Constructor: sets up all its LArVoxelReadout instances.

Definition at line 108 of file LArVoxelReadoutGeometry.cxx.

109  : G4VUserParallelWorld(name), fReadoutSetupData(setupData.readoutSetup)
110  {
112  fStepLimit = std::make_unique<G4UserLimits>(ios->StepSizeLimit());
113  }
larg4::LArVoxelReadout::Setup_t fReadoutSetupData
std::unique_ptr< G4UserLimits > fStepLimit
static IonizationAndScintillation * Instance()
then echo fcl name

Member Function Documentation

void larg4::LArVoxelReadoutGeometry::Construct ( )
override

The key method in this class; creates a parallel world view of those volumes relevant to the LAr voxel readout. Required of any class that inherits from G4VUserParallelWorld

< sizes of the volume

< divisions in each volume

have some sorting...

Definition at line 117 of file LArVoxelReadoutGeometry.cxx.

118  {
119  // With a "parallel geometry", Geant4 has already created a clone
120  // of the world physical and logical volumes. We want to place
121  // the LAr TPC, and only the LAr TPC, within this cloned world.
122 
123  // Get the parallel world physical volume.
124  G4VPhysicalVolume* parallelPhysical = GetWorld();
125 
126  // Now we want to place a parallel LAr TPC volume within this
127  // parallel world volume. We only want to duplicate the LAr TPC
128  // volume; any other volumes in the "official" geometry are going
129  // to be ignored. Our parallel world will consist only of the
130  // world volume and the LAr TPC volume.
131 
132  G4Transform3D worldTransform;
133 
134  // first get the volDetEnclosure
135  std::string daughterName("volDetEnclosure");
136  G4Transform3D detEnclosureTransform;
137  G4VPhysicalVolume* detEnclosureVolume =
138  this->FindNestedVolume(g4b::DetectorConstruction::GetWorld(),
139  worldTransform,
140  detEnclosureTransform,
141  daughterName,
142  0);
143 
144  // we keep track of all the voxel boxes with different geometries we create
145  struct VoxelSpecs_t {
146  double w, h, d; ///< sizes of the volume
147  unsigned int nw, nh, nd; ///< divisions in each volume
148 
149  /// have some sorting...
150  bool
151  operator<(const VoxelSpecs_t& vs) const
152  {
153  if (w < vs.w) return true;
154  if (w > vs.w) return false;
155  if (h < vs.h) return true;
156  if (h > vs.h) return false;
157  if (d < vs.d) return true;
158  if (d > vs.d) return false;
159  if (nw < vs.nw) return true;
160  if (nw > vs.nw) return false;
161  if (nh < vs.nh) return true;
162  if (nh > vs.nh) return false;
163  if (nd < vs.nd) return true;
164  return false;
165  } // operator<
166  }; // VoxelSpecs_t
167 
168  // TODO we don't need to cache the cell any more: a simple map will suffice
169  struct VoxelVolumes_t {
170  G4LogicalVolume* pBox;
171  G4LogicalVolume* pCell;
172 
173  VoxelVolumes_t(G4LogicalVolume* box = nullptr, G4LogicalVolume* cell = nullptr)
174  : pBox(box), pCell(cell)
175  {}
176  }; // VoxelVolumes_t
177 
178  // Define the sensitive detector for the voxel readout. This class
179  // routines will be called every time a particle deposits energy in
180  // a voxel that overlaps the LAr TPC.
181  flarVoxelReadout = new LArVoxelReadout("LArVoxelSD");
183  if ((fGeo->Ncryostats() == 1) && (fGeo->Cryostat(0).NTPC() == 1))
184  flarVoxelReadout->SetSingleTPC(0, 0); // just one TPC in the detector...
185 
186  // Tell Geant4's sensitive-detector manager that the voxel SD
187  // class exists.
188  G4SDManager* sdManager = G4SDManager::GetSDMpointer();
189  sdManager->AddNewDetector(flarVoxelReadout);
190 
191  // hope hashing doubles is reliable...
192  using VoxelCache_t = std::map<VoxelSpecs_t, VoxelVolumes_t>;
193  VoxelCache_t VoxelCache;
194 
195  for (unsigned int c = 0; c < fGeo->Ncryostats(); ++c) {
196 
197  // next get the cryostat
198  daughterName = "volCryostat";
199  G4Transform3D cryostatTransform;
200  G4VPhysicalVolume* cryostatVolume = this->FindNestedVolume(
201  detEnclosureVolume, detEnclosureTransform, cryostatTransform, daughterName, c);
202 
203  for (unsigned int t = 0; t < fGeo->Cryostat(c).NTPC(); ++t) {
204 
205  // now for the TPC
206  daughterName = "volTPC";
207  G4Transform3D tpcTransform;
208  G4VPhysicalVolume* tpcVolume =
209  this->FindNestedVolume(cryostatVolume, cryostatTransform, tpcTransform, daughterName, t);
210 
211  daughterName = "volTPCActive";
212  G4Transform3D transform;
213  G4VPhysicalVolume* larTPCPhysical =
214  this->FindNestedVolume(tpcVolume, tpcTransform, transform, daughterName, t);
215 
216  // Get the LAr TPC volume, and its shape.
217  G4LogicalVolume* larTPCLogical = larTPCPhysical->GetLogicalVolume();
218  larTPCLogical->SetUserLimits(fStepLimit.get());
219 
220  G4VSolid* larTPCShape = larTPCLogical->GetSolid();
221 
222  // We're not going to exactly duplicate the LAr TPC in our
223  // parallel world. We're going to construct a box of voxels.
224  // What should the size and position of that box be?
225 
226  // To get our first hints, we need the overall dimensions of a
227  // "bounding box" that contains the shape. For now, we'll allow
228  // two possible shapes: a box (by the far the most likely) and a
229  // cylinder (for bizarre future detectors that I know nothing
230  // about).
231 
232  G4double larTPCHalfXLength = 0;
233  G4double larTPCHalfYLength = 0;
234  G4double larTPCHalfZLength = 0;
235  G4Box* tpcBox = dynamic_cast<G4Box*>(larTPCShape);
236  if (tpcBox != 0) {
237  larTPCHalfXLength = tpcBox->GetXHalfLength();
238  larTPCHalfYLength = tpcBox->GetYHalfLength();
239  larTPCHalfZLength = tpcBox->GetZHalfLength();
240  }
241  else {
242  // It's not a box. Try a cylinder.
243  G4Tubs* tube = dynamic_cast<G4Tubs*>(larTPCShape);
244  if (tube != 0) {
245  larTPCHalfXLength = tube->GetOuterRadius();
246  larTPCHalfYLength = tube->GetOuterRadius();
247  larTPCHalfZLength = tube->GetZHalfLength();
248  }
249  else {
250  throw cet::exception("LArVoxelReadoutGeometry")
251  << "Unknown shape in readout geometry"
252  << "The LAr TPC volume is not a box or a tube. "
253  << "This routine can't convert any other shapes.\n";
254  }
255  }
256 
257  MF_LOG_DEBUG("LArVoxelReadoutGeometry") << ": larTPCHalfXLength=" << larTPCHalfXLength
258  << ": larTPCHalfYLength=" << larTPCHalfYLength
259  << ": larTPCHalfZLength=" << larTPCHalfZLength;
260 
261  // Get some constants from the LAr voxel information object.
262  // Remember, ROOT uses cm.
263  art::ServiceHandle<sim::LArVoxelCalculator const> lvc;
264  G4double voxelSizeX = lvc->VoxelSizeX() * CLHEP::cm;
265  G4double voxelSizeY = lvc->VoxelSizeY() * CLHEP::cm;
266  G4double voxelSizeZ = lvc->VoxelSizeZ() * CLHEP::cm;
267  G4double voxelOffsetX = lvc->VoxelOffsetX() * CLHEP::cm;
268  G4double voxelOffsetY = lvc->VoxelOffsetY() * CLHEP::cm;
269  G4double voxelOffsetZ = lvc->VoxelOffsetZ() * CLHEP::cm;
270 
271  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
272  << ": voxelSizeX=" << voxelSizeX << ", voxelSizeY=" << voxelSizeY
273  << ", voxelSizeZ=" << voxelSizeZ;
274 
275  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
276  << ": voxelOffsetX=" << voxelOffsetX << ", voxelOffsetY=" << voxelOffsetY
277  << ", voxelOffsetZ=" << voxelOffsetZ;
278 
279  // We want our voxelization region to be an integer multiple of
280  // the voxel sizes in all directions; if we didn't do this, we
281  // might get into trouble when we start playing with replicas.
282  // Compute the the dimensions of our voxelization to be about the
283  // size of the LAr TPC region, adjusted to be an integer number of
284  // voxels in all directions.
285 
286  G4double numberXvoxels = 2. * larTPCHalfXLength / voxelSizeX;
287  G4double numberYvoxels = 2. * larTPCHalfYLength / voxelSizeY;
288  G4double numberZvoxels = 2. * larTPCHalfZLength / voxelSizeZ;
289  numberXvoxels = trunc(numberXvoxels) + 1.;
290  numberYvoxels = trunc(numberYvoxels) + 1.;
291  numberZvoxels = trunc(numberZvoxels) + 1.;
292 
293  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
294  << "Active volume of cryo #" << c << " TPC #" << t << " will be split in "
295  << numberXvoxels << " x " << numberYvoxels << " x " << numberZvoxels << " = "
296  << (numberXvoxels * numberYvoxels * numberZvoxels) << " voxels of size " << voxelSizeX
297  << " x " << voxelSizeY << " x " << voxelSizeZ << " cm";
298 
299  VoxelSpecs_t VoxelSpecs{/* w */ 2. * larTPCHalfXLength,
300  /* h */ 2. * larTPCHalfYLength,
301  /* d */ 2. * larTPCHalfZLength,
302  /* nw */ (unsigned int)numberXvoxels,
303  /* nh */ (unsigned int)numberYvoxels,
304  /* nd */ (unsigned int)numberZvoxels};
305 
306  VoxelCache_t::iterator iVoxelVol =
307  DisableVoxelCaching ? VoxelCache.end() : VoxelCache.find(VoxelSpecs);
308  if (iVoxelVol == VoxelCache.end()) {
309  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
310  << "Creating a new voxel volume " << VoxelSpecs.w << " x " << VoxelSpecs.h << " x "
311  << VoxelSpecs.d << " cm, " << VoxelSpecs.nw << " x " << VoxelSpecs.nh << " x "
312  << VoxelSpecs.nd << " voxels";
313 
314  G4double voxelBoxHalfX = numberXvoxels * voxelSizeX / 2.;
315  G4double voxelBoxHalfY = numberYvoxels * voxelSizeY / 2.;
316  G4double voxelBoxHalfZ = numberZvoxels * voxelSizeZ / 2.;
317 
318  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
319  << ": voxelBoxHalfX=" << voxelBoxHalfX << ", voxelBoxHalfY=" << voxelBoxHalfY
320  << ", voxelBoxHalfZ=" << voxelBoxHalfZ;
321 
322  // Now we have a box that will include an integer number of voxels
323  // in each direction. Note that the material is irrelevant for a
324  // "parallel world."
325  auto voxelBox = new G4Box("VoxelBox", voxelBoxHalfX, voxelBoxHalfY, voxelBoxHalfZ);
326  auto voxelBoxLogical = new G4LogicalVolume(voxelBox, 0, "VoxelizationLogicalVolume");
327 
328  // If we generate an event display within Geant4, we won't want to
329  // see this box.
330  auto invisible = new G4VisAttributes();
331  invisible->SetVisibility(false);
332  voxelBoxLogical->SetVisAttributes(invisible);
333 
334  // Now we've fill our "box of voxels" with the voxels themselves.
335  // We'll do this by sub-dividing the volume in x, then y, then z.
336 
337  // Create an "x-slice".
338  G4Box* xSlice = new G4Box("xSlice", voxelSizeX / 2., voxelBoxHalfY, voxelBoxHalfZ);
339  G4LogicalVolume* xSliceLogical = new G4LogicalVolume(xSlice, 0, "xLArVoxelSlice");
340  xSliceLogical->SetVisAttributes(invisible);
341 
342  // Use replication to slice up the "box of voxels" along the x-axis.
343  new G4PVReplica("VoxelSlicesInX",
344  xSliceLogical,
345  voxelBoxLogical,
346  kXAxis,
347  G4int(numberXvoxels),
348  voxelSizeX);
349 
350  // Now do the same thing, dividing that x-slice along the y-axis.
351  auto ySlice = new G4Box("ySlice", voxelSizeX / 2., voxelSizeY / 2., voxelBoxHalfZ);
352  auto ySliceLogical = new G4LogicalVolume(ySlice, 0, "yLArVoxelSlice");
353  ySliceLogical->SetVisAttributes(invisible);
354  new G4PVReplica("VoxelSlicesInY",
355  ySliceLogical,
356  xSliceLogical,
357  kYAxis,
358  G4int(numberYvoxels),
359  voxelSizeY);
360 
361  // Now divide the y-slice along the z-axis, giving us our actual
362  // voxels.
363  auto zSlice = new G4Box("zSlice", voxelSizeX / 2., voxelSizeY / 2., voxelSizeZ / 2.);
364  auto voxelLogical = new G4LogicalVolume(zSlice, 0, "LArVoxel");
365  voxelLogical->SetVisAttributes(invisible);
366  new G4PVReplica(
367  "LArVoxel", voxelLogical, ySliceLogical, kZAxis, G4int(numberZvoxels), voxelSizeZ);
368 
369  iVoxelVol = VoxelCache.emplace(VoxelSpecs, VoxelVolumes_t(voxelBoxLogical, voxelLogical))
370  .first; // iterator to the inserted element
371 
372  // Set the sensitive detector of this LAr TPC (logical) volume
373  // to be the voxel readout.
374  voxelLogical->SetSensitiveDetector(flarVoxelReadout);
375 
376  } // if not cached yet
377  G4LogicalVolume* voxelBoxLogical = iVoxelVol->second.pBox;
378 
379  // We have a "box of voxels" that's the right dimensions, but we
380  // have to know exactly where to put it. The user has the option
381  // to offset the voxel co-ordinate system. We want to place our
382  // box so the edges of our voxels align with that co-ordinate
383  // system. In effect, we want to offset our "box of voxels" by
384  // the user's offsets, modulo the size of the voxel in each
385  // direction.
386 
387  G4double offsetInVoxelsX = voxelOffsetX / voxelSizeX;
388  G4double offsetInVoxelsY = voxelOffsetY / voxelSizeY;
389  G4double offsetInVoxelsZ = voxelOffsetZ / voxelSizeZ;
390  G4double fractionOffsetX = offsetInVoxelsX - trunc(offsetInVoxelsX);
391  G4double fractionOffsetY = offsetInVoxelsY - trunc(offsetInVoxelsY);
392  G4double fractionOffsetZ = offsetInVoxelsZ - trunc(offsetInVoxelsZ);
393  G4double offsetX = fractionOffsetX * voxelSizeX;
394  G4double offsetY = fractionOffsetY * voxelSizeY;
395  G4double offsetZ = fractionOffsetZ * voxelSizeZ;
396 
397  // Now we know how much to offset the "box of voxels". Include
398  // that in the transformation of the co-ordinates from world
399  // volume to LAr TPC volume.
400  transform = G4Translate3D(offsetX, offsetY, offsetZ) * transform;
401 
402  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
403  << ": offsetX=" << offsetX << ", offsetY=" << offsetY << ", offsetZ=" << offsetZ;
404 
405  // suffix representing this TPC
406  std::ostringstream nums;
407  nums << "_Cryostat" << c << "_TPC" << t;
408 
409  // Place the box of voxels, with the accumulated transformations
410  // computed above.
411  new G4PVPlacementInTPC(transform,
412  "VoxelizationPhysicalVolume" + nums.str(),
413  voxelBoxLogical,
414  parallelPhysical,
415  false, // Only one volume
416  0, // Copy number
417  false, // No surface check
418  {(unsigned short int)c, (unsigned short int)t} // TPC ID
419  );
420 
421  } // end loop over tpcs
422  } // end loop over cryostats
423 
424  if (mf::isDebugEnabled()) {
425  MF_LOG_DEBUG("LArVoxelReadoutGeometryDump") << "Dump of voxelized volume";
426  {
427  mf::LogDebug log("LArVoxelReadoutGeometryDump");
428  DumpPhysicalVolume(log, *parallelPhysical);
429  }
430  MF_LOG_DEBUG("LArVoxelReadoutGeometryDump") << "End of dump of voxelized volume";
431  }
432  }
static constexpr Sample_t transform(Sample_t sample)
larg4::LArVoxelReadout::Setup_t fReadoutSetupData
constexpr bool operator<(CryostatID const &a, CryostatID const &b)
Order cryostats with increasing ID.
Definition: geo_types.h:706
void Setup(Setup_t const &setupData)
Reads all the configuration elements from setupData
larg4::LArVoxelReadout * flarVoxelReadout
Data for LArVoxelReadout setup.
while getopts h
void SetSingleTPC(unsigned int cryostat, unsigned int tpc)
Associates this readout to one specific TPC.
std::unique_ptr< G4UserLimits > fStepLimit
constexpr bool DisableVoxelCaching
G4PVPlacementWithID< TPCID_t > G4PVPlacementInTPC
A physical volume with a TPC ID.
art::ServiceHandle< geo::Geometry const > fGeo
Handle to the geometry.
G4VPhysicalVolume * FindNestedVolume(G4VPhysicalVolume *mother, G4Transform3D &motherTransform, G4Transform3D &daughterTransform, std::string &daughterName, unsigned int expectedNum)
G4VPhysicalVolume * larg4::LArVoxelReadoutGeometry::FindNestedVolume ( G4VPhysicalVolume *  mother,
G4Transform3D &  motherTransform,
G4Transform3D &  daughterTransform,
std::string &  daughterName,
unsigned int  expectedNum 
)
private

Definition at line 442 of file LArVoxelReadoutGeometry.cxx.

447  {
448  G4LogicalVolume* logicalVolume = mother->GetLogicalVolume();
449  G4int numberDaughters = logicalVolume->GetNoDaughters();
450  for (G4int i = 0; i != numberDaughters; ++i) {
451  G4VPhysicalVolume* d = logicalVolume->GetDaughter(i);
452 
453  if (d->GetName().contains(daughterName)) {
454 
455  // check that this cryostat is the requested one using fCryostat
456  G4ThreeVector translation = d->GetObjectTranslation();
457  G4RotationMatrix rotation = d->GetObjectRotationValue();
458  G4Transform3D transform(rotation, translation);
459  daughterTransform = motherTransform * transform;
460 
461  // take the origin of the volume and transform it to
462  // world coordinated
463  G4Point3D local(0., 0., 0.);
464  G4Point3D world = daughterTransform * local;
465 
466  MF_LOG_DEBUG("LArVoxelReadoutGeometry")
467  << "current " << daughterName << " origin is at (" << world.x() / CLHEP::cm << ","
468  << world.y() / CLHEP::cm << "," << world.z() / CLHEP::cm << ")";
469 
470  // we don't bother with the cryostat number when calling Geometry::PositionToTPC
471  // because we know we have already started off with the correct cryostat volume
472  // G4 uses mm, we want cm
473  double worldPos[3] = {world.x() / CLHEP::cm, world.y() / CLHEP::cm, world.z() / CLHEP::cm};
474  unsigned int daughterNum = 0;
475  unsigned int extra = 0;
476  if (daughterName.compare("volCryostat") == 0)
477  fGeo->PositionToCryostat(worldPos, daughterNum);
478  else if (daughterName.compare("volTPC") == 0)
479  fGeo->PositionToTPC(worldPos, daughterNum, extra);
480  else if (daughterName.compare("volTPCActive") == 0 ||
481  daughterName.compare("volDetEnclosure") == 0) {
482  // for either of these volumes, we know there is only 1 in the mother volume
483  MF_LOG_DEBUG("LArVoxelReadoutGeometry") << "found the desired " << daughterName;
484  return d;
485  }
486 
487  // if we found the desired volume, stop looking
488  if (daughterNum == expectedNum) {
489  MF_LOG_DEBUG("LArVoxelReadoutGeometry") << "found the desired " << daughterName;
490  return d;
491  }
492  } // end if the volume has the right name
493  } // end loop over volumes
494 
495  throw cet::exception("LArVoxelReadoutGeometry")
496  << "could not find the desired " << daughterName << " to make LArVoxelReadoutGeometry\n";
497 
498  return 0;
499  }
static constexpr Sample_t transform(Sample_t sample)
then local
art::ServiceHandle< geo::Geometry const > fGeo
Handle to the geometry.

Member Data Documentation

art::ServiceHandle<geo::Geometry const> larg4::LArVoxelReadoutGeometry::fGeo
private

Handle to the geometry.

Definition at line 71 of file LArVoxelReadoutGeometry.h.

larg4::LArVoxelReadout* larg4::LArVoxelReadoutGeometry::flarVoxelReadout {nullptr}
private

Data for LArVoxelReadout setup.

Definition at line 76 of file LArVoxelReadoutGeometry.h.

larg4::LArVoxelReadout::Setup_t larg4::LArVoxelReadoutGeometry::fReadoutSetupData
private

Definition at line 77 of file LArVoxelReadoutGeometry.h.

std::unique_ptr<G4UserLimits> larg4::LArVoxelReadoutGeometry::fStepLimit
private

G4 doesn't handle memory management, so we have to

Definition at line 72 of file LArVoxelReadoutGeometry.h.


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