LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
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 (GeoNodePath_t const &path, size_t depth)
 
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 = DefaultVector_t>
Vector Direction () const
 Returns the wire direction as a norm-one vector. More...
 
template<typename Stream >
void PrintWireInfo (Stream &&out, std::string indent="", unsigned int verbosity=1) const
 Prints information about this wire. More...
 
const TGeoNode * Node () const
 
double ComputeZatY0 () const
 
double DistanceFrom (geo::WireGeo const &wire) const
 Returns 3D distance from the specified wire. More...
 
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 = DefaultPoint_t>
Point GetPositionFromCenter (double localz) const
 Returns the position (world coordinate) of a point on the wire. More...
 
template<typename Point = DefaultPoint_t>
Point GetPositionFromCenterUnbounded (double localz) const
 Returns the position (world coordinate) of a point on the wire. More...
 
template<typename Point = DefaultPoint_t>
Point GetCenter () const
 Returns the world coordinate of the center of the wire [cm]. More...
 
template<typename Point = DefaultPoint_t>
Point GetStart () const
 Returns the world coordinate of one end of the wire [cm]. More...
 
template<typename Point = DefaultPoint_t>
Point GetEnd () const
 Returns the world coordinate of one end of the wire [cm]. More...
 
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
 Returns angle of wire with respect to z axis in the Y-Z plane in radians. More...
 
double TanThetaZ () const
 Returns angle of wire with respect to z axis in the Y-Z plane in radians. More...
 
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...
 

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...
 

Static Public Attributes

static constexpr unsigned int MaxVerbosity = 4
 Maximum verbosity supported by PrintWireInfo(). More...
 

Private Types

using DefaultVector_t = TVector3
 
using DefaultPoint_t = TVector3
 
using LocalTransformation_t = geo::LocalTransformationGeo< TGeoHMatrix, 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...
 

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 61 of file WireGeo.h.

Member Typedef Documentation

using geo::WireGeo::DefaultPoint_t = TVector3
private

Definition at line 64 of file WireGeo.h.

using geo::WireGeo::DefaultVector_t = TVector3
private

Definition at line 63 of file WireGeo.h.

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

Definition at line 68 of file WireGeo.h.

Type of points in the local GDML wire plane frame.

Definition at line 89 of file WireGeo.h.

Definition at line 330 of file WireGeo.h.

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

Definition at line 92 of file WireGeo.h.

Constructor & Destructor Documentation

geo::WireGeo::WireGeo ( GeoNodePath_t const &  path,
size_t  depth 
)

Definition at line 29 of file WireGeo.cxx.

References evd::details::end(), fCenter, fHalfL, fThetaZ, fWireNode, util::pi(), and toWorldCoords().

30  : fTrans(path, depth)
31  , flipped(false)
32  {
33  fWireNode = path[depth];
34  fHalfL = ((TGeoTube*)fWireNode->GetVolume()->GetShape())->GetDZ();
35 
36  // uncomment the following to check the paths to the wires
37  // std::string p(base);
38  // for(int i = 0; i <= depth; ++i){
39  // p += "/";
40  // p += path[i]->GetName();
41  // }
42  // std::cout << p.c_str() << std::endl;
43 
44  // determine the orientation of the wire
45  auto lp = geo::origin<LocalPoint_t>();
46  fCenter = toWorldCoords(lp);
47 
48  lp.SetZ(fHalfL);
49  auto end = toWorldCoords(lp);
50 
51  fThetaZ = std::acos((end.Z() - fCenter.Z())/fHalfL);
52 
53  // check to see if it runs "forward" or "backwards" in z
54  // check is made looking at the y position of the end point
55  // relative to the center point because you want to know if
56  // the end point is above or below the center of the wire in
57  // the yz plane
58  if(end.Y() < fCenter.Y()) fThetaZ *= -1.;
59 
60  //This ensures we are looking at the angle between 0 and Pi
61  //as if the wire runs at one angle it also runs at that angle +-Pi
62  if(fThetaZ < 0) fThetaZ += util::pi();
63 
64  } // geo::WireGeo::WireGeo()
double fThetaZ
angle of the wire with respect to the z direction
Definition: WireGeo.h:333
geo::Point_t fCenter
Center of the wire in world coordinates.
Definition: WireGeo.h:335
geo::Point_t toWorldCoords(LocalPoint_t const &local) const
Transform point from local wire frame to world frame.
Definition: WireGeo.h:273
double fHalfL
half length of the wire
Definition: WireGeo.h:334
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:336
constexpr T pi()
Returns the constant pi (up to 35 decimal digits of precision)
std::vector< evd::details::RawDigitInfo_t >::const_iterator end(RawDigitCacheDataClass const &cache)
const TGeoNode * fWireNode
Pointer to the wire node.
Definition: WireGeo.h:332
bool flipped
whether start and end are reversed
Definition: WireGeo.h:337

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 347 of file WireGeo.h.

