Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
geo::WireGeo Class Reference

Geometry description of a TPC wireThe wire is a single straight segment on a wire plane. Different wires may be connected to the same readout channel. That is of no relevance for the geometry description. More...

#include <WireGeo.h>

Classes

struct  WireGeoCoordinatesTag
 Tag for vectors in the "local" GDML coordinate frame of the plane. More...
 

Public Types

using GeoNodePath_t = std::vector< TGeoNode const * >
 
Types for geometry-local reference vectors.

These types represents points and displacement vectors in the reference frame defined in the wire geometry "box" from the GDML geometry description.

No alias is explicitly defined for the LArSoft global vector types, geo::Point_t and geo::Vector_t.

Remember the LocalPoint_t and LocalVector_t vectors from different instances of geo::WireGeo have the same type but are not compatible.

using LocalPoint_t = geo::Point3DBase_t< WireGeoCoordinatesTag >
 Type of points in the local GDML wire plane frame. More...
 
using LocalVector_t = geo::Vector3DBase_t< WireGeoCoordinatesTag >
 Type of displacement vectors in the local GDML wire plane frame. More...
 

Public Member Functions

 WireGeo (TGeoNode const &node, geo::TransformationMatrix &&trans)
 Constructor from a ROOT geometry node and a transformation. More...
 
const TGeoNode * Node () const
 
void UpdateAfterSorting (geo::WireID const &, bool flip)
 
Size and coordinates
double RMax () const
 Returns the outer half-size of the wire [cm]. More...
 
double HalfL () const
 Returns half the length of the wire [cm]. More...
 
double RMin () const
 Returns the inner radius of the wire (usually 0) [cm]. More...
 
void GetCenter (double *xyz, double localz=0.0) const
 Fills the world coordinate of a point on the wire. More...
 
void GetStart (double *xyz) const
 
void GetEnd (double *xyz) const
 
template<typename Point >
Point GetPositionFromCenter (double localz) const
 Returns the position (world coordinate) of a point on the wire. More...
 
DefaultPoint_t GetPositionFromCenter (double localz) const
 
template<typename Point >
Point GetPositionFromCenterUnbounded (double localz) const
 Returns the position (world coordinate) of a point on the wire. More...
 
DefaultPoint_t GetPositionFromCenterUnbounded (double localz) const
 
template<typename Point >
Point GetCenter () const
 
DefaultPoint_t GetCenter () const
 
template<typename Point >
Point GetStart () const
 
DefaultPoint_t GetStart () const
 
template<typename Point >
Point GetEnd () const
 
DefaultPoint_t GetEnd () const
 
double Length () const
 Returns the wire length in centimeters. More...
 
Orientation and angles
double ThetaZ () const
 Returns angle of wire with respect to z axis in the Y-Z plane in radians. More...
 
double ThetaZ (bool degrees) const
 
double CosThetaZ () const
 Returns trigonometric operations on ThetaZ() More...
 
double SinThetaZ () const
 
double TanThetaZ () const
 
bool isHorizontal () const
 Returns if this wire is horizontal (theta_z ~ 0) More...
 
bool isVertical () const
 Returns if this wire is vertical (theta_z ~ pi/2) More...
 
bool isParallelTo (geo::WireGeo const &wire) const
 Returns if this wire is parallel to another. More...
 
template<typename Vector >
Vector Direction () const
 
DefaultVector_t Direction () const
 
Coordinate transformation

Local points and displacement vectors are described by the types geo::WireGeo::LocalPoint_t and geo::WireGeo::LocalVector_t, respectively.

void LocalToWorld (const double *wire, double *world) const
 Transform point from local wire frame to world frame. More...
 
geo::Point_t toWorldCoords (LocalPoint_t const &local) const
 Transform point from local wire frame to world frame. More...
 
void LocalToWorldVect (const double *wire, double *world) const
 Transform direction vector from local to world. More...
 
geo::Vector_t toWorldCoords (LocalVector_t const &local) const
 Transform direction vector from local to world. More...
 
void WorldToLocal (const double *world, double *wire) const
 Transform point from world frame to local wire frame. More...
 
LocalPoint_t toLocalCoords (geo::Point_t const &world) const
 Transform point from world frame to local wire frame. More...
 
void WorldToLocalVect (const double *world, double *wire) const
 Transform direction vector from world to local. More...
 
LocalVector_t toLocalCoords (geo::Vector_t const &world) const
 Transform direction vector from world to local. More...
 
Geometric properties and algorithms
double ComputeZatY0 () const
 
double DistanceFrom (geo::WireGeo const &wire) const
 Returns 3D distance from the specified wire. More...
 
template<typename Point = DefaultPoint_t>
Point IntersectionWith (geo::WireGeo const &other) const
 Returns the point of this wire that is closest to other wire. More...
 
