LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
geo::BoxBoundedGeo Class Reference

A base class aware of world box coordinatesAn object describing a simple shape can inherit from this one to gain access to a common boundary checking interface. More...

#include "BoxBoundedGeo.h"

Inheritance diagram for geo::BoxBoundedGeo:
geo::CryostatGeo geo::TPCGeo

Public Types

using Coords_t = geo::Point_t
 Type of the coordinate triplet. More...
 
using Coord_t = Coords_t::Scalar
 Type of the coordinate. More...
 

Public Member Functions

 BoxBoundedGeo ()=default
 Default constructor: sets an empty volume. More...
 
 BoxBoundedGeo (Coord_t x_min, Coord_t x_max, Coord_t y_min, Coord_t y_max, Coord_t z_min, Coord_t z_max)
 Constructor: sets the boundaries in world coordinates as specified. More...
 
 BoxBoundedGeo (Coords_t lower, Coords_t upper)
 Constructor: sets the boundaries in world coordinates as specified. More...
 
Dimension queries
double MinX () const
 Returns the world x coordinate of the start of the box. More...
 
double MaxX () const
 Returns the world x coordinate of the end of the box. More...
 
double CenterX () const
 Returns the world x coordinate of the center of the box. More...
 
double SizeX () const
 Returns the full size in the X dimension. More...
 
double HalfSizeX () const
 Returns the size from the center to the border on X dimension. More...
 
double MinY () const
 Returns the world y coordinate of the start of the box. More...
 
double MaxY () const
 Returns the world y coordinate of the end of the box. More...
 
double CenterY () const
 Returns the world y coordinate of the center of the box. More...
 
double SizeY () const
 Returns the full size in the Y dimension. More...
 
double HalfSizeY () const
 Returns the size from the center to the border on Y dimension. More...
 
double MinZ () const
 Returns the world z coordinate of the start of the box. More...
 
double MaxZ () const
 Returns the world z coordinate of the end of the box. More...
 
double CenterZ () const
 Returns the world z coordinate of the center of the box. More...
 
double SizeZ () const
 Returns the full size in the Z dimension. More...
 
double HalfSizeZ () const
 Returns the size from the center to the border on Z dimension. More...
 
geo::Point_t Min () const
 Returns the corner point with the smallest coordinates. More...
 
geo::Point_t Max () const
 Returns the corner point with the largest coordinates. More...
 
geo::Point_t Center () const
 Returns the center point of the box. More...
 
Containment in the full volume
bool ContainsX (double x, double const wiggle=1) const
 Returns whether this TPC contains the specified world x coordinate. More...
 
bool ContainsY (double y, double const wiggle=1) const
 Returns whether this TPC contains the specified world y coordinate. More...
 
bool ContainsZ (double z, double const wiggle=1) const
 Returns whether this TPC contains the specified world z coordinate. More...
 
bool ContainsYZ (double y, double z, double const wiggle=1) const
 Returns if TPC contains the specified world y and z coordinates. More...
 
bool ContainsPosition (geo::Point_t const &point, double wiggle=1.0) const
 Returns whether this volume contains the specified point. More...
 
bool ContainsPosition (TVector3 const &point, double wiggle=1.0) const
 
bool ContainsPosition (double const *point, double wiggle=1.0) const
 
Containment in a fiducial volume
bool InFiducialX (double x, double neg_margin, double pos_margin) const
 Returns whether TPC fiducial volume contains world x coordinate. More...
 
bool InFiducialX (double x, double margin) const
 Returns whether TPC fiducial volume contains world x coordinate. More...
 
bool InFiducialY (double y, double neg_margin, double pos_margin) const
 Returns whether TPC fiducial volume contains world y coordinate. More...
 
bool InFiducialY (double y, double margin) const
 Returns whether TPC fiducial volume contains world y coordinate. More...
 
bool InFiducialZ (double z, double neg_margin, double pos_margin) const
 Returns whether TPC fiducial volume contains world z coordinate. More...
 
bool InFiducialZ (double z, double margin) const
 Returns whether TPC fiducial volume contains world z coordinate. More...
 
Setting dimensions
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. More...
 
void SetBoundaries (Coords_t lower, Coords_t upper)
 Sets the boundaries in world coordinates as specified. More...
 
void ExtendToInclude (Coord_t x, Coord_t y, Coord_t z)
 Extends the current box to also include the specified point. More...
 
void ExtendToInclude (geo::Point_t const &point)
 Extends the current box to also include the specified point. More...
 
void ExtendToInclude (BoxBoundedGeo const &box)
 Extends the current box to also include the specified one. More...
 
std::vector< TVector3 > GetIntersections (TVector3 const &TrajectoryStart, TVector3 const &TrajectoryDirect) const
 Calculates the entry and exit points of a trajectory on the box surface. More...
 
std::vector< geo::Point_tGetIntersections (geo::Point_t const &TrajectoryStart, geo::Vector_t const &TrajectoryDirect) const
 Calculates the entry and exit points of a trajectory on the box surface. More...
 

Static Public Member Functions

static bool CoordinateContained (double c, double min, double max, double wiggle=1.)
 Returns whether the specified coordinate is in a range. More...
 
static bool CoordinateContained (double c, double const *range, double wiggle=1.)
 Returns whether the specified coordinate is in a range. More...
 
static void set_min (Coord_t &var, Coord_t value)
 Sets var to value if value is smaller than the current var value. More...
 
