EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Acts::ConeSurface Class Reference

#include <acts/blob/sPHENIX/Core/include/Acts/Surfaces/ConeSurface.hpp>

+ Inheritance diagram for Acts::ConeSurface:
+ Collaboration diagram for Acts::ConeSurface:

Public Member Functions

 ~ConeSurface () override=default
 
 ConeSurface ()=delete
 
ConeSurfaceoperator= (const ConeSurface &other)
 
Vector3D binningPosition (const GeometryContext &gctx, BinningValue bValue) const final
 
SurfaceType type () const override
 Return the surface type.
 
RotationMatrix3D referenceFrame (const GeometryContext &gctx, const Vector3D &position, const Vector3D &momentum) const final
 
Vector3D normal (const GeometryContext &gctx, const Vector2D &lposition) const final
 
Vector3D normal (const GeometryContext &gctx, const Vector3D &position) const final
 
virtual Vector3D rotSymmetryAxis (const GeometryContext &gctx) const
 
const ConeBoundsbounds () const final
 This method returns the ConeBounds by reference.
 
Vector3D localToGlobal (const GeometryContext &gctx, const Vector2D &lposition, const Vector3D &momentum) const final
 
Result< Vector2DglobalToLocal (const GeometryContext &gctx, const Vector3D &position, const Vector3D &momentum) const final
 
SurfaceIntersection intersect (const GeometryContext &gctx, const Vector3D &position, const Vector3D &direction, const BoundaryCheck &bcheck) const final
 
double pathCorrection (const GeometryContext &gctx, const Vector3D &position, const Vector3D &direction) const final
 
Polyhedron polyhedronRepresentation (const GeometryContext &gctx, size_t lseg) const override
 
std::string name () const override
 Return properly formatted class name for screen output.
 
LocalCartesianToBoundLocalMatrix localCartesianToBoundLocalDerivative (const GeometryContext &gctx, const Vector3D &position) const final
 
- Public Member Functions inherited from Acts::Surface
virtual ~Surface ()
 
std::shared_ptr< SurfacegetSharedPtr ()
 
std::shared_ptr< const SurfacegetSharedPtr () const
 
Surfaceoperator= (const Surface &other)
 
virtual bool operator== (const Surface &other) const
 
virtual bool operator!= (const Surface &sf) const
 
virtual const Transform3Dtransform (const GeometryContext &gctx) const
 
virtual Vector3D center (const GeometryContext &gctx) const
 
virtual Vector3D normal (const GeometryContext &gctx) const
 
const DetectorElementBaseassociatedDetectorElement () const
 
const LayerassociatedLayer () const
 
void associateLayer (const Layer &lay)
 
const ISurfaceMaterialsurfaceMaterial () const
 
const std::shared_ptr< const
ISurfaceMaterial > & 
surfaceMaterialSharedPtr () const
 
void assignSurfaceMaterial (std::shared_ptr< const ISurfaceMaterial > material)
 
bool isOnSurface (const GeometryContext &gctx, const Vector3D &position, const Vector3D &momentum, const BoundaryCheck &bcheck=true) const
 
virtual bool insideBounds (const Vector2D &lposition, const BoundaryCheck &bcheck=true) const
 
virtual void initJacobianToGlobal (const GeometryContext &gctx, BoundToFreeMatrix &jacobian, const Vector3D &position, const Vector3D &direction, const BoundVector &pars) const
 
virtual RotationMatrix3D initJacobianToLocal (const GeometryContext &gctx, FreeToBoundMatrix &jacobian, const Vector3D &position, const Vector3D &direction) const
 
virtual BoundRowVector derivativeFactors (const GeometryContext &gctx, const Vector3D &position, const Vector3D &direction, const RotationMatrix3D &rft, const BoundToFreeMatrix &jacobian) const
 
virtual std::ostream & toStream (const GeometryContext &gctx, std::ostream &sl) const
 
AlignmentToBoundMatrix alignmentToBoundDerivative (const GeometryContext &gctx, const FreeVector &derivatives, const Vector3D &position, const Vector3D &direction) const
 