template<typename Point = DefaultPoint_t>
geo::IntersectionPointAndOffsets< PointIntersectionAndOffsetsWith (geo::WireGeo const &other) const
 Returns the point of this wire that is closest to other wire. More...
 

Static Public Member Functions

static double WirePitch (geo::WireGeo const &w1, geo::WireGeo const &w2)
 Returns the pitch (distance on y/z plane) between two wires, in cm. More...
 

Private Types

using DefaultVector_t = TVector3
 
using DefaultPoint_t = TVector3
 
using LocalTransformation_t = geo::LocalTransformationGeo< ROOT::Math::Transform3D, LocalPoint_t, LocalVector_t >
 

Private Member Functions

bool isFlipped () const
 
double relLength (double local) const
 Returns the relative length from center to be used when transforming. More...
 
double capLength (double local) const
 Caps the specified local length coordinate to lay on the wire. More...
 
double capRelLength (double local) const
 Stacked capLength() and relLength(). More...
 
void Flip ()
 Set to swap the start and end wire. More...
 

Static Private Member Functions

static double gausSum (double a, double b)
 

Private Attributes

const TGeoNode * fWireNode
 Pointer to the wire node. More...
 
double fThetaZ
 angle of the wire with respect to the z direction More...
 
double fHalfL
 half length of the wire More...
 
geo::Point_t fCenter
 Center of the wire in world coordinates. More...
 
LocalTransformation_t fTrans
 Wire to world transform. More...
 
bool flipped
 whether start and end are reversed More...
 

Printing

static constexpr unsigned int MaxVerbosity = 4
 Maximum verbosity supported by PrintWireInfo(). More...
 
template<typename Stream >
void PrintWireInfo (Stream &&out, std::string indent="", unsigned int verbosity=1) const
 Prints information about this wire. More...
 
std::string WireInfo (std::string indent="", unsigned int verbosity=1) const
 Returns a string with all the information of the wire. More...
 

Detailed Description

Geometry description of a TPC wire

The wire is a single straight segment on a wire plane. Different wires may be connected to the same readout channel. That is of no relevance for the geometry description.


The wire has a start and an end point. Their definition of them is related to the other wires in the plane and to the TPC itself.

The direction of increasing wire coordinate, defined in the wire plane, is orthogonal to the wire direction and of course points to the direction where the wire number within the plane increases. This direction is indirectly defined when sorting the wires in the plane, which is done by the plane (geo::PlaneGeo). This direction lies by definition on the wire plane. The direction normal to the wire plane is defined by the TPC so that it points inward the TPC rather than outward. Finally, the wire direction is defined so that the triplet of unit vectors direction of the wire $ \hat{l} $, direction of increasing wire number $ \hat{w} $, and normal to the plane $ \hat{n} $ is positively defined ( $ \hat{l} \times \hat{w} \cdot \hat{n} = +1 $). The start $ \vec{a}_{w} $ and the end of the wire $ \vec{b}_{w} $ are defined so that their difference $ \vec{b}_{w} - \vec{a}_{w} $ points in the same direction as $ \hat{l} $.

Definition at line 65 of file WireGeo.h.

Member Typedef Documentation

using geo::WireGeo::DefaultPoint_t = TVector3
private

Definition at line 68 of file WireGeo.h.

using geo::WireGeo::DefaultVector_t = TVector3
private

Definition at line 67 of file WireGeo.h.

using geo::WireGeo::GeoNodePath_t = std::vector<TGeoNode const*>

Definition at line 72 of file WireGeo.h.

Type of points in the local GDML wire plane frame.

Definition at line 94 of file WireGeo.h.

Definition at line 481 of file WireGeo.h.

Type of displacement vectors in the local GDML wire plane frame.

Definition at line 97 of file WireGeo.h.

Constructor & Destructor Documentation

geo::WireGeo::WireGeo ( TGeoNode const &  node,
geo::TransformationMatrix &&  trans 
)

Constructor from a ROOT geometry node and a transformation.

Parameters
nodeROOT geometry node
transtransformation matrix (local to world)

The node describes the shape of the wire (the only relevant information is in fact the length), while the transformation described its positioning in the world (both position and orientation).

A pointer to the node and a copy of the transformation matrix are kept in the WireGeo object.

Definition at line 32 of file WireGeo.cxx.

