Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
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...
 
Overlaps
bool OverlapsX (geo::BoxBoundedGeo const &other) const
 Returns if the x coordinates covered by this and other box overlap. More...
 
bool OverlapsY (geo::BoxBoundedGeo const &other) const
 Returns if the y coordinates covered by this and other box overlap. More...
 
bool OverlapsZ (geo::BoxBoundedGeo const &other) const
 Returns if the z coordinates covered by this and other box overlap. More...
 
bool Overlaps (geo::BoxBoundedGeo const &other) const
 Returns if this and other box overlap. 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
 

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 33 of file BoxBoundedGeo.h.

Member Typedef Documentation

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

Type of the coordinate.

Definition at line 36 of file BoxBoundedGeo.h.

Type of the coordinate triplet.

Definition at line 35 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 62 of file BoxBoundedGeo.h.

66  :
67  c_min{ x_min, y_min, z_min }, c_max{ x_max, y_max, z_max }
68  { 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.
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 79 of file BoxBoundedGeo.h.

79  :
80  c_min(lower), c_max(upper)
81  { 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 139 of file BoxBoundedGeo.h.

140  { return { CenterX(), CenterY(), CenterZ() }; }
double CenterX() const
Returns the world x coordinate of the center of the box.
Definition: BoxBoundedGeo.h:94
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 94 of file BoxBoundedGeo.h.

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

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

Definition at line 109 of file BoxBoundedGeo.h.

109 { 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 124 of file BoxBoundedGeo.h.

124 { 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 205 of file BoxBoundedGeo.h.

206  {
207  return ContainsX(point.X(), wiggle)
208  && ContainsYZ(point.Y(), point.Z(), wiggle);
209  } // 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 26 of file BoxBoundedGeo.cxx.

27  { 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 31 of file BoxBoundedGeo.cxx.

32  { 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 159 of file BoxBoundedGeo.h.

160  { return CoordinateContained(x, MinX(), MaxX(), wiggle); }
double MinX() const
Returns the world x coordinate of the start of the box.
Definition: BoxBoundedGeo.h:88
double MaxX() const
Returns the world x coordinate of the end of the box.
Definition: BoxBoundedGeo.h:91
static bool CoordinateContained(double c, double min, double max, double wiggle=1.)
Returns whether the specified coordinate is in a range.
list x
Definition: train.py:276
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 169 of file BoxBoundedGeo.h.

170  { return CoordinateContained(y, MinY(), MaxY(), wiggle); }
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 190 of file BoxBoundedGeo.h.

191  { return ContainsY(y, wiggle) && ContainsZ(z, wiggle); }
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 179 of file BoxBoundedGeo.h.

180  { return CoordinateContained(z, MinZ(), MaxZ(), wiggle); }
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 355 of file BoxBoundedGeo.h.

356  {
357  return (c >= (min > 0? min / wiggle: min * wiggle))
358  && (c <= (max < 0? max / wiggle: max * wiggle));
359  } // CoordinateContained()
static int max(int a, int b)
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
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 373 of file BoxBoundedGeo.h.

374  { 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 414 of file BoxBoundedGeo.h.

415  { ExtendToInclude(geo::Point_t(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:184
void ExtendToInclude(Coord_t x, Coord_t y, Coord_t z)
Extends the current box to also include the specified point.
list x
Definition: train.py:276
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 421 of file BoxBoundedGeo.h.

422  {
423  set_min(c_min, point);
424  set_max(c_max, point);
425  }
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 433 of file BoxBoundedGeo.h.

434  {
435  set_min(c_min, box.Min());
436  set_max(c_max, box.Max());
437  } // 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 126 of file BoxBoundedGeo.cxx.

127  {
128  std::vector<TVector3> intersections;
129  for (auto const& point: GetIntersections(geo::Point_t(TrajectoryStart), geo::Vector_t(TrajectoryDirect)))
130  intersections.emplace_back(point.X(), point.Y(), point.Z());
131  return intersections;
132  } // 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:164
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:184
std::vector< TVector3 > GetIntersections(TVector3 const &TrajectoryStart, TVector3 const &TrajectoryDirect) const
Calculates the entry and exit points of a trajectory on the box surface.
std::vector< geo::Point_t > geo::BoxBoundedGeo::GetIntersections ( geo::Point_t const &  TrajectoryStart,
geo::Vector_t const &  TrajectoryDirect 
) const

Definition at line 35 of file BoxBoundedGeo.cxx.

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

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

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

Definition at line 115 of file BoxBoundedGeo.h.

115 { 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 130 of file BoxBoundedGeo.h.

130 { 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 241 of file BoxBoundedGeo.h.

242  {
243  return CoordinateContained(x, MinX() + neg_margin, MaxX() - pos_margin);
244  }
double MinX() const
Returns the world x coordinate of the start of the box.
Definition: BoxBoundedGeo.h:88
double MaxX() const
Returns the world x coordinate of the end of the box.
Definition: BoxBoundedGeo.h:91
static bool CoordinateContained(double c, double min, double max, double wiggle=1.)
Returns whether the specified coordinate is in a range.
list x
Definition: train.py:276
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 251 of file BoxBoundedGeo.h.

252  { return InFiducialX(x, margin, margin); }
list x
Definition: train.py:276
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 272 of file BoxBoundedGeo.h.

273  {
274  return CoordinateContained(y, MinY() + neg_margin, MaxY() - pos_margin);
275  }
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 282 of file BoxBoundedGeo.h.

283  { return InFiducialY(y, margin, margin); }
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 303 of file BoxBoundedGeo.h.

304  {
305  return CoordinateContained(z, MinZ() + neg_margin, MaxZ() - pos_margin);
306  }
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 313 of file BoxBoundedGeo.h.

314  { return InFiducialZ(z, margin, margin); }
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 136 of file BoxBoundedGeo.h.

136 { return c_max; }
Coords_t c_max
maximum coordinates (x, y, z)
double geo::BoxBoundedGeo::MaxX ( ) const
inline

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

Definition at line 91 of file BoxBoundedGeo.h.

91 { return c_max.X(); }
Coords_t c_max
maximum coordinates (x, y, z)
double geo::BoxBoundedGeo::MaxY ( ) const
inline

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

Definition at line 106 of file BoxBoundedGeo.h.

106 { return c_max.Y(); }
Coords_t c_max
maximum coordinates (x, y, z)
double geo::BoxBoundedGeo::MaxZ ( ) const
inline

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

Definition at line 121 of file BoxBoundedGeo.h.

121 { return c_max.Z(); }
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 133 of file BoxBoundedGeo.h.

133 { return c_min; }
Coords_t c_min
minimum coordinates (x, y, z)
double geo::BoxBoundedGeo::MinX ( ) const
inline

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

Definition at line 88 of file BoxBoundedGeo.h.

88 { return c_min.X(); }
Coords_t c_min
minimum coordinates (x, y, z)
double geo::BoxBoundedGeo::MinY ( ) const
inline

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

Definition at line 103 of file BoxBoundedGeo.h.

103 { return c_min.Y(); }
Coords_t c_min
minimum coordinates (x, y, z)
double geo::BoxBoundedGeo::MinZ ( ) const
inline

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

Definition at line 118 of file BoxBoundedGeo.h.

118 { return c_min.Z(); }
Coords_t c_min
minimum coordinates (x, y, z)
bool geo::BoxBoundedGeo::Overlaps ( geo::BoxBoundedGeo const &  other) const
inline

Returns if this and other box overlap.

Definition at line 336 of file BoxBoundedGeo.h.

337  { return OverlapsX(other) && OverlapsY(other) && OverlapsZ(other); }
bool OverlapsX(geo::BoxBoundedGeo const &other) const
Returns if the x coordinates covered by this and other box overlap.
bool OverlapsZ(geo::BoxBoundedGeo const &other) const
Returns if the z coordinates covered by this and other box overlap.
bool OverlapsY(geo::BoxBoundedGeo const &other) const
Returns if the y coordinates covered by this and other box overlap.
bool geo::BoxBoundedGeo::OverlapsX ( geo::BoxBoundedGeo const &  other) const
inline

Returns if the x coordinates covered by this and other box overlap.

Definition at line 324 of file BoxBoundedGeo.h.

325  { return std::min(MaxX(), other.MaxX()) > std::max(MinX(), other.MinX()); }
double MinX() const
Returns the world x coordinate of the start of the box.
Definition: BoxBoundedGeo.h:88
double MaxX() const
Returns the world x coordinate of the end of the box.
Definition: BoxBoundedGeo.h:91
static int max(int a, int b)
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
bool geo::BoxBoundedGeo::OverlapsY ( geo::BoxBoundedGeo const &  other) const
inline

Returns if the y coordinates covered by this and other box overlap.

Definition at line 328 of file BoxBoundedGeo.h.

329  { return std::min(MaxY(), other.MaxY()) > std::max(MinY(), other.MinY()); }
static int max(int a, int b)
double MaxY() const
Returns the world y coordinate of the end of the box.
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
double MinY() const
Returns the world y coordinate of the start of the box.
bool geo::BoxBoundedGeo::OverlapsZ ( geo::BoxBoundedGeo const &  other) const
inline

Returns if the z coordinates covered by this and other box overlap.

Definition at line 332 of file BoxBoundedGeo.h.

333  { return std::min(MaxZ(), other.MaxZ()) > std::max(MinZ(), other.MinZ()); }
double MinZ() const
Returns the world z coordinate of the start of the box.
static int max(int a, int b)
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
double MaxZ() const
Returns the world z coordinate of the end of the box.
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 468 of file BoxBoundedGeo.h.

469  { if (value > var) var = value; }
int var
Definition: 018_def.c:9
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 480 of file BoxBoundedGeo.h.

481  {
482  if (value.X() > var.X()) var.SetX(value.X());
483  if (value.Y() > var.Y()) var.SetY(value.Y());
484  if (value.Z() > var.Z()) var.SetZ(value.Z());
485  }
int var
Definition: 018_def.c:9
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 464 of file BoxBoundedGeo.h.

465  { if (value < var) var = value; }
int var
Definition: 018_def.c:9
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 472 of file BoxBoundedGeo.h.

473  {
474  if (value.X() < var.X()) var.SetX(value.X());
475  if (value.Y() < var.Y()) var.SetY(value.Y());
476  if (value.Z() < var.Z()) var.SetZ(value.Z());
477  }
int var
Definition: 018_def.c:9
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 389 of file BoxBoundedGeo.h.

394  {
395  c_min.SetXYZ(x_min, y_min, z_min);
396  c_max.SetXYZ(x_max, y_max, z_max);
397  SortCoordinates();
398  }
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.
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 405 of file BoxBoundedGeo.h.

406  { 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 97 of file BoxBoundedGeo.h.

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

Returns the full size in the Y dimension.

Definition at line 112 of file BoxBoundedGeo.h.

112 { 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 127 of file BoxBoundedGeo.h.

127 { 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 135 of file BoxBoundedGeo.cxx.

135  {
136 
137  for (auto coordMan: geo::vect::coordManagers<geo::Point_t>()) {
138  auto min = geo::vect::bindCoord(c_min, coordMan);
139  auto max = geo::vect::bindCoord(c_max, coordMan);
140  if (min() > max()) {
141  auto temp = min();
142  min = max();
143  max = temp;
144  }
145  } // for
146 
147  } // BoxBoundedGeo::SortCoordinates()
constexpr auto bindCoord(Vector const &v, CoordReader_t< Vector > helper)
Binds the specified constant vector to the coordinate reader.
static int max(int a, int b)
Coords_t c_max
maximum coordinates (x, y, z)
Coords_t c_min
minimum coordinates (x, y, z)
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55

Member Data Documentation

Coords_t geo::BoxBoundedGeo::c_max
private

maximum coordinates (x, y, z)

Definition at line 491 of file BoxBoundedGeo.h.

Coords_t geo::BoxBoundedGeo::c_min
private

minimum coordinates (x, y, z)

Definition at line 490 of file BoxBoundedGeo.h.


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