static void set_max (Coord_t &var, Coord_t value)
 Sets var to value if value is larger than the current var value. More...
 
static void set_min (Coords_t &var, geo::Point_t const &value)
 Sets each coordinate of var to the one in value if the latter is smaller. More...
 
static void set_max (Coords_t &var, geo::Point_t const &value)
 Sets each coordinate of var to the one in value if the latter is larger. More...
 

Private Member Functions

void SortCoordinates ()
 Makes sure each coordinate of the minimum point is smaller than maximum. More...
 

Private Attributes

Coords_t c_min
 minimum coordinates (x, y, z) More...
 
Coords_t c_max
 maximum coordinates (x, y, z) More...
 

Detailed Description

A base class aware of world box coordinates

An object describing a simple shape can inherit from this one to gain access to a common boundary checking interface.

The boundary is a simple box with axes parallel to the coordinate system.

Definition at line 35 of file BoxBoundedGeo.h.

Member Typedef Documentation

using geo::BoxBoundedGeo::Coord_t = Coords_t::Scalar

Type of the coordinate.

Definition at line 38 of file BoxBoundedGeo.h.

Type of the coordinate triplet.

Definition at line 37 of file BoxBoundedGeo.h.

Constructor & Destructor Documentation

geo::BoxBoundedGeo::BoxBoundedGeo ( )
default

Default constructor: sets an empty volume.

See also
SetBoundaries

We allow for a default construction since the derived object might need some specific action before being aware of its boundaries. In that case, SetBoundaries() will set the boundaries.

geo::BoxBoundedGeo::BoxBoundedGeo ( Coord_t  x_min,
Coord_t  x_max,
Coord_t  y_min,
Coord_t  y_max,
Coord_t  z_min,
Coord_t  z_max 
)
inline

Constructor: sets the boundaries in world coordinates as specified.

Parameters
x_minlower x coordinate
x_maxupper x coordinate
y_minlower y coordinate
y_maxupper y coordinate
z_minlower z coordinate
z_maxupper z coordinate
See also
SetBoundaries()

Note that is assumed that each minimum is larger than its maximum, and no check is performed.

Definition at line 64 of file BoxBoundedGeo.h.

References c_max, SortCoordinates(), x_max, and x_min.

68  :
69  c_min{ x_min, y_min, z_min }, c_max{ x_max, y_max, z_max }
70  { SortCoordinates(); }
double x_min
Definition: berger.C:15
Coords_t c_max
maximum coordinates (x, y, z)
double x_max
Definition: berger.C:16
Coords_t c_min
minimum coordinates (x, y, z)
void SortCoordinates()
Makes sure each coordinate of the minimum point is smaller than maximum.
geo::BoxBoundedGeo::BoxBoundedGeo ( Coords_t  lower,
Coords_t  upper 
)
inline

Constructor: sets the boundaries in world coordinates as specified.

Parameters
lowerlower coordinates (x, y, z)
upperupper coordinates (x, y, z)
See also
SetBoundaries()

Note that is assumed that each minimum is larger than its maximum, and no check is performed.

Definition at line 81 of file BoxBoundedGeo.h.

References SortCoordinates().

81  :
82  c_min(lower), c_max(upper)
83  { SortCoordinates(); }
Coords_t c_max
maximum coordinates (x, y, z)
Coords_t c_min
minimum coordinates (x, y, z)
void SortCoordinates()
Makes sure each coordinate of the minimum point is smaller than maximum.

Member Function Documentation

geo::Point_t geo::BoxBoundedGeo::Center ( ) const
inline

Returns the center point of the box.

Definition at line 141 of file BoxBoundedGeo.h.

References CenterX(), CenterY(), and CenterZ().

Referenced by geo::buildDriftVolumes(), geo::CryostatGeo::GetCenter(), and geo::PlaneGeo::UpdatePlaneNormal().

142  { return { CenterX(), CenterY(), CenterZ() }; }
double CenterX() const
Returns the world x coordinate of the center of the box.
Definition: BoxBoundedGeo.h:96
double CenterZ() const
Returns the world z coordinate of the center of the box.
double CenterY() const
Returns the world y coordinate of the center of the box.
double geo::BoxBoundedGeo::CenterX ( ) const
inline

Returns the world x coordinate of the center of the box.

Definition at line 96 of file BoxBoundedGeo.h.

References MaxX(), and MinX().

Referenced by Center().

96 { return (MinX() + MaxX()) / 2.; }
double MinX() const
Returns the world x coordinate of the start of the box.
Definition: BoxBoundedGeo.h:90
double MaxX() const
Returns the world x coordinate of the end of the box.
Definition: BoxBoundedGeo.h:93
double geo::BoxBoundedGeo::CenterY ( ) const
inline

Returns the world y coordinate of the center of the box.

Definition at line 111 of file BoxBoundedGeo.h.

References MaxY(), and MinY().

Referenced by Center().

111 { return (MinY() + MaxY()) / 2.; }
double MaxY() const
Returns the world y coordinate of the end of the box.
double MinY() const
Returns the world y coordinate of the start of the box.
double geo::BoxBoundedGeo::CenterZ ( ) const
inline

Returns the world z coordinate of the center of the box.

Definition at line 126 of file BoxBoundedGeo.h.

References MaxZ(), and MinZ().

Referenced by Center().