33  : fWireNode(&node)
34  , fTrans(std::move(trans))
35  , flipped(false)
36  {
37  fHalfL = ((TGeoTube*)fWireNode->GetVolume()->GetShape())->GetDZ();
38 
39  // uncomment the following to check the paths to the wires
40  // std::string p(base);
41  // for(int i = 0; i <= depth; ++i){
42  // p += "/";
43  // p += path[i]->GetName();
44  // }
45  // std::cout << p.c_str() << std::endl;
46 
47  // determine the orientation of the wire
48  auto lp = geo::origin<LocalPoint_t>();
49  fCenter = toWorldCoords(lp);
50 
51  lp.SetZ(fHalfL);
52  auto end = toWorldCoords(lp);
53 
54  fThetaZ = std::acos(std::clamp((end.Z() - fCenter.Z())/fHalfL, -1.0, +1.0));
55 
56  // check to see if it runs "forward" or "backwards" in z
57  // check is made looking at the y position of the end point
58  // relative to the center point because you want to know if
59  // the end point is above or below the center of the wire in
60  // the yz plane
61  if(end.Y() < fCenter.Y()) fThetaZ *= -1.;
62 
63  //This ensures we are looking at the angle between 0 and Pi
64  //as if the wire runs at one angle it also runs at that angle +-Pi
65  if(fThetaZ < 0) fThetaZ += util::pi();
66 
67  assert(std::isfinite(fThetaZ));
68 
69  } // geo::WireGeo::WireGeo()
end
while True: pbar.update(maxval-len(onlies[E][S])) #print iS, "/", len(onlies[E][S]) found = False for...
double fThetaZ
angle of the wire with respect to the z direction
Definition: WireGeo.h:484
geo::Point_t fCenter
Center of the wire in world coordinates.
Definition: WireGeo.h:486
geo::Point_t toWorldCoords(LocalPoint_t const &local) const
Transform point from local wire frame to world frame.
Definition: WireGeo.h:359
double fHalfL
half length of the wire
Definition: WireGeo.h:485
bool isfinite(Vector const &v)
Returns whether all components of the vector are finite.
def move(depos, offset)
Definition: depos.py:107
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:487
constexpr T pi()
Returns the constant pi (up to 35 decimal digits of precision)
const TGeoNode * fWireNode
Pointer to the wire node.
Definition: WireGeo.h:483
bool flipped
whether start and end are reversed
Definition: WireGeo.h:488

Member Function Documentation

double geo::WireGeo::capLength ( double  local) const
inlineprivate

Caps the specified local length coordinate to lay on the wire.

Definition at line 498 of file WireGeo.h.

499  { return std::min(+HalfL(), std::max(-HalfL(), local)); }
static int max(int a, int b)
double HalfL() const
Returns half the length of the wire [cm].
Definition: WireGeo.h:128
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
double geo::WireGeo::capRelLength ( double  local) const
inlineprivate

Stacked capLength() and relLength().

Definition at line 502 of file WireGeo.h.

503  { return capLength(relLength(local)); }
double relLength(double local) const
Returns the relative length from center to be used when transforming.
Definition: WireGeo.h:495
double capLength(double local) const
Caps the specified local length coordinate to lay on the wire.
Definition: WireGeo.h:498
double geo::WireGeo::ComputeZatY0 ( ) const
inline

Returns the z coordinate, in centimetres, at the point where y = 0. Assumes the wire orthogonal to x axis and the wire not parallel to z.

Definition at line 399 of file WireGeo.h.

400  { return fCenter.Z() - fCenter.Y() / TanThetaZ(); }
geo::Point_t fCenter
Center of the wire in world coordinates.
Definition: WireGeo.h:486
double TanThetaZ() const
Definition: WireGeo.h:266
double geo::WireGeo::CosThetaZ ( ) const
inline

Returns trigonometric operations on ThetaZ()

Definition at line 264 of file WireGeo.h.

264 { return std::cos(ThetaZ()); }
double ThetaZ() const
Returns angle of wire with respect to z axis in the Y-Z plane in radians.
Definition: WireGeo.h:250
template<typename Vector >
Vector geo::WireGeo::Direction ( ) const

Returns the wire direction as a norm-one vector.

Template Parameters
Vectortype of the vector being returned

Definition at line 587 of file WireGeo.h.

587  {
588  // maybe (GetCenter() - GetStart()) / HalfL() would be faster;
589  // strangely, TVector3 does not implement operator/ (double).
590  return geo::vect::convertTo<Vector>(GetEnd<geo::Point_t>() - GetStart<geo::Point_t>()) * (1.0 / Length());
591 } // geo::WireGeo::Direction()
double Length() const
Returns the wire length in centimeters.
Definition: WireGeo.h:236
DefaultVector_t geo::WireGeo::Direction ( ) const
inline

Definition at line 293 of file WireGeo.h.

293 { return Direction<DefaultVector_t>(); }
double geo::WireGeo::DistanceFrom ( geo::WireGeo const &  wire) const

Returns 3D distance from the specified wire.

Returns
the signed distance in centimetres (0 if wires are not parallel)

If the specified wire is "ahead" in z respect to this, the distance is returned negative.

Definition at line 113 of file WireGeo.cxx.