References HalfL(), max, and min.

Referenced by capRelLength(), and GetPositionFromCenter().

348  { return std::min(+HalfL(), std::max(-HalfL(), local)); }
Int_t max
Definition: plot.C:27
double HalfL() const
Returns half the length of the wire [cm].
Definition: WireGeo.h:106
Int_t min
Definition: plot.C:26
double geo::WireGeo::capRelLength ( double  local) const
inlineprivate

Stacked capLength() and relLength().

Definition at line 351 of file WireGeo.h.

References capLength(), Flip(), and relLength().

352  { return capLength(relLength(local)); }
double relLength(double local) const
Returns the relative length from center to be used when transforming.
Definition: WireGeo.h:344
double capLength(double local) const
Caps the specified local length coordinate to lay on the wire.
Definition: WireGeo.h:347
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 307 of file WireGeo.h.

References DistanceFrom(), fCenter, TanThetaZ(), and UpdateAfterSorting().

308  { return fCenter.Z() - fCenter.Y() / TanThetaZ(); }
geo::Point_t fCenter
Center of the wire in world coordinates.
Definition: WireGeo.h:335
double TanThetaZ() const
Returns angle of wire with respect to z axis in the Y-Z plane in radians.
Definition: WireGeo.h:205
double geo::WireGeo::CosThetaZ ( ) const
inline

Returns trigonometric operations on ThetaZ()

Definition at line 203 of file WireGeo.h.

Referenced by geo::ChannelMapStandardAlg::Initialize(), isVertical(), and ThetaZ().