126 { return (MinZ() + MaxZ()) / 2.; }
double MinZ() const
Returns the world z coordinate of the start of the box.
double MaxZ() const
Returns the world z coordinate of the end of the box.
bool geo::BoxBoundedGeo::ContainsPosition ( geo::Point_t const &  point,
double  wiggle = 1.0 
) const
inline

Returns whether this volume contains the specified point.

Parameters
pointthe point [cm]
wiggleexpansion factor for the range
Returns
whether the specified coordinate is in this volume

If the wiggle is larger than 1, each size of the volume is expanded by the wiggle factor. If the wiggle is less than 1, each size is shrinked.

Definition at line 207 of file BoxBoundedGeo.h.

References ContainsX(), and ContainsYZ().

Referenced by ContainsPosition(), ContainsYZ(), and geo::GeometryCore::MaterialName().

208  {
209  return ContainsX(point.X(), wiggle)
210  && ContainsYZ(point.Y(), point.Z(), wiggle);
211  } // ContainsPosition()
bool ContainsX(double x, double const wiggle=1) const
Returns whether this TPC contains the specified world x coordinate.
bool ContainsYZ(double y, double z, double const wiggle=1) const
Returns if TPC contains the specified world y and z coordinates.
bool geo::BoxBoundedGeo::ContainsPosition ( TVector3 const &  point,
double  wiggle = 1.0 
) const
See also
ContainsPosition(geo::Point_t const&, double) const.

Definition at line 18 of file BoxBoundedGeo.cxx.

References ContainsPosition(), and geo::vect::toPoint().

19  { return ContainsPosition(geo::vect::toPoint(point), wiggle); }
::geo::Point_t toPoint(Point const &p)
Convert the specified point into a geo::Point_t.
bool ContainsPosition(geo::Point_t const &point, double wiggle=1.0) const
Returns whether this volume contains the specified point.
bool geo::BoxBoundedGeo::ContainsPosition ( double const *  point,
double  wiggle = 1.0 
) const
See also
ContainsPosition(geo::Point_t const&, double) const.

Definition at line 23 of file BoxBoundedGeo.cxx.

References geo::vect::makePointFromCoords().

24  { return ContainsPosition(geo::vect::makePointFromCoords(point), wiggle); }
GENVECTOR_CONSTEXPR::geo::Point_t makePointFromCoords(Coords &&coords)
Creates a geo::Point_t from its coordinates (see makeFromCoords()).
bool ContainsPosition(geo::Point_t const &point, double wiggle=1.0) const
Returns whether this volume contains the specified point.
bool geo::BoxBoundedGeo::ContainsX ( double  x,
double const  wiggle = 1 
) const
inline

Returns whether this TPC contains the specified world x coordinate.

Parameters
xthe absolute ("world") coordinate x
wiggleexpansion factor for the range (see ContainsPosition())
Returns
whether the specified coordinate is in this TPC
See also
ContainsPosition()

Note that x is by definition the drift direction, and a reconstructed x typically depends on an assumption respect to the event time.

Definition at line 161 of file BoxBoundedGeo.h.

References CoordinateContained(), MaxX(), and MinX().

Referenced by ContainsPosition().

162  { return CoordinateContained(x, MinX(), MaxX(), wiggle); }
Float_t x
Definition: compare.C:6
double MinX() const
Returns the world x coordinate of the start of the box.
Definition: BoxBoundedGeo.h:90
double MaxX() const
Returns the world x coordinate of the end of the box.
Definition: BoxBoundedGeo.h:93
static bool CoordinateContained(double c, double min, double max, double wiggle=1.)
Returns whether the specified coordinate is in a range.
bool geo::BoxBoundedGeo::ContainsY ( double  y,
double const  wiggle = 1 
) const
inline

Returns whether this TPC contains the specified world y coordinate.

Parameters
ythe absolute ("world") coordinate y
wiggleexpansion factor for the range (see ContainsPosition())
Returns
whether the specified coordinate is in this TPC
See also
ContainsPosition()

Definition at line 171 of file BoxBoundedGeo.h.

References CoordinateContained(), MaxY(), and MinY().

Referenced by ContainsYZ().

172  { return CoordinateContained(y, MinY(), MaxY(), wiggle); }
Float_t y
Definition: compare.C:6
double MaxY() const
Returns the world y coordinate of the end of the box.
static bool CoordinateContained(double c, double min, double max, double wiggle=1.)
Returns whether the specified coordinate is in a range.
double MinY() const
Returns the world y coordinate of the start of the box.
bool geo::BoxBoundedGeo::ContainsYZ ( double  y,
double  z,
double const  wiggle = 1 
) const
inline

Returns if TPC contains the specified world y and z coordinates.

Parameters
ythe absolute ("world") coordinate y
zthe absolute ("world") coordinate z
wiggleexpansion factor for the range (see ContainsPosition())
Returns
whether the specified coordinate is in this TPC
See also
ContainsPosition()

Definition at line 192 of file BoxBoundedGeo.h.

References ContainsPosition(), ContainsY(), and ContainsZ().

Referenced by ContainsPosition().

193  { return ContainsY(y, wiggle) && ContainsZ(z, wiggle); }
Float_t y
Definition: compare.C:6
Double_t z
Definition: plot.C:279
bool ContainsY(double y, double const wiggle=1) const
Returns whether this TPC contains the specified world y coordinate.
bool ContainsZ(double z, double const wiggle=1) const
Returns whether this TPC contains the specified world z coordinate.
bool geo::BoxBoundedGeo::ContainsZ ( double  z,
double const  wiggle = 1 
) const
inline