113  {
114  //
115  // The algorithm assumes that picking any point on the wire will do,
116  // that is, that the wires are parallel.
117  //
118 
119  if (!isParallelTo(wire)) return 0;
120 
121  // a vector connecting to the other wire
122  auto const toWire = wire.GetCenter() - GetCenter();
123 
124  // the distance is that vector, times the sine of the angle with the wire
125  // direction; we get that in a generic way with a cross product.
126  // We don't even care about the sign here
127  // (if we did, we would do a dot-product with the normal to the plane,
128  // and we should get a positive distance if the other wire has larger wire
129  // coordinate than this one).
130  return toWire.Cross(Direction()).Mag();
131 
132  } // WireGeo::DistanceFrom()
Point GetCenter() const
Definition: WireGeo.h:212
Vector Direction() const
Definition: WireGeo.h:587
bool isParallelTo(geo::WireGeo const &wire) const
Returns if this wire is parallel to another.
Definition: WireGeo.h:281
void geo::WireGeo::Flip ( )
private

Set to swap the start and end wire.

Definition at line 144 of file WireGeo.cxx.

144  {
145  // we don't need to do much to flip so far:
146  // - ThetaZ() is defined in [0, pi], invariant to flipping
147  // - we don't change the transformation matrices, that we want to be
148  // untouched and coherent with the original geometry source
149  // - center is invariant for flipping
150  // - start and end are computed on the fly (taking flipping into account)
151  // - ... and we chose to leave half length unsigned and independent
152 
153  // change the flipping bit
154  flipped = !flipped;
155 
156  } // WireGeo::Flip()
bool flipped
whether start and end are reversed
Definition: WireGeo.h:488
static double geo::WireGeo::gausSum ( double  a,
double  b 
)
inlinestaticprivate

Definition at line 509 of file WireGeo.h.

509 { return std::sqrt(a*a + b*b); }
const double a
static bool * b
Definition: config.cpp:1043
void geo::WireGeo::GetCenter ( double *  xyz,
double  localz = 0.0 
) const

Fills the world coordinate of a point on the wire.

Parameters
xyz_(output)_ the position to be filled, as [ x, y, z ] (in cm)
localzdistance of the requested point from the middle of the wire
See also
GetCenter(), GetStart(), GetEnd(), GetPositionFromCenter()

The center of the wires corresponds to localz equal to 0; negative positions head toward the start of the wire, positive toward the end.

If the localz position would put the point outside the wire, the returned position is the wire end closest to the requested position.

Deprecated:
Use the version returning a vector instead.

Definition at line 73 of file WireGeo.cxx.

74  {
75  if (localz == 0.) { // if no dislocation is requested, we already have it
77  return;
78  }
79 
80  double locz = relLength(localz);
81  if (std::abs(locz) > fHalfL) {
82  mf::LogWarning("WireGeo") << "asked for position along wire that"
83  " extends beyond the wire, returning position at end point";
84  locz = relLength((locz < 0)? -fHalfL: fHalfL);
85  }
86  const double local[3] = { 0., 0., locz };
87  LocalToWorld(local, xyz);
88  }
geo::Point_t fCenter
Center of the wire in world coordinates.
Definition: WireGeo.h:486
unsigned int fillCoords(Coords &dest, Vector const &src)
Fills a coordinate array with the coordinates of a vector.
double fHalfL
half length of the wire
Definition: WireGeo.h:485
T abs(T value)
void LocalToWorld(const double *wire, double *world) const
Transform point from local wire frame to world frame.
Definition: WireGeo.h:355
double relLength(double local) const
Returns the relative length from center to be used when transforming.
Definition: WireGeo.h:495
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
template<typename Point >
Point geo::WireGeo::GetCenter ( ) const
inline

Returns the world coordinate of the center of the wire [cm]

Template Parameters
Pointtype of the point being returned

Definition at line 212 of file WireGeo.h.

212 { return geo::vect::convertTo<Point>(fCenter); }
geo::Point_t fCenter
Center of the wire in world coordinates.
Definition: WireGeo.h:486
DefaultPoint_t geo::WireGeo::GetCenter ( ) const
inline

Definition at line 213 of file WireGeo.h.

213 { return GetCenter<DefaultPoint_t>(); }
void geo::WireGeo::GetEnd ( double *  xyz) const
inline

Fills the world coordinate of one end of the wire

Deprecated:
Use the version returning a vector instead.

Definition at line 163 of file WireGeo.h.

163 { GetCenter(xyz, +fHalfL); }
Point GetCenter() const
Definition: WireGeo.h:212
double fHalfL
half length of the wire
Definition: WireGeo.h:485
template<typename Point >
Point geo::WireGeo::GetEnd ( ) const
inline

