LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
TPCGeo.cxx
Go to the documentation of this file.
1 
7 // class header
9 
10 // LArSoft includes
14 
15 // Framework includes
17 #include "cetlib_except/exception.h"
18 
19 // ROOT includes
20 #include "TGeoNode.h"
21 #include "TGeoMatrix.h"
22 #include "TGeoBBox.h"
23 
24 // C/C++ standard libraries
25 #include <cmath>
26 #include <cassert>
27 #include <map>
28 #include <algorithm> // std::max(), std::copy()
29 #include <iterator> // std::inserter()
30 #include <functional> // std::mem_fn()
31 
32 
33 namespace geo{
34 
35 
36  //......................................................................
37  TPCGeo::TPCGeo(GeoNodePath_t& path, size_t depth)
38  : BoxBoundedGeo() // we initialize boundaries at the end of construction
39  , fTrans(path, depth)
40  , fActiveVolume(0)
41  , fTotalVolume(0)
42  , fDriftDirection(geo::kUnknownDrift)
43  , fWidthDir (geo::Xaxis())
44  , fHeightDir(geo::Yaxis())
45  , fLengthDir(geo::Zaxis())
46  , fDriftDir() // null until known
47  {
48 
49  // all planes are going to be contained in the volume named volTPC
50  // now get the total volume of the TPC
51  TGeoVolume *vc = path[depth]->GetVolume();
52  if(!vc){
53  throw cet::exception("Geometry") << "cannot find detector outline volume - bail ungracefully\n";
54  }
55 
56  fTotalVolume = vc;
57 
58  // loop over the daughters of this node and look for the active volume
59  int nd = vc->GetNdaughters();
60  TGeoNode const* pActiveVolNode = nullptr;
61  for(int i = 0; i < nd; ++i){
62  if(strncmp(vc->GetNode(i)->GetName(), "volTPCActive", 12) != 0) continue;
63 
64  pActiveVolNode = vc->GetNode(i);
65  TGeoVolume *vca = pActiveVolNode->GetVolume();
66  if(vca) fActiveVolume = vca;
67  break;
68 
69  }// end loop over daughters of the volume
70 
72 
73  LOG_DEBUG("Geometry") << "detector total volume is " << fTotalVolume->GetName()
74  << "\ndetector active volume is " << fActiveVolume->GetName();
75 
76  // compute the active volume transformation too
77  TGeoHMatrix ActiveHMatrix(fTrans.Matrix());
78  if (pActiveVolNode) ActiveHMatrix.Multiply(pActiveVolNode->GetMatrix());
79  // we don't keep the active volume information... just store its center:
80  std::array<double, 3> localActiveCenter, worldActiveCenter;
81  localActiveCenter.fill(0.0);
82  ActiveHMatrix.LocalToMaster
83  (localActiveCenter.data(), worldActiveCenter.data());
84  fActiveCenter = geo::vect::makeFromCoords<geo::Point_t>(worldActiveCenter);
85 
86 
87  // find the wires for the plane so that you can use them later
88  this->FindPlane(path, depth);
89 
90  // set the width, height, and lengths
91  fActiveHalfWidth = ((TGeoBBox*)fActiveVolume->GetShape())->GetDX();
92  fActiveHalfHeight = ((TGeoBBox*)fActiveVolume->GetShape())->GetDY();
93  fActiveLength = 2.0*((TGeoBBox*)fActiveVolume->GetShape())->GetDZ();
94 
95  fHalfWidth = ((TGeoBBox*)fTotalVolume->GetShape())->GetDX();
96  fHalfHeight = ((TGeoBBox*)fTotalVolume->GetShape())->GetDY();
97  fLength = 2.0*((TGeoBBox*)fTotalVolume->GetShape())->GetDZ();
98 
99  // check that the rotation matrix to the world is the identity, if not
100  // we need to change the width, height and length values;
101  // the correspondence of these to x, y and z are not guaranteed to be
102  // trivial, so we store the two independently (cartesian dimensions in the
103  // bounding boxes, the sizes in data members directly)
104  double const* rotMatrix = fTrans.Matrix().GetRotationMatrix();
105  if(rotMatrix[0] != 1){
106  if(std::abs(rotMatrix[2]) == 1){
107  fActiveHalfWidth = ((TGeoBBox*)fActiveVolume->GetShape())->GetDZ();
108  fHalfWidth = ((TGeoBBox*)fTotalVolume->GetShape())->GetDZ();
109  fWidthDir = Zaxis();
110  }
111  if(std::abs(rotMatrix[1]) == 1){
112  fActiveHalfWidth = ((TGeoBBox*)fActiveVolume->GetShape())->GetDY();
113  fHalfWidth = ((TGeoBBox*)fTotalVolume->GetShape())->GetDY();
114  fWidthDir = Yaxis();
115  }
116  }
117  if(rotMatrix[4] != 1){
118  if(std::abs(rotMatrix[3]) == 1){
119  fActiveHalfHeight = ((TGeoBBox*)fActiveVolume->GetShape())->GetDX();
120  fHalfHeight = ((TGeoBBox*)fTotalVolume->GetShape())->GetDX();
121  fHeightDir = Xaxis();
122  }
123  if(std::abs(rotMatrix[5]) == 1){
124  fActiveHalfHeight = ((TGeoBBox*)fActiveVolume->GetShape())->GetDZ();
125  fHalfHeight = ((TGeoBBox*)fTotalVolume->GetShape())->GetDZ();
126  fHeightDir = Zaxis();
127  }
128  }
129  if(rotMatrix[8] != 1){
130  if(std::abs(rotMatrix[6]) == 1){
131  fActiveLength = 2.*((TGeoBBox*)fActiveVolume->GetShape())->GetDX();
132  fLength = 2.*((TGeoBBox*)fTotalVolume->GetShape())->GetDX();
133  fLengthDir = Xaxis();
134  }
135  if(std::abs(rotMatrix[7]) == 1){
136  fActiveLength = 2.*((TGeoBBox*)fActiveVolume->GetShape())->GetDY();
137  fLength = 2.*((TGeoBBox*)fTotalVolume->GetShape())->GetDY();
138  fLengthDir = Yaxis();
139  }
140  }
141 
144 
145  } // TPCGeo::TPCGeo()
146 
147  //......................................................................
148  void TPCGeo::FindPlane(GeoNodePath_t& path, size_t depth)
149  {
150 
151  const char* nm = path[depth]->GetName();
152  if( (strncmp(nm, "volTPCPlane", 11) == 0) ){
153  this->MakePlane(path,depth);
154  return;
155  }
156 
157  //explore the next layer down
158  unsigned int deeper = depth+1;
159  if(deeper >= path.size()){
160  throw cet::exception("BadTGeoNode") << "exceeded maximum TGeoNode depth\n";
161  }
162 
163  const TGeoVolume *v = path[depth]->GetVolume();
164  int nd = v->GetNdaughters();
165  for(int i = 0; i < nd; ++i){
166  path[deeper] = v->GetNode(i);
167  this->FindPlane(path, deeper);
168  }
169 
170  }
171 
172 
173  //......................................................................
174  void TPCGeo::MakePlane(GeoNodePath_t& path, size_t depth)
175  {
176  fPlanes.emplace_back(path, depth);
177  }
178 
179 
180  //......................................................................
181  short int TPCGeo::DetectDriftDirection() const {
182 
183  //
184  // 1. determine the drift axis
185  // 2. determine the drift direction on it
186  //
187  // We assume that all the planes cover most of the TPC face; therefore,
188  // the centre of the plane and the one of the TPC should be very close
189  // to each other, when projected on the same drift plane.
190  // Here we find which is the largest coordinate difference.
191 
192  if (Nplanes() == 0) {
193  // chances are that we get this because stuff is not initialised yet,
194  // and then even the ID might be wrong
195  throw cet::exception("TPCGeo")
196  << "DetectDriftDirection(): no planes in TPC " << std::string(ID())
197  << "\n";
198  }
199 
200  auto const TPCcenter = GetCenter();
201  auto const PlaneCenter = Plane(0).GetBoxCenter(); // any will do
202 
203  auto const driftVector = PlaneCenter - TPCcenter; // approximation!
204 
205  if ((std::abs(driftVector.X()) > std::abs(driftVector.Y()))
206  && (std::abs(driftVector.X()) > std::abs(driftVector.Z())))
207  {
208  // x is the solution
209  return (driftVector.X() > 0)? +1: -1;
210  }
211  else if (std::abs(driftVector.Y()) > std::abs(driftVector.Z()))
212  {
213  // y is the man
214  return (driftVector.Y() > 0)? +2: -2;
215  }
216  else {
217  // z is the winner
218  return (driftVector.Z() > 0)? +3: -3;
219  }
220 
221  } // TPCGeo::DetectDriftDirection()
222 
223  //......................................................................
224  // sort the PlaneGeo objects and the WireGeo objects inside
226  {
228 
229  double origin[3] = {0.};
230 
231  // set the plane pitch for this TPC
232  double xyz[3] = {0.};
233  fPlanes[0].LocalToWorld(origin,xyz);
234  double xyz1[3] = {0.};
235  fPlaneLocation.clear();
236  fPlaneLocation.resize(fPlanes.size());
237  for(unsigned int i = 0; i < fPlaneLocation.size(); ++i) fPlaneLocation[i].resize(3);
238  fPlane0Pitch.clear();
239  fPlane0Pitch.resize(this->Nplanes(), 0.);
240  for(size_t p = 0; p < this->Nplanes(); ++p){
241  fPlanes[p].LocalToWorld(origin,xyz1);
242  if(p > 0) fPlane0Pitch[p] = fPlane0Pitch[p-1] + std::abs(xyz1[0]-xyz[0]);
243  else fPlane0Pitch[p] = 0.;
244  xyz[0] = xyz1[0];
245  fPlaneLocation[p][0] = xyz1[0];
246  fPlaneLocation[p][1] = xyz1[1];
247  fPlaneLocation[p][2] = xyz1[2];
248  }
249 
250  // the PlaneID_t cast convert InvalidID into a rvalue (non-reference);
251  // leaving it a reference would cause C++ to treat it as such,
252  // that can't be because InvalidID is a static member constant without an address
253  // (it is not defined in any translation unit, just declared in header)
254  fViewToPlaneNumber.resize
256  for(size_t p = 0; p < this->Nplanes(); ++p)
257  fViewToPlaneNumber[(size_t) fPlanes[p].View()] = p;
258 
259  for(size_t p = 0; p < fPlanes.size(); ++p) fPlanes[p].SortWires(sorter);
260 
261  }
262 
263 
264  //......................................................................
266 
267  // reset the ID
268  fID = tpcid;
269 
270  // ask the planes to update; also check
271 
272  for (unsigned int plane = 0; plane < Nplanes(); ++plane) {
273  fPlanes[plane].UpdateAfterSorting(geo::PlaneID(fID, plane), *this);
274 
275  // check that the plane normal is opposite to the TPC drift direction
277  .equal(-(fPlanes[plane].GetNormalDirection()), DriftDir()));
278 
279  } // for
280 
282 
283  } // TPCGeo::UpdateAfterSorting()
284 
285 
286  //......................................................................
287  const PlaneGeo& TPCGeo::Plane(unsigned int iplane) const
288  {
289  geo::PlaneGeo const* pPlane = PlanePtr(iplane);
290  if (!pPlane){
291  throw cet::exception("PlaneOutOfRange") << "Request for non-existant plane " << iplane << "\n";
292  }
293 
294  return *pPlane;
295  }
296 
297  //......................................................................
298  const PlaneGeo& TPCGeo::Plane(geo::View_t view) const
299  {
300  geo::PlaneID::PlaneID_t const p = fViewToPlaneNumber[size_t(view)];
301  if (p == geo::PlaneID::InvalidID) {
302  throw cet::exception("TPCGeo")
303  << "TPCGeo[" << ((void*) this) << "]::Plane(): no plane for view #"
304  << (size_t) view << "\n";
305  }
306  return fPlanes[p];
307  } // TPCGeo::Plane(geo::View_t)
308 
309 
310  //......................................................................
312 
313  //
314  // Returns the plane with the smallest width x depth. No nonsense here.
315  //
316 
317  auto iPlane = fPlanes.begin(), pend = fPlanes.end();
318  auto smallestPlane = iPlane;
319  double smallestSurface = smallestPlane->Width() * smallestPlane->Depth();
320  while (++iPlane != pend) {
321  double const surface = iPlane->Width() * iPlane->Depth();
322  if (surface > smallestSurface) continue;
323  smallestSurface = surface;
324  smallestPlane = iPlane;
325  } // while
326  return *smallestPlane;
327 
328  } // TPCGeo::SmallestPlane()
329 
330 
331  //......................................................................
332  unsigned int TPCGeo::MaxWires() const {
333  unsigned int maxWires = 0;
334  for (geo::PlaneGeo const& plane: fPlanes) {
335  unsigned int maxWiresInPlane = plane.Nwires();
336  if (maxWiresInPlane > maxWires) maxWires = maxWiresInPlane;
337  } // for
338  return maxWires;
339  } // TPCGeo::MaxWires()
340 
341 
342  //......................................................................
343  std::set<geo::View_t> TPCGeo::Views() const {
344  std::set<geo::View_t> views;
345  std::transform(fPlanes.cbegin(), fPlanes.cend(),
346  std::inserter(views, views.begin()), std::mem_fn(&geo::PlaneGeo::View));
347  return views;
348  } // TPCGeo::Views()
349 
350 
351  //......................................................................
352  // returns distance between plane 0 to each of the remaining planes
353  // not the distance between two consecutive planes
354  double TPCGeo::Plane0Pitch(unsigned int p) const
355  {
356  return fPlane0Pitch[p];
357  }
358 
359  //......................................................................
361 
362  //
363  // 1. find the center of the face of the TPC opposite to the anode
364  // 2. compute the distance of it from the last wire plane
365  //
366 
367  //
368  // find the cathode center
369  //
370  geo::Point_t cathodeCenter = GetActiveVolumeCenter<geo::Point_t>();
371  switch (DetectDriftDirection()) {
372  case -1:
373  geo::vect::Xcoord(cathodeCenter) += ActiveHalfWidth();
374  // cathodeCenter.SetX(cathodeCenter.X() + ActiveHalfWidth());
375  break;
376  case +1:
377  cathodeCenter.SetX(cathodeCenter.X() - ActiveHalfWidth());
378  break;
379  case -2:
380  cathodeCenter.SetY(cathodeCenter.Y() + ActiveHalfHeight());
381  break;
382  case +2:
383  cathodeCenter.SetY(cathodeCenter.Y() - ActiveHalfHeight());
384  break;
385  case -3:
386  cathodeCenter.SetZ(cathodeCenter.Z() + ActiveLength() / 2.0);
387  break;
388  case +3:
389  cathodeCenter.SetZ(cathodeCenter.Z() - ActiveLength() / 2.0);
390  break;
391  case 0:
392  default:
393  // in this case, a better algorithm is probably needed
394  throw cet::exception("TPCGeo")
395  << "CathodeCenter(): Can't determine the cathode plane (code="
396  << DetectDriftDirection() << ")\n";
397  } // switch
398  return cathodeCenter;
399 
400  } // TPCGeo::GetCathodeCenterImpl()
401 
402 
403  //......................................................................
405  auto const& activeBox = ActiveBoundingBox();
406  return { activeBox.CenterX(), activeBox.CenterY(), activeBox.MinZ() };
407  } // TPCGeo::GetFrontFaceCenterImpl()
408 
409 
410  //......................................................................
411  // returns xyz location of planes in TPC
412  const double* TPCGeo::PlaneLocation(unsigned int p) const
413  {
414  return &fPlaneLocation[p][0];
415  }
416 
417  //......................................................................
418  double TPCGeo::PlanePitch(unsigned int p1,
419  unsigned int p2) const
420  {
421  return std::abs(fPlane0Pitch[p2] - fPlane0Pitch[p1]);
422  }
423 
424  //......................................................................
425  // This method returns the distance between wires in the given plane.
426  double TPCGeo::WirePitch(unsigned plane) const
427  {
428  return this->Plane(plane).WirePitch();
429  }
430 
431  //......................................................................
433 
434  auto const driftDirCode = DetectDriftDirection();
435  switch (driftDirCode) {
436  case +1:
437  fDriftDirection = geo::kPosX; // this is the same as kPos!
438  fDriftDir = geo::Xaxis();
439  break;
440  case -1:
441  fDriftDirection = geo::kNegX; // this is the same as kNeg!
442  fDriftDir = -geo::Xaxis();
443  break;
444  case +2:
445  fDriftDir = geo::Yaxis();
447  break;
448  case -2:
449  fDriftDir = -geo::Yaxis();
451  break;
452  case +3:
453  fDriftDir = geo::Zaxis();
455  break;
456  case -3:
457  fDriftDir = -geo::Zaxis();
459  break;
460  default:
461  // TPC ID is likely not yet set
463 
464  // we estimate the drift direction roughly from the geometry
466 
467  mf::LogError("TPCGeo")
468  << "Unable to detect drift direction (result: " << driftDirCode
469  << ", drift: ( " << fDriftDir.X() << " ; " << fDriftDir.Y() << " ; "
470  << fDriftDir.Z() << " )";
471  break;
472  } // switch
473 
475 
476  } // TPCGeo::ResetDriftDirection()
477 
478 
479  //......................................................................
481 
482  //
483  // 1. find the center of the face of the TPC opposite to the anode
484  // 2. compute the distance of it from the last wire plane
485  //
486 
487  geo::PlaneGeo const& plane = fPlanes.back();
488  return std::abs(plane.DistanceFromPlane(GetCathodeCenter()));
489 
490  } // TPCGeo::ComputeDriftDistance()
491 
492 
493  //......................................................................
495  // note that this assumes no rotations of the TPC
496  // (except for rotations of a flat angle around one of the three main axes);
497  // to avoid this, we should transform the six vertices
498  // rather than just the centre
499 
500  // we rely on the asumption that the center of TPC is at the local origin
504  );
505 
506  // the center of the active volume may be elsewhere than the local origin:
507  auto const& activeCenter = GetActiveVolumeCenter<geo::Point_t>();
509  activeCenter.X() - ActiveHalfWidth(), activeCenter.X() + ActiveHalfWidth(),
510  activeCenter.Y() - ActiveHalfHeight(), activeCenter.Y() + ActiveHalfHeight(),
511  activeCenter.Z() - ActiveHalfLength(), activeCenter.Z() + ActiveHalfLength()
512  );
513 
514 
515  } // CryostatGeo::InitTPCBoundaries()
516 
517  //......................................................................
518 
520 
521  // the PlaneID_t cast convert InvalidID into a rvalue (non-reference);
522  // leaving it a reference would cause C++ to treat it as such,
523  // that can't be because InvalidID is a static member constant without an address
524  // (it is not defined in any translation unit, just declared in header)
525  fViewToPlaneNumber.clear();
526  fViewToPlaneNumber.resize
528  for(size_t p = 0; p < Nplanes(); ++p)
529  fViewToPlaneNumber[(size_t) fPlanes[p].View()] = p;
530 
531  } // TPCGeo::UpdatePlaneViewCache()
532 
533 
534  //......................................................................
535  void TPCGeo::SortPlanes(std::vector<geo::PlaneGeo>& planes) const {
536  //
537  // Sort planes by increasing drift distance.
538  //
539  // This function should work in bootstrap mode, relying on least things as
540  // possible. Therefore we compute here a proxy of the drift axis.
541  //
542 
543  //
544  // determine the drift axis (or close to): from TPC center to plane center
545  //
546 
547  // Instead of using the plane center, which might be not available yet,
548  // we use the plane box center, which only needs the geometry description
549  // to be available.
550  // We use the first plane -- it does not make any difference.
551  decltype(auto) TPCcenter = GetCenter();
552  auto driftAxis
553  = geo::vect::normalize(planes[0].GetBoxCenter() - TPCcenter);
554 
555  //
556  // associate each plane with its distance from the center of TPC
557  //
558  decltype(auto) center = GetCenter();
559  // pair: <plane pointer,drift distance>
560  std::vector<std::pair<geo::PlaneGeo*, double>> planesWithDistance;
561  planesWithDistance.reserve(planes.size());
562  for (geo::PlaneGeo& plane: planes) {
563  double const driftDistance
564  = geo::vect::dot(plane.GetBoxCenter() - TPCcenter, driftAxis);
565  planesWithDistance.emplace_back(&plane, driftDistance);
566  } // for
567 
568  //
569  // sort by distance
570  //
571  std::sort(planesWithDistance.begin(), planesWithDistance.end(),
572  [](auto const& a, auto const& b){ return a.second < b.second; }
573  );
574 
575  //
576  // extract the result
577  //
578  std::vector<geo::PlaneGeo> sortedPlanes;
579  sortedPlanes.reserve(planesWithDistance.size());
580  for (auto const& pair: planesWithDistance)
581  sortedPlanes.push_back(std::move(*(pair.first)));
582 
583  planes = std::move(sortedPlanes);
584 
585  } // TPCGeo::SortPlanes()
586 
587  //......................................................................
588 
589 }
geo::TPCID const & ID() const
Returns the identifier of this TPC.
Definition: TPCGeo.h:278
void InitTPCBoundaries()
Recomputes the TPC boundary.
Definition: TPCGeo.cxx:494
geo::Point_t GetCathodeCenterImpl() const
Definition: TPCGeo.cxx:360
void round01(Vector &v, Scalar tol)
Returns a vector with all components rounded if close to 0, -1 or +1.
Vector DriftDir() const
Returns the direction of the drift (vector pointing toward the planes).
Definition: TPCGeo.h:687
geo::BoxBoundedGeo fActiveBox
Box of the active volume.
Definition: TPCGeo.h:652
std::vector< double > fPlane0Pitch
Pitch between planes.
Definition: TPCGeo.h:636
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:418
void MakePlane(GeoNodePath_t &path, size_t depth)
Definition: TPCGeo.cxx:174
constexpr auto dot(Vector const &a, Vector const &b)
Return cross product of two vectors.
Drift direction is unknown.
Definition: geo_types.h:105
void FindPlane(GeoNodePath_t &path, size_t depth)
Definition: TPCGeo.cxx:148
double fLength
Length of total volume.
Definition: TPCGeo.h:645
TPCGeo(GeoNodePath_t &path, size_t depth)
Definition: TPCGeo.cxx:37
geo::Point3DBase_t< TPCGeoCoordinatesTag > LocalPoint_t
Type of points in the local GDML TPC frame.
Definition: TPCGeo.h:65
double ComputeDriftDistance() const
Definition: TPCGeo.cxx:480
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
double ActiveHalfHeight() const
Half height (associated with y coordinate) of active TPC volume [cm].
Definition: TPCGeo.h:91
Drift towards positive values.
Definition: geo_types.h:106
Unknown view.
Definition: geo_types.h:83
unsigned int Nplanes() const
Number of planes in this tpc.
Definition: TPCGeo.h:145
double fActiveHalfWidth
Half width of active volume.
Definition: TPCGeo.h:640
Point GetCathodeCenter() const
Definition: TPCGeo.h:253
TransformationMatrix_t const & Matrix() const
Direct access to the transformation matrix.
geo::Point_t GetFrontFaceCenterImpl() const
Definition: TPCGeo.cxx:404
The data type to uniquely identify a Plane.
Definition: geo_types.h:250
constexpr Vector Yaxis()
Returns a y axis vector of the specified type.
Definition: geo_vectors.h:222
Point GetBoxCenter() const
Returns the centre of the box representing the plane.
Definition: PlaneGeo.h:441
geo::BoxBoundedGeo const & ActiveBoundingBox() const
Returns the box of the active volume of this TPC.
Definition: TPCGeo.h:266
std::vector< geo::PlaneID::PlaneID_t > fViewToPlaneNumber
Index of the plane for each view (InvalidID if none).
Definition: TPCGeo.h:657
std::vector< PlaneGeo > fPlanes
List of planes in this plane.
Definition: TPCGeo.h:632
double HalfLength() const
Length is associated with z coordinate [cm].
Definition: TPCGeo.h:109
void UpdateAfterSorting(geo::TPCID tpcid)
Performs all updates after cryostat has sorted TPCs.
Definition: TPCGeo.cxx:265
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
Drift towards negative X values.
Definition: geo_types.h:109
void SortSubVolumes(geo::GeoObjectSorter const &sorter)
Apply sorting to the PlaneGeo objects.
Definition: TPCGeo.cxx:225
geo::Vector_t fHeightDir
Direction height refers to.
Definition: TPCGeo.h:648
geo::Point_t toWorldCoords(LocalPoint_t const &local) const
Transform point from local TPC frame to world frame.
Definition: TPCGeo.h:498
geo::TPCID fID
ID of this TPC.
Definition: TPCGeo.h:654
Class for approximate comparisons.
unsigned int PlaneID_t
Type for the ID number.
Definition: geo_types.h:251
View_t View() const
Which coordinate does this plane measure.
Definition: PlaneGeo.h:171
geo::Vector_t fWidthDir
Direction width refers to.
Definition: TPCGeo.h:647
double ActiveHalfLength() const
Length (associated with z coordinate) of active TPC volume [cm].
Definition: TPCGeo.h:97
double fHalfWidth
Half width of total volume.
Definition: TPCGeo.h:643
auto makeVector3DComparison(RealType threshold)
Creates a Vector3DComparison from a RealComparisons object.
geo::PlaneGeo const & SmallestPlane() const
Returns the wire plane with the smallest surface.
Definition: TPCGeo.cxx:311
double fActiveLength
Length of active volume.
Definition: TPCGeo.h:642
double ActiveHalfWidth() const
Half width (associated with x coordinate) of active TPC volume [cm].
Definition: TPCGeo.h:87
Geometry information for a single wire plane.The plane is represented in the geometry by a solid whic...
Definition: PlaneGeo.h:78
void UpdatePlaneViewCache()
Refills the plane vs. view cache of the TPC.
Definition: TPCGeo.cxx:519
constexpr Vector Xaxis()
Returns a x axis vector of the specified type.
Definition: geo_vectors.h:218
TGeoVolume * fTotalVolume
Total volume of TPC, called volTPC in GDML file.
Definition: TPCGeo.h:634
DriftDirection_t fDriftDirection
Direction of the electron drift in the TPC.
Definition: TPCGeo.h:635
double WirePitch(unsigned plane=0) const
Returns the center of the TPC volume in world coordinates [cm].
Definition: TPCGeo.cxx:426
The data type to uniquely identify a TPC.
Definition: geo_types.h:195
unsigned int MaxWires() const
Returns the largest number of wires among the planes in this TPC.
Definition: TPCGeo.cxx:332
double Plane0Pitch(unsigned int p) const
Returns the center of the TPC volume in world coordinates [cm].
Definition: TPCGeo.cxx:354
void SortPlanes(std::vector< geo::PlaneGeo > &) const
Sorts (in place) the specified PlaneGeo objects by drift distance.
Definition: TPCGeo.cxx:535
double ActiveLength() const
Length (associated with z coordinate) of active TPC volume [cm].
Definition: TPCGeo.h:95
std::vector< std::vector< double > > fPlaneLocation
xyz locations of planes in the TPC.
Definition: TPCGeo.h:637
Encapsulate the geometry of a wire.
auto Xcoord(Vector &v)
Returns an object to manage the coordinate X of the vector v.
double HalfHeight() const
Height is associated with y coordinate [cm].
Definition: TPCGeo.h:103
constexpr Vector Zaxis()
Returns a z axis vector of the specified type.
Definition: geo_vectors.h:226
TGeoVolume * fActiveVolume
Active volume of LAr, called volTPCActive in GDML file.
Definition: TPCGeo.h:633
static const PlaneID_t InvalidID
Special code for an invalid ID.
Definition: geo_types.h:255
A base class aware of world box coordinatesAn object describing a simple shape can inherit from this ...
Definition: BoxBoundedGeo.h:35
Drift towards positive X values.
Definition: geo_types.h:108
Encapsulate the construction of a single detector plane.
double DistanceFromPlane(geo::Point_t const &point) const
Returns the distance of the specified point from the wire plane.
Definition: PlaneGeo.h:557
PlaneGeo const * PlanePtr(unsigned int iplane) const
Returns the plane number iplane from this TPC.
Definition: TPCGeo.h:203
void SetBoundaries(Coord_t x_min, Coord_t x_max, Coord_t y_min, Coord_t y_max, Coord_t z_min, Coord_t z_max)
Sets the boundaries in world coordinates as specified.
short int DetectDriftDirection() const
Returns the expected drift direction based on geometry.
Definition: TPCGeo.cxx:181
geo::Vector_t fLengthDir
Direction length refers to.
Definition: TPCGeo.h:649
LocalTransformation_t fTrans
TPC-to-world transformation.
Definition: TPCGeo.h:630
geo::WireGeo::GeoNodePath_t GeoNodePath_t
Definition: TPCGeo.h:44
void ResetDriftDirection()
Recomputes the drift direction; needs planes to have been initialised.
Definition: TPCGeo.cxx:432
Drift towards negative values.
Definition: geo_types.h:107
#define LOG_DEBUG(id)
geo::Vector_t fDriftDir
Direction electrons drift along.
Definition: TPCGeo.h:650
std::set< geo::View_t > Views() const
Returns a set of all views covered in this TPC.
Definition: TPCGeo.cxx:343
PlaneGeo const & Plane(geo::View_t view) const
Return the plane in the tpc with View_t view.
Definition: TPCGeo.cxx:298
Vector normalize(Vector const &v)
Returns a vector parallel to v and with norm 1.
Float_t e
Definition: plot.C:34
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:187
Namespace collecting geometry-related classes utilities.
geo::Point_t fActiveCenter
Center of the active volume, in world coordinates [cm].
Definition: TPCGeo.h:638
const double * PlaneLocation(unsigned int p) const
Returns the coordinates of the center of the specified plane [cm].
Definition: TPCGeo.cxx:412
double fHalfHeight
Half height of total volume.
Definition: TPCGeo.h:644
double HalfWidth() const
Width is associated with x coordinate [cm].
Definition: TPCGeo.h:99
constexpr Point origin()
Returns a origin position with a point of the specified type.
Definition: geo_vectors.h:230
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
double WirePitch() const
Return the wire pitch (in centimeters). It is assumed constant.
Definition: PlaneGeo.h:367
Encapsulate the construction of a single detector plane.
double fActiveHalfHeight
Half height of active volume.
Definition: TPCGeo.h:641
Point GetCenter() const
Returns the center of the TPC volume in world coordinates [cm].
Definition: TPCGeo.h:693