Returns whether this TPC contains the specified world z coordinate.

Parameters
zthe absolute ("world") coordinate z
wiggleexpansion factor for the range (see ContainsPosition())
Returns
whether the specified coordinate is in this TPC
See also
ContainsPosition()

Definition at line 181 of file BoxBoundedGeo.h.

References CoordinateContained(), MaxZ(), and MinZ().

Referenced by ContainsYZ().

182  { return CoordinateContained(z, MinZ(), MaxZ(), wiggle); }
Double_t z
Definition: plot.C:279
double MinZ() const
Returns the world z coordinate of the start of the box.
double MaxZ() const
Returns the world z coordinate of the end of the box.
static bool CoordinateContained(double c, double min, double max, double wiggle=1.)
Returns whether the specified coordinate is in a range.
static bool geo::BoxBoundedGeo::CoordinateContained ( double  c,
double  min,
double  max,
double  wiggle = 1. 
)
inlinestatic

Returns whether the specified coordinate is in a range.

Parameters
cthe coordinate
minlower boundary of the range
maxupper boundary of the range
wiggleexpansion factor for the range
Returns
whether the specified coordinate is in a range

If the wiggle is larger than 1, the range is expanded by the wiggle factor. If the wiggle is less than 1, the range is shrinked.

Definition at line 333 of file BoxBoundedGeo.h.

Referenced by ContainsX(), ContainsY(), ContainsZ(), CoordinateContained(), InFiducialX(), InFiducialY(), and InFiducialZ().

334  {
335  return (c >= (min > 0? min / wiggle: min * wiggle))
336  && (c <= (max < 0? max / wiggle: max * wiggle));
337  } // CoordinateContained()
Int_t max
Definition: plot.C:27
Int_t min
Definition: plot.C:26
static bool geo::BoxBoundedGeo::CoordinateContained ( double  c,
double const *  range,
double  wiggle = 1. 
)
inlinestatic

Returns whether the specified coordinate is in a range.

Parameters
cthe coordinate
rangepointer to [ lower boundary, upper boundary ] of the range
wiggleexpansion factor for the range
Returns
whether the specified coordinate is in a range
See also
CoordinateContained(double, double, double, double)

If the wiggle is larger than 1, the range is expanded by the wiggle factor. If the wiggle is less than 1, the range is shrinked.

Definition at line 351 of file BoxBoundedGeo.h.

References CoordinateContained().

352  { return CoordinateContained(c, range[0], range[1], wiggle); }
static bool CoordinateContained(double c, double min, double max, double wiggle=1.)
Returns whether the specified coordinate is in a range.
void geo::BoxBoundedGeo::ExtendToInclude ( Coord_t  x,
Coord_t  y,
Coord_t  z 
)
inline

Extends the current box to also include the specified point.

Parameters
xx coordinate of the point to include
yy coordinate of the point to include
zz coordinate of the point to include

Definition at line 392 of file BoxBoundedGeo.h.

Referenced by geo::PlaneGeo::BoundingBox(), evd_tool::StandardDrawer::DetOutline3D(), and lar::example::SpacePointIsolationAlg::fillAlgConfigFromGeometry().

393  { ExtendToInclude(geo::Point_t(x, y, z)); }
Float_t x
Definition: compare.C:6
Float_t y
Definition: compare.C:6
Double_t z
Definition: plot.C:279
void ExtendToInclude(Coord_t x, Coord_t y, Coord_t z)
Extends the current box to also include the specified point.
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
void geo::BoxBoundedGeo::ExtendToInclude ( geo::Point_t const &  point)
inline

Extends the current box to also include the specified point.

Parameters
pointcoordinates of the point to include

Definition at line 399 of file BoxBoundedGeo.h.

References c_max, c_min, set_max(), and set_min().

400  {
401  set_min(c_min, point);
402  set_max(c_max, point);
403  }
static void set_min(Coord_t &var, Coord_t value)
Sets var to value if value is smaller than the current var value.
static void set_max(Coord_t &var, Coord_t value)
Sets var to value if value is larger than the current var value.
Coords_t c_max
maximum coordinates (x, y, z)
Coords_t c_min
minimum coordinates (x, y, z)
void geo::BoxBoundedGeo::ExtendToInclude ( BoxBoundedGeo const &  box)
inline

Extends the current box to also include the specified one.

Parameters
boxthe box to include

It is assumed that the box has its boundaries properly sorted.

Definition at line 411 of file BoxBoundedGeo.h.

References c_max, c_min, GetIntersections(), Max(), Min(), set_max(), and set_min().

412  {
413  set_min(c_min, box.Min());
414  set_max(c_max, box.Max());
415  } // ExtendToInclude()
static void set_min(Coord_t &var, Coord_t value)
Sets var to value if value is smaller than the current var value.
static void set_max(Coord_t &var, Coord_t value)
Sets var to value if value is larger than the current var value.
Coords_t c_max
maximum coordinates (x, y, z)
Coords_t c_min
minimum coordinates (x, y, z)
std::vector< TVector3 > geo::BoxBoundedGeo::GetIntersections ( TVector3 const &  TrajectoryStart,
TVector3 const &  TrajectoryDirect 
) const

Calculates the entry and exit points of a trajectory on the box surface.