Returns the world coordinate of one end of the wire [cm]

Template Parameters
Pointtype of the point being returned

Definition at line 229 of file WireGeo.h.

230  { return GetPositionFromCenterUnbounded<Point>(+HalfL()); }
double HalfL() const
Returns half the length of the wire [cm].
Definition: WireGeo.h:128
DefaultPoint_t geo::WireGeo::GetEnd ( ) const
inline

Definition at line 231 of file WireGeo.h.

231 { return GetEnd<DefaultPoint_t>(); }
template<typename Point >
Point geo::WireGeo::GetPositionFromCenter ( double  localz) const
inline

Returns the position (world coordinate) of a point on the wire.

Template Parameters
Pointtype of vector to be returned (current default: TVector3)
Parameters
localzdistance of the requested point from the middle of the wire
Returns
the position of the requested point (in cm)
See also
GetCenter(), GetStart(), GetEnd(), GetPositionFromCenterUnbounded()

The center of the wires corresponds to localz equal to 0; negative positions head toward the start of the wire, positive toward the end.

If the localz position would put the point outside the wire, the returned position is the wire end closest to the requested position.

Definition at line 182 of file WireGeo.h.

183  { return GetPositionFromCenterUnbounded<Point>(capLength(localz)); }
double capLength(double local) const
Caps the specified local length coordinate to lay on the wire.
Definition: WireGeo.h:498
DefaultPoint_t geo::WireGeo::GetPositionFromCenter ( double  localz) const
inline

Definition at line 184 of file WireGeo.h.

185  { return GetPositionFromCenter<DefaultPoint_t>(localz); }
template<typename Point >
Point geo::WireGeo::GetPositionFromCenterUnbounded ( double  localz) const

Returns the position (world coordinate) of a point on the wire.

Template Parameters
Pointtype of vector to be returned (current default: TVector3)
Parameters
localzdistance of the requested point from the middle of the wire
Returns
the position of the requested point (in cm)
See also
GetCenter(), GetStart(), GetEnd(), GetPositionFromCenter()

The center of the wires corresponds to localz equal to 0; negative positions head toward the start of the wire, positive toward the end.

If the localz position would put the point outside the wire, the returned position will lie beyond the end of the wire.

Definition at line 579 of file WireGeo.h.

579  {
580  return geo::vect::convertTo<Point>
581  (toWorldCoords(LocalPoint_t{ 0.0, 0.0, relLength(localz) }));
582 } // geo::WireGeo::GetPositionFromCenterImpl()
geo::Point_t toWorldCoords(LocalPoint_t const &local) const
Transform point from local wire frame to world frame.
Definition: WireGeo.h:359
double relLength(double local) const
Returns the relative length from center to be used when transforming.
Definition: WireGeo.h:495
geo::Point3DBase_t< WireGeoCoordinatesTag > LocalPoint_t
Type of points in the local GDML wire plane frame.
Definition: WireGeo.h:94
DefaultPoint_t geo::WireGeo::GetPositionFromCenterUnbounded ( double  localz) const
inline

Definition at line 204 of file WireGeo.h.

205  { return GetPositionFromCenterUnbounded<DefaultPoint_t>(localz); }
void geo::WireGeo::GetStart ( double *  xyz) const
inline

Fills the world coordinate of one end of the wire

Deprecated:
Use the version returning a vector instead.

Definition at line 157 of file WireGeo.h.

157 { GetCenter(xyz, -fHalfL); }
Point GetCenter() const
Definition: WireGeo.h:212
double fHalfL
half length of the wire
Definition: WireGeo.h:485
template<typename Point >
Point geo::WireGeo::GetStart ( ) const
inline

Returns the world coordinate of one end of the wire [cm]

Template Parameters
Pointtype of the point being returned

Definition at line 220 of file WireGeo.h.

221  { return GetPositionFromCenterUnbounded<Point>(-HalfL()); }
double HalfL() const
Returns half the length of the wire [cm].
Definition: WireGeo.h:128
DefaultPoint_t geo::WireGeo::GetStart ( ) const
inline

Definition at line 222 of file WireGeo.h.

222 { return GetStart<DefaultPoint_t>(); }
double geo::WireGeo::HalfL ( ) const
inline

Returns half the length of the wire [cm].

Definition at line 128 of file WireGeo.h.

128 { return fHalfL; }
double fHalfL
half length of the wire
Definition: WireGeo.h:485
template<typename Point >
geo::IntersectionPointAndOffsets< Point > geo::WireGeo::IntersectionAndOffsetsWith ( geo::WireGeo const &  other) const

Returns the point of this wire that is closest to other wire.

Template Parameters
Pointthe type of point returned
Parameters
otherthe other wire
Returns
a data structure with three fields: point: the point of this wire closest to other; offset1: its offset on this wire [cm]; offset2: its offset on the other wire [cm]
See also
IntersectionWith()