203 { 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:192
template<typename Vector >
Vector geo::WireGeo::Direction ( ) const

Returns the wire direction as a norm-one vector.

Definition at line 377 of file WireGeo.h.

References Length().

Referenced by lar_cluster3d::MinSpanTreeAlg::configure(), DistanceFrom(), lar_cluster3d::PrincipalComponentsAlg::getHit2DPocaToAxis(), geo::PlaneGeo::GetNormalAxis(), geo::GeometryCore::IntersectSegments(), isParallelTo(), lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_2D(), lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_calc2DDocas(), and geo::PlaneGeo::shouldFlipWire().

377  {
378  // maybe (GetCenter() - GetStart()) / HalfL() would be faster;
379  // strangely, TVector3 does not implement operator/ (double).
380  return geo::vect::convertTo<Vector>(GetEnd<geo::Point_t>() - GetStart<geo::Point_t>()) * (1.0 / Length());
381 } // geo::WireGeo::Direction()
double Length() const
Returns the wire length in centimeters.
Definition: WireGeo.h:184
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 100 of file WireGeo.cxx.

References Direction(), GetCenter(), and isParallelTo().

Referenced by ComputeZatY0(), and WirePitch().

100  {
101  //
102  // The algorithm assumes that picking any point on the wire will do,
103  // that is, that the wires are parallel.
104  //
105 
106  if (!isParallelTo(wire)) return 0;
107 
108  // a vector connecting to the other wire
109  auto const toWire = wire.GetCenter() - GetCenter();
110 
111  // the distance is that vector, times the sine of the angle with the wire
112  // direction; we get that in a generic way with a cross product.
113  // We don't even care about the sign here
114  // (if we did, we would do a dot-product with the normal to the plane,
115  // and we should get a positive distance if the other wire has larger wire
116  // coordinate than this one).
117  return toWire.Cross(Direction()).Mag();
118 
119  } // WireGeo::DistanceFrom()
Point GetCenter() const
Returns the world coordinate of the center of the wire [cm].
Definition: WireGeo.h:171
Vector Direction() const
Returns the wire direction as a norm-one vector.
Definition: WireGeo.h:377
bool isParallelTo(geo::WireGeo const &wire) const
Returns if this wire is parallel to another.
Definition: WireGeo.h:215
void geo::WireGeo::Flip ( )
private

Set to swap the start and end wire.

Definition at line 131 of file WireGeo.cxx.

References flipped.

Referenced by capRelLength(), and UpdateAfterSorting().

131  {
132  // we don't need to do much to flip so far:
133  // - ThetaZ() is defined in [0, pi], invariant to flipping
134  // - we don't change the transformation matrices, that we want to be
135  // untouched and coherent with the original geometry source
136  // - center is invariant for flipping
137  // - start and end are computed on the fly (taking flipping into account)
138  // - ... and we chose to leave half length unsigned and independent
139 
140  // change the flipping bit
141  flipped = !flipped;
142 
143  } // WireGeo::Flip()
bool flipped
whether start and end are reversed
Definition: WireGeo.h:337
static double geo::WireGeo::gausSum ( double  a,
double  b 
)
inlinestaticprivate

Definition at line 358 of file WireGeo.h.

358 { return std::sqrt(a*a + b*b); }
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 68 of file WireGeo.cxx.

References fCenter, fHalfL, LocalToWorld(), and relLength().

Referenced by hit::MagDriftAna::analyze(), trkf::SpacePointAlg::compatible(), lar_pandora::LArPandoraInput::CreatePandoraHits2D(), lar_pandora::LArPandoraInput::CreatePandoraReadoutGaps(), hit::DisambigCheater::DisambigCheater(), DistanceFrom(), trkf::SpacePointAlg::fillComplexSpacePoint(), trkf::SpacePointAlg::fillSpacePoint(), lar_cluster3d::PrincipalComponentsAlg::getHit2DPocaToAxis(), geo::PlaneGeo::GetNormalAxis(), cluster::MergeClusterAlg::GlobalWire(), cluster::BlurredClusteringAlg::GlobalWire(), shower::EMShowerAlg::GlobalWire(), geo::ChannelMapStandardAlg::Initialize(), geo::GeometryCore::IntersectSegments(), recob::tracking::makePlane(), trkf::Track3DKalmanHitAlg::makeSeed(), trkf::SpacePointAlg::makeSpacePoints(), lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_2D(), lar_cluster3d::PrincipalComponentsAlg::PCAAnalysis_calc2DDocas(), geo::PlaneGeo::PlaneCoordinateFrom(), trkf::SpacePointAlg::separation(), geo::sortWireStandard(), trkf::SurfWireLine::SurfWireLine(), trkf::SurfWireX::SurfWireX(), and apa::APAGeometryAlg::ThreeChanPos().

69  {
70  if (localz == 0.) { // if no dislocation is requested, we already have it
71  xyz[0] = fCenter.X();
72  xyz[1] = fCenter.Y();
73  xyz[2] = fCenter.Z();
74  return;
75  }
76 
77  double locz = relLength(localz);
78  if (std::abs(locz) > fHalfL) {
79  mf::LogWarning("WireGeo") << "asked for position along wire that"
80  " extends beyond the wire, returning position at end point";
81  locz = relLength((locz < 0)? -fHalfL: fHalfL);
82  }
83  const double local[3] = { 0., 0., locz };
84  LocalToWorld(local, xyz);
85  }
geo::Point_t fCenter
Center of the wire in world coordinates.
Definition: WireGeo.h:335
double fHalfL
half length of the wire
Definition: WireGeo.h:334
void LocalToWorld(const double *wire, double *world) const
Transform point from local wire frame to world frame.
Definition: WireGeo.h:269
double relLength(double local) const
Returns the relative length from center to be used when transforming.
Definition: WireGeo.h:344
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
template<typename Point = DefaultPoint_t>
Point geo::WireGeo::GetCenter ( ) const
inline

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

Definition at line 171 of file WireGeo.h.

References fCenter.

Referenced by DistanceFrom(), GetEnd(), GetStart(), and HalfL().

171 { return geo::vect::convertTo<Point>(fCenter); }
geo::Point_t fCenter
Center of the wire in world coordinates.
Definition: WireGeo.h:335
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 133 of file WireGeo.h.

References fHalfL, and GetCenter().

Referenced by geo::PlaneGeo::Coverage(), evd_tool::ICARUSDrawer::DrawBadChannels(), evd_tool::MicroBooNEDrawer::DrawBadChannels(), lar_pandora::PFParticleHitDumper::GetUVW(), geo::details::ActiveAreaCalculator::initializeWireEnds(), and geo::GeometryCore::WireEndPoints().

133 { GetCenter(xyz, +fHalfL); }
Point GetCenter() const
Returns the world coordinate of the center of the wire [cm].
Definition: WireGeo.h:171
double fHalfL
half length of the wire
Definition: WireGeo.h:334
template<typename Point = DefaultPoint_t>
Point geo::WireGeo::GetEnd ( ) const
inline

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

Definition at line 180 of file WireGeo.h.

References HalfL().

181  { return GetPositionFromCenterUnbounded<Point>(+HalfL()); }
double HalfL() const
Returns half the length of the wire [cm].
Definition: WireGeo.h:106
template<typename Point = DefaultPoint_t>
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 150 of file WireGeo.h.

References capLength(), and GetPositionFromCenterUnbounded().

151  { return GetPositionFromCenterUnbounded<Point>(capLength(localz)); }
double capLength(double local) const
Caps the specified local length coordinate to lay on the wire.
Definition: WireGeo.h:347
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 369 of file WireGeo.h.

References relLength(), and toWorldCoords().

Referenced by GetPositionFromCenter().

369  {
370  return geo::vect::convertTo<Point>
371  (toWorldCoords(LocalPoint_t{ 0.0, 0.0, relLength(localz) }));
372 } // geo::WireGeo::GetPositionFromCenterImpl()
geo::Point_t toWorldCoords(LocalPoint_t const &local) const
Transform point from local wire frame to world frame.
Definition: WireGeo.h:273
double relLength(double local) const
Returns the relative length from center to be used when transforming.
Definition: WireGeo.h:344
geo::Point3DBase_t< WireGeoCoordinatesTag > LocalPoint_t
Type of points in the local GDML wire plane frame.
Definition: WireGeo.h:89
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 129 of file WireGeo.h.

References fHalfL, and GetCenter().

Referenced by geo::PlaneGeo::Coverage(), evd_tool::ICARUSDrawer::DrawBadChannels(), evd_tool::MicroBooNEDrawer::DrawBadChannels(), lar_pandora::PFParticleHitDumper::GetUVW(), geo::details::ActiveAreaCalculator::initializeWireEnds(), and geo::GeometryCore::WireEndPoints().

129 { GetCenter(xyz, -fHalfL); }
Point GetCenter() const
Returns the world coordinate of the center of the wire [cm].
Definition: WireGeo.h:171
double fHalfL
half length of the wire
Definition: WireGeo.h:334
template<typename Point = DefaultPoint_t>
Point geo::WireGeo::GetStart ( ) const
inline

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

Definition at line 175 of file WireGeo.h.

References HalfL().

176  { return GetPositionFromCenterUnbounded<Point>(-HalfL()); }
double HalfL() const
Returns half the length of the wire [cm].
Definition: WireGeo.h:106
bool geo::WireGeo::isFlipped ( ) const
inlineprivate

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

Definition at line 341 of file WireGeo.h.

References flipped.

Referenced by relLength().

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

Returns if this wire is horizontal (theta_z ~ 0)

Definition at line 209 of file WireGeo.h.

References e, and SinThetaZ().

Referenced by PrintWireInfo().

209 { return std::abs(SinThetaZ()) < 1e-5; }
double SinThetaZ() const
Returns angle of wire with respect to z axis in the Y-Z plane in radians.
Definition: WireGeo.h:204
Float_t e
Definition: plot.C:34
bool geo::WireGeo::isParallelTo ( geo::WireGeo const &  wire) const
inline

Returns if this wire is parallel to another.

Definition at line 215 of file WireGeo.h.

References Direction(), e, art::detail::indent(), and PrintWireInfo().

Referenced by DistanceFrom().

216  {
217  return // parallel if the dot product of the directions is about +/- 1
218  std::abs(std::abs(Direction<geo::Vector_t>().Dot(wire.Direction<geo::Vector_t>())) - 1.) < 1e-5;
219  }
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
Float_t e
Definition: plot.C:34
bool geo::WireGeo::isVertical ( ) const
inline

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

Definition at line 212 of file WireGeo.h.

References CosThetaZ(), and e.

Referenced by PrintWireInfo().

212 { return std::abs(CosThetaZ()) < 1e-5; }
double CosThetaZ() const
Returns trigonometric operations on ThetaZ()
Definition: WireGeo.h:203
Float_t e
Definition: plot.C:34
double geo::WireGeo::Length ( ) const
inline

Returns the wire length in centimeters.

Definition at line 184 of file WireGeo.h.

References HalfL().

Referenced by Direction(), and PrintWireInfo().

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

Transform point from local wire frame to world frame.

Definition at line 269 of file WireGeo.h.

Referenced by GetCenter().

270  { 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:336
void geo::WireGeo::LocalToWorldVect ( const double *  wire,
double *  world 
) const
inline

Transform direction vector from local to world.

Definition at line 277 of file WireGeo.h.

References fTrans, and geo::LocalTransformation< StoredMatrix >::LocalToWorldVect().

278  { 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:336
const TGeoNode* geo::WireGeo::Node ( ) const
inline

Definition at line 303 of file WireGeo.h.

References fWireNode.

303 { return fWireNode; }
const TGeoNode * fWireNode
Pointer to the wire node.
Definition: WireGeo.h:332
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 386 of file WireGeo.h.

References isHorizontal(), isVertical(), Length(), and ThetaZ().

Referenced by isParallelTo(), and geo::GeometryCore::Print().

390  {
391 
392  //----------------------------------------------------------------------------
393  out << "wire from " << GetStart<geo::Point_t>()
394  << " to " << GetEnd<geo::Point_t>();
395 
396  if (verbosity-- <= 0) return; // 0
397 
398  //----------------------------------------------------------------------------
399  out << " (" << Length() << " cm long)";
400 
401  if (verbosity-- <= 0) return; // 1
402 
403  //----------------------------------------------------------------------------
404  out << ", theta(z)=" << ThetaZ() << " rad";
405 
406  if (verbosity-- <= 0) return; // 2
407 
408  //----------------------------------------------------------------------------
409  out << "\n" << indent
410  << " center at " << GetCenter<geo::Point_t>() << " cm";
411 
412  if (verbosity-- <= 0) return; // 3
413 
414  //----------------------------------------------------------------------------
415  out << ", direction: " << Direction<geo::Vector_t>();
416  if (isHorizontal()) out << " (horizontal)";
417  if (isVertical()) out << " (vertical)";
418 
419 // if (verbosity-- <= 0) return; // 4
420 
421  //----------------------------------------------------------------------------
422 } // geo::WireGeo::PrintWireInfo()
bool isVertical() const
Returns if this wire is vertical (theta_z ~ pi/2)
Definition: WireGeo.h:212
double Length() const
Returns the wire length in centimeters.
Definition: WireGeo.h:184
double ThetaZ() const
Returns angle of wire with respect to z axis in the Y-Z plane in radians.
Definition: WireGeo.h:192
bool isHorizontal() const
Returns if this wire is horizontal (theta_z ~ 0)
Definition: WireGeo.h:209
std::string indent(std::size_t const i)
double geo::WireGeo::relLength ( double  local) const
inlineprivate

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

Definition at line 344 of file WireGeo.h.

References isFlipped().

Referenced by capRelLength(), GetCenter(), and GetPositionFromCenterUnbounded().

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

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

Definition at line 88 of file WireGeo.cxx.

References fWireNode.

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

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

Definition at line 92 of file WireGeo.cxx.

References fWireNode.

Referenced by HalfL().

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

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

Definition at line 204 of file WireGeo.h.

References ThetaZ().

Referenced by geo::ChannelMapStandardAlg::Initialize(), and isHorizontal().

204 { 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:192
double geo::WireGeo::TanThetaZ ( ) const
inline

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

Definition at line 205 of file WireGeo.h.

References ThetaZ().

Referenced by ComputeZatY0().

205 { 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:192
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 96 of file WireGeo.cxx.

References fThetaZ, and util::RadiansToDegrees().

97  { return degrees? util::RadiansToDegrees(fThetaZ): fThetaZ; }
double fThetaZ
angle of the wire with respect to the z direction
Definition: WireGeo.h:333
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 289 of file WireGeo.h.

References fTrans, and geo::LocalTransformationGeo< StoredMatrix, LocalPoint, LocalVector >::toLocalCoords().

290  { 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:336
LocalVector_t geo::WireGeo::toLocalCoords ( geo::Vector_t const &  world) const
inline

Transform direction vector from world to local.

Definition at line 297 of file WireGeo.h.

References fTrans, and geo::LocalTransformationGeo< StoredMatrix, LocalPoint, LocalVector >::toLocalCoords().

298  { 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:336
geo::Point_t geo::WireGeo::toWorldCoords ( LocalPoint_t const &  local) const
inline

Transform point from local wire frame to world frame.

Definition at line 273 of file WireGeo.h.

References fTrans, and geo::LocalTransformationGeo< StoredMatrix, LocalPoint, LocalVector >::toWorldCoords().

Referenced by GetPositionFromCenterUnbounded(), and WireGeo().

274  { return fTrans.toWorldCoords(local); }
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:336
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 281 of file WireGeo.h.

References fTrans, and geo::LocalTransformationGeo< StoredMatrix, LocalPoint, LocalVector >::toWorldCoords().

282  { return fTrans.toWorldCoords(local); }
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:336
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 123 of file WireGeo.cxx.

References Flip().

Referenced by ComputeZatY0().

123  {
124 
125  // flip, if asked
126  if (flip) Flip();
127 
128  } // WireGeo::UpdateAfterSorting()
void Flip()
Set to swap the start and end wire.
Definition: WireGeo.cxx:131
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 325 of file WireGeo.h.

References DistanceFrom().

Referenced by geo::PlaneGeo::UpdateWirePitch(), and geo::PlaneGeo::UpdateWirePitchSlow().

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

Transform point from world frame to local wire frame.

Definition at line 285 of file WireGeo.h.

References fTrans, and geo::LocalTransformation< StoredMatrix >::WorldToLocal().

286  { 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:336
void geo::WireGeo::WorldToLocalVect ( const double *  world,
double *  wire 
) const
inline

Transform direction vector from world to local.

Definition at line 293 of file WireGeo.h.

References fTrans, and geo::LocalTransformation< StoredMatrix >::WorldToLocalVect().

294  { fTrans.WorldToLocalVect(world, wire); }
LocalTransformation_t fTrans
Wire to world transform.
Definition: WireGeo.h:336
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 335 of file WireGeo.h.

Referenced by ComputeZatY0(), GetCenter(), and WireGeo().

double geo::WireGeo::fHalfL
private

half length of the wire

Definition at line 334 of file WireGeo.h.

Referenced by GetCenter(), GetEnd(), GetStart(), HalfL(), and WireGeo().

bool geo::WireGeo::flipped
private

whether start and end are reversed

Definition at line 337 of file WireGeo.h.

Referenced by Flip(), and isFlipped().

double geo::WireGeo::fThetaZ
private

angle of the wire with respect to the z direction

Definition at line 333 of file WireGeo.h.

Referenced by ThetaZ(), and WireGeo().

LocalTransformation_t geo::WireGeo::fTrans
private

Wire to world transform.

Definition at line 336 of file WireGeo.h.

Referenced by LocalToWorldVect(), toLocalCoords(), toWorldCoords(), WorldToLocal(), and WorldToLocalVect().

const TGeoNode* geo::WireGeo::fWireNode
private

Pointer to the wire node.

Definition at line 332 of file WireGeo.h.

Referenced by Node(), RMax(), RMin(), and WireGeo().

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

Maximum verbosity supported by PrintWireInfo().

Definition at line 251 of file WireGeo.h.

Referenced by geo::GeometryCore::Print().


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