Author
Christoph Rudolf von Rohr (crohr.nosp@m.@fna.nosp@m.l.gov)
Date
July 27th, 2015
Parameters
TrajectoryStartposition of the trajectory source
TrajectoryDirectdirection vector of the trajectory

This member is public since it just gives an output and does not change any member. The algorithm works only for a box shaped active volume with facing walls parallel to axis. If the return std::vector is empty the trajectory does not intersect with the box. Normally the return value should have one (if the trajectory originates in the box) or two (else) entries. If the return value has two entries the first represents the entry point and the second the exit point

Definition at line 118 of file BoxBoundedGeo.cxx.

Referenced by ExtendToInclude(), and GetIntersections().

119  {
120  std::vector<TVector3> intersections;
121  for (auto const& point: GetIntersections(geo::Point_t(TrajectoryStart), geo::Vector_t(TrajectoryDirect)))
122  intersections.emplace_back(point.X(), point.Y(), point.Z());
123  return intersections;
124  } // GetIntersections(TVector3)
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Vector_t
Type for representation of momenta in 3D space.
Definition: geo_vectors.h:167
std::vector< TVector3 > GetIntersections(TVector3 const &TrajectoryStart, TVector3 const &TrajectoryDirect) const
Calculates the entry and exit points of a trajectory on the box surface.
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
std::vector< geo::Point_t > geo::BoxBoundedGeo::GetIntersections ( geo::Point_t const &  TrajectoryStart,
geo::Vector_t const &  TrajectoryDirect 
) const

Calculates the entry and exit points of a trajectory on the box surface.

Author
Christoph Rudolf von Rohr (crohr.nosp@m.@fna.nosp@m.l.gov)
Date
July 27th, 2015
Parameters
TrajectoryStartposition of the trajectory source
TrajectoryDirectdirection vector of the trajectory

This member is public since it just gives an output and does not change any member. The algorithm works only for a box shaped active volume with facing walls parallel to axis. If the return std::vector is empty the trajectory does not intersect with the box. Normally the return value should have one (if the trajectory originates in the box) or two (else) entries. If the return value has two entries the first represents the entry point and the second the exit point

Definition at line 27 of file BoxBoundedGeo.cxx.

References geo::vect::bindCoord(), c_max, c_min, geo::vect::coord(), GetIntersections(), Max(), and Min().

30  {
31 
32  std::vector<geo::Point_t> IntersectionPoints;
33  std::vector<double> LineParameters;
34 
35  // Generate normal vectors and offsets for every plane of the box
36  // All normal vectors are headed outwards
37  // BUG the double brace syntax is required to work around clang bug 21629
38  // (https://bugs.llvm.org/show_bug.cgi?id=21629)
39  static std::array<geo::Vector_t, 6U> const NormalVectors = {{
40  -geo::Xaxis<geo::Vector_t>(), geo::Xaxis<geo::Vector_t>(), // anode, cathode,
41  -geo::Yaxis<geo::Vector_t>(), geo::Yaxis<geo::Vector_t>(), // bottom, top,
42  -geo::Zaxis<geo::Vector_t>(), geo::Zaxis<geo::Vector_t>() // upstream, downstream
43  }};
44  // BUG the double brace syntax is required to work around clang bug 21629
45  // (https://bugs.llvm.org/show_bug.cgi?id=21629)
46  std::array<geo::Point_t, 6U> const NormalVectorOffsets = {{
47  geo::Point_t{ Min().X(), Min().Y(), Min().Z() }, // Anode offset
48  geo::Point_t{ Max().X(), Min().Y(), Min().Z() }, // Cathode
49  geo::Point_t{ Min().X(), Min().Y(), Min().Z() }, // Bottom
50  geo::Point_t{ Min().X(), Max().Y(), Min().Z() }, // Top
51  geo::Point_t{ Min().X(), Min().Y(), Min().Z() }, // upstream
52  geo::Point_t{ Min().X(), Min().Y(), Max().Z() } // downstream
53  }};
54 
55  // Loop over all surfaces of the box
56  for(unsigned int face_no = 0; face_no < NormalVectors.size(); face_no++)
57  {
58  // Check if trajectory and surface are not parallel
59  if(NormalVectors[face_no].Dot(TrajectoryDirect))
60  {
61  // Calculate the line parameter for the intersection points
62  LineParameters.push_back( NormalVectors[face_no].Dot(NormalVectorOffsets.at(face_no) - TrajectoryStart)
63  / NormalVectors[face_no].Dot(TrajectoryDirect) );
64  }
65  else continue;
66 
67  // Calculate intersection point using the line parameter
68  IntersectionPoints.push_back( TrajectoryStart + LineParameters.back()*TrajectoryDirect );
69 
70  // Coordinate which should be ignored when checking for limits added by Christoph Rudolf von Rohr 05/21/2016
71  unsigned int NoCheckCoord;
72 
73  // Calculate NoCheckCoord out of the face_no
74  if(face_no % 2)
75  {
76  // Convert odd face number to coordinate
77  NoCheckCoord = (face_no - 1)/2;
78  }
79  else
80  {
81  // Convert even face number to coordinate
82  NoCheckCoord = face_no/2;
83  }
84 
85  // Loop over all three space coordinates
86  unsigned int coord = 0;
87  for(auto extractCoord: geo::vect::coordReaders<geo::Point_t>())
88  {
89  auto const lastPointCoord = geo::vect::bindCoord(IntersectionPoints.back(), extractCoord);
90  auto const minCoord = geo::vect::bindCoord(c_min, extractCoord);
91  auto const maxCoord = geo::vect::bindCoord(c_max, extractCoord);
92 
93  // Changed by Christoph Rudolf von Rohr 05/21/2016
94  // Then check if point is not within the surface limits at this coordinate, without looking
95  // at the plane normal vector coordinate. We can assume, that our algorithm already found this coordinate correctily.
96  // In rare cases, looking for boundaries in this coordinate causes unexpected behavior due to floating point inaccuracies.
97  if( coord++ != NoCheckCoord && ((lastPointCoord() > maxCoord()) || (lastPointCoord() < minCoord)) )
98  {
99  // if off limits, get rid of the useless data and break the coordinate loop
100  LineParameters.pop_back();
101  IntersectionPoints.pop_back();
102  break;
103  }
104  }// coordinate loop
105  }// Surcaces loop
106 
107  // sort points according to their parameter value (first is entry, second is exit)
108  if(LineParameters.size() == 2 && LineParameters.front() > LineParameters.back())
109  {
110  std::swap(IntersectionPoints.front(),IntersectionPoints.back());
111  }
112 
113  return IntersectionPoints;
114  } // GetIntersections()
auto coord(Vector &v, unsigned int n) noexcept
Returns an object to manage the coordinate n of a vector.
constexpr auto bindCoord(Vector const &v, CoordReader_t< Vector > helper)
Binds the specified constant vector to the coordinate reader.
geo::Point_t Min() const
Returns the corner point with the smallest coordinates.
Coords_t c_max
maximum coordinates (x, y, z)
void swap(art::HLTGlobalStatus &lhs, art::HLTGlobalStatus &rhs)
Coords_t c_min
minimum coordinates (x, y, z)
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
geo::Point_t Max() const
Returns the corner point with the largest coordinates.
double geo::BoxBoundedGeo::HalfSizeX ( ) const
inline