The point of this wire that is closest to any point of the other wire is returned. The returned intersection point is the same as for IntersectionWith(other). The return value is actually triplet, though, which is most easily unpacked immediately:

auto [ point, offset, otherOffset ]
= wire.IntersectionAndOffsetsWith(otherWire);

The two other elements of the triplets are the distances of the intersection point from the center of this wire (offset in the example) and from the center of the other wire (otherOffset), in centimeters. The sign of the offsets are positive if the intersection points lie on the side pointed by the Direction() of the respective wires.

To reassign the variables after they have been defined, instead:

std::tie(point, otherOffset, offset)
= otherWire.IntersectionAndOffsetsWith(wire);

Definition at line 645 of file WireGeo.h.

646 {
647  auto const& [ point, ofsA, ofsB ] = WiresIntersectionAndOffsets(*this, other);
648  return { geo::vect::convertTo<Point>(point), ofsA, ofsB };
649 } // geo::WireGeo::IntersectionAndOffsetsWith()
geo::IntersectionPointAndOffsets< geo::Point_t > WiresIntersectionAndOffsets(geo::WireGeo const &wireA, geo::WireGeo const &wireB)
Returns the point of wireA that is closest to wireB.
Definition: WireGeo.h:668
template<typename Point >
Point geo::WireGeo::IntersectionWith ( geo::WireGeo const &  other) const

Returns the point of this wire that is closest to other wire.

Template Parameters
Pointthe type of point returned
Parameters
otherthe other wire
Returns
the point of this wire closest to other
See also
IntersectionAndOffsetsWith()

The point of this wire that is closest to any point of the other wire is returned.

The other wire is assumed not to be parallel to this one, and when this prerequisite is not met the behaviour is undefined.

Another method, IntersectionAndOffsetsWith(), also returns the offset of the intersection from the two wire centers.

Definition at line 637 of file WireGeo.h.

637  {
638  return geo::vect::convertTo<Point>(WiresIntersection(*this, other));
639 } // geo::WireGeo::IntersectionWith()
geo::Point_t WiresIntersection(geo::WireGeo const &wireA, geo::WireGeo const &wireB)
Returns the point of wireA that is closest to wireB.
Definition: WireGeo.h:654
bool geo::WireGeo::isFlipped ( ) const
inlineprivate

Returns whether ( 0, 0, fHalfL ) identifies end (false) or start (true) of the wire.

Definition at line 492 of file WireGeo.h.

492 { return flipped; }
bool flipped
whether start and end are reversed
Definition: WireGeo.h:488
bool geo::WireGeo::isHorizontal ( ) const
inline

Returns if this wire is horizontal (theta_z ~ 0)

Definition at line 271 of file WireGeo.h.

271 { return std::abs(SinThetaZ()) < 1e-5; }
double SinThetaZ() const
Definition: WireGeo.h:265
T abs(T value)
const double e
bool geo::WireGeo::isParallelTo ( geo::WireGeo const &  wire) const
inline

Returns if this wire is parallel to another.

Definition at line 281 of file WireGeo.h.

282  {
283  return // parallel if the dot product of the directions is about +/- 1
284  std::abs(std::abs(Direction<geo::Vector_t>().Dot(wire.Direction<geo::Vector_t>())) - 1.) < 1e-5;
285  }
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
T abs(T value)
const double e
bool geo::WireGeo::isVertical ( ) const
inline

Returns if this wire is vertical (theta_z ~ pi/2)

Definition at line 276 of file WireGeo.h.

276 { return std::abs(CosThetaZ()) < 1e-5; }
T abs(T value)
const double e
double CosThetaZ() const
Returns trigonometric operations on ThetaZ()
Definition: WireGeo.h:264
double geo::WireGeo::Length ( ) const
inline

Returns the wire length in centimeters.

Definition at line 236 of file WireGeo.h.

236 { return 2. * HalfL(); }
double HalfL() const
Returns half the length of the wire [cm].
Definition: WireGeo.h:128
void geo::WireGeo::LocalToWorld ( const double *  wire,
double *  world 
) const
inline

Transform point from local wire frame to world frame.

Definition at line 355 of file WireGeo.h.

356  { fTrans.LocalToWorld(wire, world); }
void LocalToWorld(double const *local, double *world) const
Transforms a point from local frame to world frame.
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:487
void geo::WireGeo::LocalToWorldVect ( const double *  wire,
double *  world 
) const
inline

Transform direction vector from local to world.

Definition at line 363 of file WireGeo.h.

364  { fTrans.LocalToWorldVect(wire, world); }
void LocalToWorldVect(double const *local, double *world) const
Transforms a vector from local frame to world frame.
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:487
const TGeoNode* geo::WireGeo::Node ( ) const
inline