virtual AlignmentRowVector alignmentToPathDerivative (const GeometryContext &gctx, const RotationMatrix3D &rotToLocalZAxis, const Vector3D &position, const Vector3D &direction) const
 
- Public Member Functions inherited from Acts::GeometryObject
 GeometryObject ()=default
 Defaulted construrctor.
 
 GeometryObject (const GeometryObject &)=default
 Defaulted copy constructor.
 
 GeometryObject (const GeometryIdentifier &geometryId)
 
GeometryObjectoperator= (const GeometryObject &geometryId)
 
const GeometryIdentifiergeometryId () const
 
virtual double binningPositionValue (const GeometryContext &gctx, BinningValue bValue) const
 
void assignGeometryId (const GeometryIdentifier &geometryId)
 

Protected Member Functions

 ConeSurface (const Transform3D &transform, double alpha, bool symmetric=false)
 
 ConeSurface (const Transform3D &transform, double alpha, double zmin, double zmax, double halfPhi=M_PI)
 
 ConeSurface (const Transform3D &transform, const std::shared_ptr< const ConeBounds > &cbounds)
 
 ConeSurface (const ConeSurface &other)
 
 ConeSurface (const GeometryContext &gctx, const ConeSurface &other, const Transform3D &shift)
 
- Protected Member Functions inherited from Acts::Surface
 Surface (const Transform3D &transform=Transform3D::Identity())
 
 Surface (const Surface &other)
 
 Surface (const DetectorElementBase &detelement)
 
 Surface (const GeometryContext &gctx, const Surface &other, const Transform3D &shift)
 

Protected Attributes

std::shared_ptr< const ConeBoundsm_bounds
 bounds (shared)
 
- Protected Attributes inherited from Acts::Surface
Transform3D m_transform = Transform3D::Identity()
 
const DetectorElementBasem_associatedDetElement {nullptr}
 Pointer to the a DetectorElementBase.
 
const Layerm_associatedLayer {nullptr}
 
const TrackingVolumem_associatedTrackingVolume {nullptr}
 
std::shared_ptr< const
ISurfaceMaterial
m_surfaceMaterial
 Possibility to attach a material descrption.
 
- Protected Attributes inherited from Acts::GeometryObject
GeometryIdentifier m_geometryId
 

Private Member Functions

detail::RealQuadraticEquation intersectionSolver (const GeometryContext &gctx, const Vector3D &position, const Vector3D &direction) const
 

Private Attributes

friend Surface
 

Additional Inherited Members

- Public Types inherited from Acts::Surface
enum  SurfaceType {
  Cone = 0, Cylinder = 1, Disc = 2, Perigee = 3,
  Plane = 4, Straw = 5, Curvilinear = 6, Other = 7
}
 
- Static Public Member Functions inherited from Acts::Surface
template<class T , typename... Args>
static std::shared_ptr< TmakeShared (Args &&...args)
 

Detailed Description

Class for a conical surface in the Tracking geometry. It inherits from Surface.

The ConeSurface is special since no corresponding Track parameters exist since they're numerical instable at the tip of the cone. Propagations to a cone surface will be returned in curvilinear coordinates.

Definition at line 32 of file ConeSurface.hpp.

View newest version in sPHENIX GitHub at line 32 of file ConeSurface.hpp

Constructor & Destructor Documentation

Acts::ConeSurface::ConeSurface ( const Transform3D transform,
double  alpha,
bool  symmetric = false 
)
protected

Constructor form HepTransform and an opening angle

Parameters
transformis the transform to place to cone in a 3D frame
alphais the opening angle of the cone
symmetricindicates if the cones are built to +/1 z

Definition at line 31 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 31 of file ConeSurface.cpp

Acts::ConeSurface::ConeSurface ( const Transform3D transform,
double  alpha,
double  zmin,
double  zmax,
double  halfPhi = M_PI 
)
protected

Constructor form HepTransform and an opening angle

Parameters
transformis the transform that places the cone in the global frame
alphais the opening angle of the cone
zminis the z range over which the cone spans
zmaxis the z range over which the cone spans
halfPhiis the openen angle for cone ssectors