Returns the size from the center to the border on X dimension.

Definition at line 102 of file BoxBoundedGeo.h.

References SizeX().

102 { return SizeX() / 2.0; }
double SizeX() const
Returns the full size in the X dimension.
Definition: BoxBoundedGeo.h:99
double geo::BoxBoundedGeo::HalfSizeY ( ) const
inline

Returns the size from the center to the border on Y dimension.

Definition at line 117 of file BoxBoundedGeo.h.

References SizeY().

117 { return SizeY() / 2.0; }
double SizeY() const
Returns the full size in the Y dimension.
double geo::BoxBoundedGeo::HalfSizeZ ( ) const
inline

Returns the size from the center to the border on Z dimension.

Definition at line 132 of file BoxBoundedGeo.h.

References SizeZ().

132 { return SizeZ() / 2.0; }
double SizeZ() const
Returns the full size in the Z dimension.
bool geo::BoxBoundedGeo::InFiducialX ( double  x,
double  neg_margin,
double  pos_margin 
) const
inline

Returns whether TPC fiducial volume contains world x coordinate.

Parameters
xthe absolute ("world") coordinate x
neg_marginhow far within the TPC the fiducial region starts
pos_marginhow far before the TPC the fiducial region ends
Returns
whether the specified coordinate is in this TPC fiducial volume

The fiducial volume is defined by the specified margins, that denote how much of the coordinate is not fiducial. For example,

bool bWithin = tpc->InFiducialX(x, 5., 8.);

on a TPC with x from 0 to 100 (cm, presumably) will return true for all x between 5 and 92. If positive margin is not specified, it's assumed to be the same as the negative one. Specifying a negative mergin effectively extends the TPC volume. Note that x is by definition the drift direction, and a reconstructed x typically depends on an assumption respect to the event time.

Definition at line 243 of file BoxBoundedGeo.h.

References CoordinateContained(), MaxX(), and MinX().

Referenced by InFiducialX().

244  {
245  return CoordinateContained(x, MinX() + neg_margin, MaxX() - pos_margin);
246  }
Float_t x
Definition: compare.C:6
double MinX() const
Returns the world x coordinate of the start of the box.
Definition: BoxBoundedGeo.h:90
double MaxX() const
Returns the world x coordinate of the end of the box.
Definition: BoxBoundedGeo.h:93
static bool CoordinateContained(double c, double min, double max, double wiggle=1.)
Returns whether the specified coordinate is in a range.
bool geo::BoxBoundedGeo::InFiducialX ( double  x,
double  margin 
) const
inline

Returns whether TPC fiducial volume contains world x coordinate.

See also
InFiducialX(double, double, double) const

Margins are symmetric.

Definition at line 253 of file BoxBoundedGeo.h.

References InFiducialX().

254  { return InFiducialX(x, margin, margin); }
Float_t x
Definition: compare.C:6
bool InFiducialX(double x, double neg_margin, double pos_margin) const
Returns whether TPC fiducial volume contains world x coordinate.
bool geo::BoxBoundedGeo::InFiducialY ( double  y,
double  neg_margin,
double  pos_margin 
) const
inline

Returns whether TPC fiducial volume contains world y coordinate.

Parameters
ythe absolute ("world") coordinate y
neg_marginhow far within the TPC the fiducial region starts
pos_marginhow far before the TPC the fiducial region ends
Returns
whether the specified coordinate is in this TPC fiducial volume

The fiducial volume is defined by the specified margins, that denote how much of the coordinate is not fiducial. For example,

bool bWithin = tpc->InFiducialY(y, 5., 8.);