Definition at line 390 of file WireGeo.h.

390 { return fWireNode; }
const TGeoNode * fWireNode
Pointer to the wire node.
Definition: WireGeo.h:483
template<typename Stream >
void geo::WireGeo::PrintWireInfo ( Stream &&  out,
std::string  indent = "",
unsigned int  verbosity = 1 
) const

Prints information about this wire.

Template Parameters
Streamtype of output stream to use
Parameters
outstream to send the information to
indentprepend each line with this string
verbosityamount of information printed

Note that the first line out the output is not indented.

Verbosity levels

  • 0: only start and end
  • 1 _(default)_: also length
  • 2: also angle with z axis
  • 3: also center
  • 4: also direction

The constant MaxVerbosity is set to the highest supported verbosity level.

Definition at line 596 of file WireGeo.h.

600  {
601 
602  //----------------------------------------------------------------------------
603  out << "wire from " << GetStart<geo::Point_t>()
604  << " to " << GetEnd<geo::Point_t>();
605 
606  if (verbosity-- <= 0) return; // 0
607 
608  //----------------------------------------------------------------------------
609  out << " (" << Length() << " cm long)";
610 
611  if (verbosity-- <= 0) return; // 1
612 
613  //----------------------------------------------------------------------------
614  out << ", theta(z)=" << ThetaZ() << " rad";
615 
616  if (verbosity-- <= 0) return; // 2
617 
618  //----------------------------------------------------------------------------
619  out << "\n" << indent
620  << " center at " << GetCenter<geo::Point_t>() << " cm";
621 
622  if (verbosity-- <= 0) return; // 3
623 
624  //----------------------------------------------------------------------------
625  out << ", direction: " << Direction<geo::Vector_t>();
626  if (isHorizontal()) out << " (horizontal)";
627  if (isVertical()) out << " (vertical)";
628 
629 // if (verbosity-- <= 0) return; // 4
630 
631  //----------------------------------------------------------------------------
632 } // geo::WireGeo::PrintWireInfo()
bool isVertical() const
Returns if this wire is vertical (theta_z ~ pi/2)
Definition: WireGeo.h:276
double Length() const
Returns the wire length in centimeters.
Definition: WireGeo.h:236
double ThetaZ() const
Returns angle of wire with respect to z axis in the Y-Z plane in radians.
Definition: WireGeo.h:250
bool isHorizontal() const
Returns if this wire is horizontal (theta_z ~ 0)
Definition: WireGeo.h:271
double geo::WireGeo::relLength ( double  local) const
inlineprivate

Returns the relative length from center to be used when transforming.

Definition at line 495 of file WireGeo.h.

495 { return isFlipped()? -local: local; }
bool isFlipped() const
Definition: WireGeo.h:492
double geo::WireGeo::RMax ( ) const

Returns the outer half-size of the wire [cm].

Definition at line 91 of file WireGeo.cxx.

92  { return ((TGeoTube*)fWireNode->GetVolume()->GetShape())->GetRmax(); }
const TGeoNode * fWireNode
Pointer to the wire node.
Definition: WireGeo.h:483
double geo::WireGeo::RMin ( ) const

Returns the inner radius of the wire (usually 0) [cm].

Definition at line 95 of file WireGeo.cxx.

96  { return ((TGeoTube*)fWireNode->GetVolume()->GetShape())->GetRmin(); }
const TGeoNode * fWireNode
Pointer to the wire node.
Definition: WireGeo.h:483
double geo::WireGeo::SinThetaZ ( ) const
inline

Definition at line 265 of file WireGeo.h.

265 { return std::sin(ThetaZ()); }
double ThetaZ() const
Returns angle of wire with respect to z axis in the Y-Z plane in radians.
Definition: WireGeo.h:250
double geo::WireGeo::TanThetaZ ( ) const
inline

Definition at line 266 of file WireGeo.h.

266 { return std::tan(ThetaZ()); }
double ThetaZ() const
Returns angle of wire with respect to z axis in the Y-Z plane in radians.
Definition: WireGeo.h:250
double geo::WireGeo::ThetaZ ( ) const
inline

Returns angle of wire with respect to z axis in the Y-Z plane in radians.

Definition at line 250 of file WireGeo.h.

250 { return fThetaZ; }
double fThetaZ
angle of the wire with respect to the z direction
Definition: WireGeo.h:484
double geo::WireGeo::ThetaZ ( bool  degrees) const

Returns angle of wire with respect to z axis in the Y-Z plane

Parameters
degreesreturn the angle in degrees rather than radians
Returns
wire angle

Definition at line 99 of file WireGeo.cxx.