Definition at line 37 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 37 of file ConeSurface.cpp

Acts::ConeSurface::ConeSurface ( const Transform3D transform,
const std::shared_ptr< const ConeBounds > &  cbounds 
)
protected

Constructor from HepTransform and ConeBounds

Parameters
transformis the transform that places the cone in the global frame
cboundsis the boundary class, the bounds must exit

Definition at line 44 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 44 of file ConeSurface.cpp

References throw_assert.

Acts::ConeSurface::ConeSurface ( const ConeSurface other)
protected

Copy constructor

Parameters
otheris the source cone surface

Definition at line 23 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 23 of file ConeSurface.cpp

Acts::ConeSurface::ConeSurface ( const GeometryContext gctx,
const ConeSurface other,
const Transform3D shift 
)
protected

Copy constructor - with shift

Parameters
gctxThe current geometry context object, e.g. alignment
otheris the source cone surface
shiftis the additional transfrom applied after copying

Definition at line 26 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 26 of file ConeSurface.cpp

Acts::ConeSurface::~ConeSurface ( )
overridedefault
Acts::ConeSurface::ConeSurface ( )
delete

Member Function Documentation

Acts::Vector3D Acts::ConeSurface::binningPosition ( const GeometryContext gctx,
Acts::BinningValue  bValue 
) const
finalvirtual

The binning position method - is overloaded for r-type binning

Parameters
gctxThe current geometry context object, e.g. alignment
bValuedefines the type of binning applied in the global frame
Returns
The return type is a vector for positioning in the global frame

Implements Acts::GeometryObject.

Definition at line 50 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 50 of file ConeSurface.cpp

References Acts::binR, and Acts::binRPhi.

const Acts::ConeBounds & Acts::ConeSurface::bounds ( ) const
finalvirtual

This method returns the ConeBounds by reference.

Implements Acts::Surface.

Definition at line 166 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 166 of file ConeSurface.cpp

Acts::Result< Acts::Vector2D > Acts::ConeSurface::globalToLocal ( const GeometryContext gctx,
const Vector3D position,
const Vector3D momentum 
) const
finalvirtual

Global to local transformation

Parameters
gctxThe current geometry context object, e.g. alignment
positionis the global position to be transformed
momentumis the global momentum (ignored in this operation)
Returns
a Result<Vector2D> which can be !ok() if the operation fails

Implements Acts::Surface.

Definition at line 112 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 112 of file ConeSurface.cpp

References kdfinder::abs(), Acts::VectorHelpers::perp(), Acts::VectorHelpers::position(), Acts::s_onSurfaceTolerance, and Acts::Test::transform.

+ Here is the call graph for this function:

SurfaceIntersection ConeSurface::intersect ( const GeometryContext gctx,
const Vector3D position,
const Vector3D direction,
const BoundaryCheck bcheck 
) const
inlinefinalvirtual

Straight line intersection schema from position/direction

Parameters
gctxThe current geometry context object, e.g. alignment
positionThe position to start from
directionThe direction at start
bcheckthe Boundary Check

If possible returns both solutions for the cylinder

Returns
SurfaceIntersection object (contains intersection & surface)

Implements Acts::Surface.

Definition at line 32 of file ConeSurface.ipp.

View newest version in sPHENIX GitHub at line 32 of file ConeSurface.ipp

References Acts::ObjectIntersection< object_t, representation_t >::alternative, Acts::ObjectIntersection< object_t, representation_t >::intersection, Acts::s_onSurfaceTolerance, and Acts::Test::transform.

detail::RealQuadraticEquation ConeSurface::intersectionSolver ( const GeometryContext gctx,
const Vector3D position,
const Vector3D direction 
) const
inlineprivate

Implementation of the intersection solver

mathematical motivation:

The calculation will be done in the 3-dim frame of the cone, i.e. the symmetry axis of the cone is the z-axis, x- and y-axis are perpendicular to the the z-axis. In this frame the cone is centered around the origin. Therefore the two points describing the line have to be first recalculated into the new frame. Suppose, this is done, the points of intersection can be obtained as follows:

The cone is described by the implicit equation $x^2 + y^2 = z^2 \tan \alpha$ where $\alpha$ is opening half-angle of the cone the and the line by the parameter equation (with $t$ the parameter and $x_1$ and $x_2$ are points on the line) $(x,y,z) = \vec x_1 + (\vec x_2 - \vec x_2) t $. The intersection is the given to the value of $t$ where the $(x,y,z)$ coordinates of the line satisfy the implicit equation of the cone. Inserting the expression for the points on the line into the equation of the cone and rearranging to the form of a gives (letting $ \vec x_d = \frac{\vec x_2 - \vec x_1}{|\vec x_2 - \vec x_1|} $): $t^2 (x_d^2 + y_d^2 - z_d^2 \tan^2 \alpha) + 2 t (x_1 x_d + y_1 y_d - z_1 z_d \tan^2 \alpha) + (x_1^2 + y_1^2 - z_1^2 \tan^2 \alpha) = 0 $ Solving the above for $t$ and putting the values into the equation of the line gives the points of intersection. $t$ is also the length of the path, since we normalized $x_d$ to be unit length.

Returns
the quadratic equation

Definition at line 9 of file ConeSurface.ipp.

View newest version in sPHENIX GitHub at line 9 of file ConeSurface.ipp

References Acts::UnitConstants::C, Acts::UnitConstants::e, Acts::VectorHelpers::position(), and Acts::Test::transform.

+ Here is the call graph for this function:

LocalCartesianToBoundLocalMatrix ConeSurface::localCartesianToBoundLocalDerivative ( const GeometryContext gctx,
const Vector3D position 
) const
inlinefinalvirtual

Calculate the derivative of bound track parameters local position w.r.t. position in local 3D Cartesian coordinates

Parameters
gctxThe current geometry context object, e.g. alignment
positionThe position of the paramters in global
Returns
Derivative of bound local position w.r.t. position in local 3D cartesian coordinates

Implements Acts::Surface.

Definition at line 91 of file ConeSurface.ipp.

View newest version in sPHENIX GitHub at line 91 of file ConeSurface.ipp

References cos(), Acts::VectorHelpers::perp(), phi, Acts::VectorHelpers::position(), Acts::IntegrationTest::R, and Acts::Test::transform.

+ Here is the call graph for this function:

Acts::Vector3D Acts::ConeSurface::localToGlobal ( const GeometryContext gctx,
const Vector2D lposition,
const Vector3D momentum 
) const
finalvirtual

Local to global transformation

Parameters
gctxThe current geometry context object, e.g. alignment
lpositionis the local position to be transformed
momentumis the global momentum (ignored in this operation)
Returns
The global position by value

Implements Acts::Surface.

Definition at line 102 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 102 of file ConeSurface.cpp

References cos(), Acts::eBoundLoc0, Acts::eBoundLoc1, phi, and Acts::Test::transform.

+ Here is the call graph for this function:

std::string Acts::ConeSurface::name ( ) const
overridevirtual

Return properly formatted class name for screen output.

Implements Acts::Surface.

Definition at line 140 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 140 of file ConeSurface.cpp

Referenced by eicpy.verify.PythiaHistograms::__init__(), and eicpy.verify.DjangohHistograms::__init__().

+ Here is the caller graph for this function:

Acts::Vector3D Acts::ConeSurface::normal ( const GeometryContext gctx,
const Vector2D lposition 
) const
finalvirtual

Return method for surface normal information

Parameters
gctxThe current geometry context object, e.g. alignment
lpositionis the local position at normal vector request
Returns
Vector3D normal vector in global frame

Implements Acts::Surface.

Definition at line 144 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 144 of file ConeSurface.cpp

References cos(), Acts::ConeBounds::eAlpha, Acts::eBoundLoc0, Acts::eBoundLoc1, phi, and Acts::Test::transform.

+ Here is the call graph for this function:

Acts::Vector3D Acts::ConeSurface::normal ( const GeometryContext gctx,
const Vector3D position 
) const
finalvirtual

Return method for surface normal information

Parameters
gctxThe current geometry context object, e.g. alignment
positionis the global position as normal vector base
Returns
Vector3D normal vector in global frame