on a TPC with y from 0 to 100 (cm, presumably) will return true for all y between 5 and 92. If positive margin is not specified, it's assumed to be the same as the negative one. Specifying a negative mergin effectively extends the TPC volume.

Definition at line 274 of file BoxBoundedGeo.h.

References CoordinateContained(), MaxY(), and MinY().

Referenced by InFiducialY().

275  {
276  return CoordinateContained(y, MinY() + neg_margin, MaxY() - pos_margin);
277  }
Float_t y
Definition: compare.C:6
double MaxY() const
Returns the world y coordinate of the end of the box.
static bool CoordinateContained(double c, double min, double max, double wiggle=1.)
Returns whether the specified coordinate is in a range.
double MinY() const
Returns the world y coordinate of the start of the box.
bool geo::BoxBoundedGeo::InFiducialY ( double  y,
double  margin 
) const
inline

Returns whether TPC fiducial volume contains world y coordinate.

See also
InFiducialY(double, double, double) const

Margins are symmetric.

Definition at line 284 of file BoxBoundedGeo.h.

References InFiducialY().

285  { return InFiducialY(y, margin, margin); }
Float_t y
Definition: compare.C:6
bool InFiducialY(double y, double neg_margin, double pos_margin) const
Returns whether TPC fiducial volume contains world y coordinate.
bool geo::BoxBoundedGeo::InFiducialZ ( double  z,
double  neg_margin,
double  pos_margin 
) const
inline

Returns whether TPC fiducial volume contains world z coordinate.

Parameters
zthe absolute ("world") coordinate y
neg_marginhow far within the TPC the fiducial region starts
pos_marginhow far before the TPC the fiducial region ends
Returns
whether the specified coordinate is in this TPC fiducial volume

The fiducial volume is defined by the specified margins, that denote how much of the coordinate is not fiducial. For example,

bool bWithin = tpc->InFiducialZ(z, 5., 8.);

on a TPC with z from 0 to 100 (cm, presumably) will return true for all z between 5 and 92. If positive margin is not specified, it's assumed to be the same as the negative one. Specifying a negative mergin effectively extends the TPC volume.

Definition at line 305 of file BoxBoundedGeo.h.

References CoordinateContained(), MaxZ(), and MinZ().

Referenced by InFiducialZ().

306  {
307  return CoordinateContained(z, MinZ() + neg_margin, MaxZ() - pos_margin);
308  }
Double_t z
Definition: plot.C:279
double MinZ() const
Returns the world z coordinate of the start of the box.
double MaxZ() const
Returns the world z coordinate of the end of the box.
static bool CoordinateContained(double c, double min, double max, double wiggle=1.)
Returns whether the specified coordinate is in a range.
bool geo::BoxBoundedGeo::InFiducialZ ( double  z,
double  margin 
) const
inline

Returns whether TPC fiducial volume contains world z coordinate.

See also
InFiducialZ(double, double, double) const

Margins are symmetric.

Definition at line 315 of file BoxBoundedGeo.h.

References CoordinateContained(), and InFiducialZ().

316  { return InFiducialZ(z, margin, margin); }
Double_t z
Definition: plot.C:279
bool InFiducialZ(double z, double neg_margin, double pos_margin) const
Returns whether TPC fiducial volume contains world z coordinate.
geo::Point_t geo::BoxBoundedGeo::Max ( ) const
inline

Returns the corner point with the largest coordinates.

Definition at line 138 of file BoxBoundedGeo.h.

References c_max.

Referenced by ExtendToInclude(), GetIntersections(), geo::GeometryCore::MaterialName(), geo::CryostatGeo::PrintCryostatInfo(), and geo::TPCGeo::PrintTPCInfo().

138 { return c_max; }
Coords_t c_max
maximum coordinates (x, y, z)
geo::Point_t geo::BoxBoundedGeo::Min ( ) const
inline

Returns the corner point with the smallest coordinates.

Definition at line 135 of file BoxBoundedGeo.h.

References c_min.

Referenced by ExtendToInclude(), GetIntersections(), geo::GeometryCore::MaterialName(), geo::CryostatGeo::PrintCryostatInfo(), and geo::TPCGeo::PrintTPCInfo().

135 { return c_min; }
Coords_t c_min
minimum coordinates (x, y, z)
static void geo::BoxBoundedGeo::set_max ( Coord_t var,
Coord_t  value 
)
inlinestatic

Sets var to value if value is larger than the current var value.

Definition at line 446 of file BoxBoundedGeo.h.

References fhicl::detail::atom::value().

Referenced by ExtendToInclude().

447  { if (value > var) var = value; }
std::string value(boost::any const &)
static void geo::BoxBoundedGeo::set_max ( Coords_t var,
geo::Point_t const &  value 
)
inlinestatic

Sets each coordinate of var to the one in value if the latter is larger.

Definition at line 458 of file BoxBoundedGeo.h.

459  {
460  if (value.X() > var.X()) var.SetX(value.X());
461  if (value.Y() > var.Y()) var.SetY(value.Y());
462  if (value.Z() > var.Z()) var.SetZ(value.Z());
463  }
std::string value(boost::any const &)
static void geo::BoxBoundedGeo::set_min ( Coord_t var,
Coord_t  value 
)
inlinestatic

Sets var to value if value is smaller than the current var value.

Definition at line 442 of file BoxBoundedGeo.h.

References fhicl::detail::atom::value().

Referenced by ExtendToInclude().