100  { return degrees? util::RadiansToDegrees(fThetaZ): fThetaZ; }
double fThetaZ
angle of the wire with respect to the z direction
Definition: WireGeo.h:484
constexpr T RadiansToDegrees(T angle)
Converts the argument angle from radians into degrees ( )
LocalPoint_t geo::WireGeo::toLocalCoords ( geo::Point_t const &  world) const
inline

Transform point from world frame to local wire frame.

Definition at line 375 of file WireGeo.h.

376  { return fTrans.toLocalCoords(world); }
LocalPoint_t toLocalCoords(GlobalPoint_t const &world) const
Transforms a point from world frame to local frame.
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:487
LocalVector_t geo::WireGeo::toLocalCoords ( geo::Vector_t const &  world) const
inline

Transform direction vector from world to local.

Definition at line 383 of file WireGeo.h.

384  { return fTrans.toLocalCoords(world); }
LocalPoint_t toLocalCoords(GlobalPoint_t const &world) const
Transforms a point from world frame to local frame.
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:487
geo::Point_t geo::WireGeo::toWorldCoords ( LocalPoint_t const &  local) const
inline

Transform point from local wire frame to world frame.

Definition at line 359 of file WireGeo.h.

360  { return fTrans.toWorldCoords(local); }
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:487
GlobalPoint_t toWorldCoords(LocalPoint_t const &local) const
Transforms a point from local frame to world frame.
geo::Vector_t geo::WireGeo::toWorldCoords ( LocalVector_t const &  local) const
inline

Transform direction vector from local to world.

Definition at line 367 of file WireGeo.h.

368  { return fTrans.toWorldCoords(local); }
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:487
GlobalPoint_t toWorldCoords(LocalPoint_t const &local) const
Transforms a point from local frame to world frame.
void geo::WireGeo::UpdateAfterSorting ( geo::WireID const &  ,
bool  flip 
)

Internal updates after the relative position of the wire is known (currently no-op)

Definition at line 136 of file WireGeo.cxx.

136  {
137 
138  // flip, if asked
139  if (flip) Flip();
140 
141  } // WireGeo::UpdateAfterSorting()
void Flip()
Set to swap the start and end wire.
Definition: WireGeo.cxx:144
std::string geo::WireGeo::WireInfo ( std::string  indent = "",
unsigned int  verbosity = 1 
) const

Returns a string with all the information of the wire.

See also
PrintWireInfo()

All arguments are equivalent to the ones of PrintWireInfo.

Definition at line 104 of file WireGeo.cxx.

105  {
106  std::ostringstream sstr;
107  PrintWireInfo(sstr, indent, verbosity);
108  return sstr.str();
109  } // WireGeo::WireInfo()
void PrintWireInfo(Stream &&out, std::string indent="", unsigned int verbosity=1) const
Prints information about this wire.
Definition: WireGeo.h:596
static double geo::WireGeo::WirePitch ( geo::WireGeo const &  w1,
geo::WireGeo const &  w2 
)
inlinestatic

Returns the pitch (distance on y/z plane) between two wires, in cm.

Definition at line 476 of file WireGeo.h.

477  { return std::abs(w2.DistanceFrom(w1)); }
T abs(T value)
void geo::WireGeo::WorldToLocal ( const double *  world,
double *  wire 
) const
inline

Transform point from world frame to local wire frame.

Definition at line 371 of file WireGeo.h.

372  { fTrans.WorldToLocal(world, wire); }
void WorldToLocal(double const *world, double *local) const
Transforms a point from world frame to local frame.
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:487
void geo::WireGeo::WorldToLocalVect ( const double *  world,
double *  wire 
) const
inline

Transform direction vector from world to local.

Definition at line 379 of file WireGeo.h.

380  { fTrans.WorldToLocalVect(world, wire); }
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:487
void WorldToLocalVect(const double *world, double *local) const
Transforms a vector from world frame to local frame.

Member Data Documentation

geo::Point_t geo::WireGeo::fCenter
private

Center of the wire in world coordinates.

Definition at line 486 of file WireGeo.h.

double geo::WireGeo::fHalfL
private

half length of the wire

Definition at line 485 of file WireGeo.h.

bool geo::WireGeo::flipped
private

whether start and end are reversed

Definition at line 488 of file WireGeo.h.

double geo::WireGeo::fThetaZ
private

angle of the wire with respect to the z direction

Definition at line 484 of file WireGeo.h.

LocalTransformation_t geo::WireGeo::fTrans
private

Wire to world transform.

Definition at line 487 of file WireGeo.h.

const TGeoNode* geo::WireGeo::fWireNode
private

Pointer to the wire node.

Definition at line 483 of file WireGeo.h.

constexpr unsigned int geo::WireGeo::MaxVerbosity = 4
static

Maximum verbosity supported by PrintWireInfo().

Definition at line 338 of file WireGeo.h.


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