Reimplemented from Acts::Surface.

Definition at line 157 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 157 of file ConeSurface.cpp

References Acts::VectorHelpers::position(), and Acts::Test::transform.

+ Here is the call graph for this function:

Acts::ConeSurface & Acts::ConeSurface::operator= ( const ConeSurface other)

Assignment operator

Parameters
otheris the source surface for the assignment

Definition at line 68 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 68 of file ConeSurface.cpp

References m_bounds, and Acts::Surface::operator=().

+ Here is the call graph for this function:

double Acts::ConeSurface::pathCorrection ( const GeometryContext gctx,
const Vector3D position,
const Vector3D direction 
) const
finalvirtual

The pathCorrection for derived classes with thickness

Parameters
gctxThe current geometry context object, e.g. alignment
positionis the global potion at the correction point
directionis the momentum direction at the correction point
Returns
is the path correction due to incident angle

Implements Acts::Surface.

Definition at line 124 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 124 of file ConeSurface.cpp

References kdfinder::abs(), cos(), Acts::ConeBounds::eAlpha, Acts::VectorHelpers::phi(), phi, Acts::VectorHelpers::position(), and Acts::Test::transform.

+ Here is the call graph for this function:

Acts::Polyhedron Acts::ConeSurface::polyhedronRepresentation ( const GeometryContext gctx,
size_t  lseg 
) const
overridevirtual

Return a Polyhedron for the surfaces

Parameters
gctxThe current geometry context object, e.g. alignment
lsegNumber of segments along curved lines, it represents the full 2*M_PI coverange, if lseg is set to 1 only the extrema are given
Note
that a surface transform can invalidate the extrema in the transformed space
Returns
A list of vertices and a face/facett description of it

Implements Acts::Surface.

Definition at line 171 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 171 of file ConeSurface.cpp

References kdfinder::abs(), Acts::detail::VerticesHelper::createSegment(), Acts::detail::FacesHelper::cylindricalFaceMesh(), Acts::ConeBounds::eAveragePhi, Acts::ConeBounds::eHalfPhiSector, Acts::ConeBounds::eMaxZ, Acts::ConeBounds::eMinZ, M_PI, one, Acts::detail::VerticesHelper::phiSegments(), Acts::s_onSurfaceTolerance, boost::swap(), three, Acts::Test::transform, two, and z.

+ Here is the call graph for this function:

Acts::RotationMatrix3D Acts::ConeSurface::referenceFrame ( const GeometryContext gctx,
const Vector3D position,
const Vector3D momentum 
) const
finalvirtual

Return the measurement frame - this is needed for alignment, in particular for StraightLine and Perigee Surface

  • the default implementation is the the RotationMatrix3D of the transform
Parameters
gctxThe current geometry context object, e.g. alignment
positionis the global position where the measurement frame is constructed
momentumis the momentum used for the measurement frame construction
Returns
matrix that indicates the measurement frame

<

Reimplemented from Acts::Surface.

Definition at line 81 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 81 of file ConeSurface.cpp

Acts::Vector3D Acts::ConeSurface::rotSymmetryAxis ( const GeometryContext gctx) const
virtual
Parameters
gctxThe current geometry context object, e.g. alignment

Definition at line 76 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 76 of file ConeSurface.cpp

References matrix(), and Acts::Test::transform.

+ Here is the call graph for this function:

Acts::Surface::SurfaceType Acts::ConeSurface::type ( ) const
overridevirtual

Return the surface type.

Implements Acts::Surface.

Definition at line 64 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 64 of file ConeSurface.cpp

References Acts::Surface::Cone.

Member Data Documentation

std::shared_ptr<const ConeBounds> Acts::ConeSurface::m_bounds
protected

bounds (shared)

Definition at line 211 of file ConeSurface.hpp.

View newest version in sPHENIX GitHub at line 211 of file ConeSurface.hpp

Referenced by operator=().

friend Acts::ConeSurface::Surface
private

Definition at line 33 of file ConeSurface.hpp.

View newest version in sPHENIX GitHub at line 33 of file ConeSurface.hpp


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