443  { if (value < var) var = value; }
std::string value(boost::any const &)
static void geo::BoxBoundedGeo::set_min ( Coords_t var,
geo::Point_t const &  value 
)
inlinestatic

Sets each coordinate of var to the one in value if the latter is smaller.

Definition at line 450 of file BoxBoundedGeo.h.

451  {
452  if (value.X() < var.X()) var.SetX(value.X());
453  if (value.Y() < var.Y()) var.SetY(value.Y());
454  if (value.Z() < var.Z()) var.SetZ(value.Z());
455  }
std::string value(boost::any const &)
void geo::BoxBoundedGeo::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 
)
inline

Sets the boundaries in world coordinates as specified.

Parameters
x_minlower x coordinate
x_maxupper x coordinate
y_minlower y coordinate
y_maxupper y coordinate
z_minlower z coordinate
z_maxupper z coordinate

Definition at line 367 of file BoxBoundedGeo.h.

References c_max, c_min, and SortCoordinates().

Referenced by geo::PlaneGeo::BoundingBox(), geo::CryostatGeo::InitCryoBoundaries(), and geo::TPCGeo::InitTPCBoundaries().

372  {
373  c_min.SetXYZ(x_min, y_min, z_min);
374  c_max.SetXYZ(x_max, y_max, z_max);
375  SortCoordinates();
376  }
double x_min
Definition: berger.C:15
Coords_t c_max
maximum coordinates (x, y, z)
double x_max
Definition: berger.C:16
Coords_t c_min
minimum coordinates (x, y, z)
void SortCoordinates()
Makes sure each coordinate of the minimum point is smaller than maximum.
void geo::BoxBoundedGeo::SetBoundaries ( Coords_t  lower,
Coords_t  upper 
)
inline

Sets the boundaries in world coordinates as specified.

Parameters
lowerlower coordinates (x, y, z)
upperupper coordinates (x, y, z)

Definition at line 383 of file BoxBoundedGeo.h.

References c_max, c_min, and SortCoordinates().

384  { c_min = lower; c_max = upper; SortCoordinates(); }
Coords_t c_max
maximum coordinates (x, y, z)
Coords_t c_min
minimum coordinates (x, y, z)
void SortCoordinates()
Makes sure each coordinate of the minimum point is smaller than maximum.
double geo::BoxBoundedGeo::SizeX ( ) const
inline

Returns the full size in the X dimension.

Definition at line 99 of file BoxBoundedGeo.h.

References MaxX(), and MinX().

Referenced by HalfSizeX().

99 { return MaxX() - MinX(); }
double MinX() const
Returns the world x coordinate of the start of the box.
Definition: BoxBoundedGeo.h:90
double MaxX() const
Returns the world x coordinate of the end of the box.
Definition: BoxBoundedGeo.h:93
double geo::BoxBoundedGeo::SizeY ( ) const
inline

Returns the full size in the Y dimension.

Definition at line 114 of file BoxBoundedGeo.h.

References MaxY(), and MinY().

Referenced by HalfSizeY().

114 { return MaxY() - MinY(); }
double MaxY() const
Returns the world y coordinate of the end of the box.
double MinY() const
Returns the world y coordinate of the start of the box.
double geo::BoxBoundedGeo::SizeZ ( ) const
inline

Returns the full size in the Z dimension.

Definition at line 129 of file BoxBoundedGeo.h.

References MaxZ(), and MinZ().

Referenced by HalfSizeZ().

129 { return MaxZ() - MinZ(); }
double MinZ() const
Returns the world z coordinate of the start of the box.
double MaxZ() const
Returns the world z coordinate of the end of the box.
void geo::BoxBoundedGeo::SortCoordinates ( )
private

Makes sure each coordinate of the minimum point is smaller than maximum.

Definition at line 127 of file BoxBoundedGeo.cxx.

References geo::vect::bindCoord(), c_max, c_min, max, and min.

Referenced by BoxBoundedGeo(), and SetBoundaries().

127  {
128 
129  for (auto coordMan: geo::vect::coordManagers<geo::Point_t>()) {
130  auto min = geo::vect::bindCoord(c_min, coordMan);
131  auto max = geo::vect::bindCoord(c_max, coordMan);
132  if (min() > max()) {
133  auto temp = min();
134  min = max();
135  max = temp;
136  }
137  } // for
138 
139  } // BoxBoundedGeo::SortCoordinates()
constexpr auto bindCoord(Vector const &v, CoordReader_t< Vector > helper)
Binds the specified constant vector to the coordinate reader.
Int_t max
Definition: plot.C:27
Coords_t c_max
maximum coordinates (x, y, z)
Coords_t c_min
minimum coordinates (x, y, z)
Int_t min
Definition: plot.C:26

Member Data Documentation

Coords_t geo::BoxBoundedGeo::c_max
private

maximum coordinates (x, y, z)

Definition at line 469 of file BoxBoundedGeo.h.

Referenced by BoxBoundedGeo(), ExtendToInclude(), GetIntersections(), Max(), MaxX(), MaxY(), MaxZ(), SetBoundaries(), and SortCoordinates().

Coords_t geo::BoxBoundedGeo::c_min
private

minimum coordinates (x, y, z)

Definition at line 468 of file BoxBoundedGeo.h.

Referenced by ExtendToInclude(), GetIntersections(), Min(), MinX(), MinY(), MinZ(), SetBoundaries(), and SortCoordinates().


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