Analysis Software
Documentation for sPHENIX simulation 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)
 
Vector3 binningPosition (const GeometryContext &gctx, BinningValue bValue) const final
 
SurfaceType type () const override
 Return the surface type.
 
RotationMatrix3 referenceFrame (const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction) const final
 
Vector3 normal (const GeometryContext &gctx, const Vector2 &lposition) const final
 
Vector3 normal (const GeometryContext &gctx, const Vector3 &position) const final
 
virtual Vector3 rotSymmetryAxis (const GeometryContext &gctx) const
 
const ConeBoundsbounds () const final
 This method returns the ConeBounds by reference.
 
Vector3 localToGlobal (const GeometryContext &gctx, const Vector2 &lposition, const Vector3 &direction) const final
 
Result< Vector2globalToLocal (const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction, double tolerance=s_onSurfaceTolerance) const final
 
SurfaceMultiIntersection intersect (const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction, const BoundaryCheck &bcheck=false, double tolerance=s_onSurfaceTolerance) const final
 
double pathCorrection (const GeometryContext &gctx, const Vector3 &position, const Vector3 &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.
 
AlignmentToPathMatrix alignmentToPathDerivative (const GeometryContext &gctx, const FreeVector &parameters) const final
 
ActsMatrix< 2, 3 > localCartesianToBoundLocalDerivative (const GeometryContext &gctx, const Vector3 &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 Transform3transform (const GeometryContext &gctx) const
 
virtual Vector3 center (const GeometryContext &gctx) const
 
virtual Vector3 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 assignDetectorElement (const DetectorElementBase &detelement)
 
void assignSurfaceMaterial (std::shared_ptr< const ISurfaceMaterial > material)
 
bool isOnSurface (const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction, const BoundaryCheck &bcheck=true) const
 
virtual bool insideBounds (const Vector2 &lposition, const BoundaryCheck &bcheck=true) const
 
virtual BoundToFreeMatrix boundToFreeJacobian (const GeometryContext &gctx, const BoundVector &boundParams) const
 
virtual FreeToBoundMatrix freeToBoundJacobian (const GeometryContext &gctx, const FreeVector &parameters) const
 
virtual FreeToPathMatrix freeToPathDerivative (const GeometryContext &gctx, const FreeVector &parameters) const
 
virtual std::ostream & toStream (const GeometryContext &gctx, std::ostream &sl) const
 
std::string toString (const GeometryContext &gctx) const
 
AlignmentToBoundMatrix alignmentToBoundDerivative (const GeometryContext &gctx, const FreeVector &parameters, const FreeVector &pathDerivative) const
 
- Public Member Functions inherited from Acts::GeometryObject
 GeometryObject ()=default
 Defaulted constructor.
 
 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 Transform3 &transform, double alpha, bool symmetric=false)
 
 ConeSurface (const Transform3 &transform, double alpha, double zmin, double zmax, double halfPhi=M_PI)
 
 ConeSurface (const Transform3 &transform, std::shared_ptr< const ConeBounds > cbounds)
 
 ConeSurface (const ConeSurface &other)
 
 ConeSurface (const GeometryContext &gctx, const ConeSurface &other, const Transform3 &shift)
 
- Protected Member Functions inherited from Acts::Surface
 Surface (const Transform3 &transform=Transform3::Identity())
 
 Surface (const Surface &other)
 
 Surface (const DetectorElementBase &detelement)
 
 Surface (const GeometryContext &gctx, const Surface &other, const Transform3 &shift)
 

Protected Attributes

std::shared_ptr< const ConeBoundsm_bounds
 bounds (shared)
 
- Protected Attributes inherited from Acts::Surface
Transform3 m_transform = Transform3::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 description.
 
- Protected Attributes inherited from Acts::GeometryObject
GeometryIdentifier m_geometryId
 

Private Member Functions

detail::RealQuadraticEquation intersectionSolver (const GeometryContext &gctx, const Vector3 &position, const Vector3 &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)
 
- Static Public Attributes inherited from Acts::Surface
static std::array< std::string,
SurfaceType::Other > 
s_surfaceTypeNames
 Helper strings for screen output.
 

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 42 of file ConeSurface.hpp.

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

Constructor & Destructor Documentation

Acts::ConeSurface::ConeSurface ( const Transform3 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 39 of file ConeSurface.cpp.

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

Acts::ConeSurface::ConeSurface ( const Transform3 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 opening angle for cone ssectors

Definition at line 45 of file ConeSurface.cpp.

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

Acts::ConeSurface::ConeSurface ( const Transform3 transform,
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 52 of file ConeSurface.cpp.

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

References m_bounds, and throw_assert.

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

Copy constructor

Parameters
otheris the source cone surface

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 GeometryContext gctx,
const ConeSurface other,
const Transform3 shift 
)
protected

Copy constructor - with shift

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

Definition at line 34 of file ConeSurface.cpp.

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

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

Member Function Documentation

Acts::AlignmentToPathMatrix Acts::ConeSurface::alignmentToPathDerivative ( const GeometryContext gctx,
const FreeVector parameters 
) const
finalvirtual

Calculate the derivative of path length at the geometry constraint or point-of-closest-approach w.r.t. alignment parameters of the surface (i.e. local frame origin in global 3D Cartesian coordinates and its rotation represented with extrinsic Euler angles)

Parameters
gctxThe current geometry context object, e.g. alignment
parametersis the free parameters
Returns
Derivative of path length w.r.t. the alignment parameters

Reimplemented from Acts::Surface.

Definition at line 330 of file ConeSurface.cpp.

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

References bounds, dy, dz, Acts::eAlignmentCenter0, Acts::eAlignmentRotation0, Acts::eFreeDir0, Acts::eFreePos0, norm, position, Acts::detail::rotationToLocalAxesDerivative(), and Acts::Test::transform.

+ Here is the call graph for this function:

Acts::Vector3 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 58 of file ConeSurface.cpp.

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

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

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

This method returns the ConeBounds by reference.

Implements Acts::Surface.

Definition at line 173 of file ConeSurface.cpp.

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

Acts::Result< Acts::Vector2 > Acts::ConeSurface::globalToLocal ( const GeometryContext gctx,
const Vector3 position,
const Vector3 direction,
double  tolerance = s_onSurfaceTolerance 
) const
finalvirtual

Global to local transformation

Parameters
gctxThe current geometry context object, e.g. alignment
positionis the global position to be transformed
directionis the global momentum direction (ignored in this operation)
toleranceoptional tolerance within which a point is considered valid on surface
Returns
a Result<Vector2> which can be !ok() if the operation fails

Implements Acts::Surface.

Definition at line 120 of file ConeSurface.cpp.

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

References bounds, Acts::VectorHelpers::perp(), position, physmon_track_finding_ttbar::r, and Acts::Test::transform.

+ Here is the call graph for this function:

Acts::SurfaceMultiIntersection Acts::ConeSurface::intersect ( const GeometryContext gctx,
const Vector3 position,
const Vector3 direction,
const BoundaryCheck bcheck = false,
double  tolerance = s_onSurfaceTolerance 
) const
finalvirtual

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
tolerancethe tolerance used for the intersection

If possible returns both solutions for the cylinder

Returns
SurfaceMultiIntersection object (contains intersection & surface)

Implements Acts::Surface.

Definition at line 288 of file ConeSurface.cpp.

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

References Acts::Intersection< 3 >::invalid(), Acts::Intersection< DIM >::pathLength(), and Acts::Test::transform.

+ Here is the call graph for this function:

Acts::detail::RealQuadraticEquation Acts::ConeSurface::intersectionSolver ( const GeometryContext gctx,
const Vector3 position,
const Vector3 direction 
) const
private

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 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 265 of file ConeSurface.cpp.

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

References A, bounds, C, Acts::UnitConstants::e, position, and Acts::Test::transform.

Acts::ActsMatrix< 2, 3 > Acts::ConeSurface::localCartesianToBoundLocalDerivative ( const GeometryContext gctx,
const Vector3 position 
) const
finalvirtual

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 parameters in global
Returns
Derivative of bound local position w.r.t. position in local 3D cartesian coordinates

Implements Acts::Surface.

Definition at line 383 of file ConeSurface.cpp.

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

References bounds, Acts::VectorHelpers::perp(), ActsTests::PropagationDatasets::phi, Acts::VectorHelpers::phi(), position, Acts::IntegrationTest::R, sTransform, and Acts::Test::transform.

+ Here is the call graph for this function:

Acts::Vector3 Acts::ConeSurface::localToGlobal ( const GeometryContext gctx,
const Vector2 lposition,
const Vector3 direction 
) const
finalvirtual

Local to global transformation

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

Implements Acts::Surface.

Definition at line 110 of file ConeSurface.cpp.

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

References bounds, Acts::eBoundLoc0, Acts::eBoundLoc1, ActsTests::PropagationDatasets::phi, physmon_track_finding_ttbar::r, and Acts::Test::transform.

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

Return properly formatted class name for screen output.

Implements Acts::Surface.

Definition at line 148 of file ConeSurface.cpp.

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

Acts::Vector3 Acts::ConeSurface::normal ( const GeometryContext gctx,
const Vector2 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
Vector3 normal vector in global frame

Implements Acts::Surface.

Definition at line 152 of file ConeSurface.cpp.

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

References bounds, Acts::ConeBounds::eAlpha, Acts::eBoundLoc0, Acts::eBoundLoc1, ActsTests::PropagationDatasets::phi, and Acts::Test::transform.

Acts::Vector3 Acts::ConeSurface::normal ( const GeometryContext gctx,
const Vector3 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
Vector3 normal vector in global frame

Reimplemented from Acts::Surface.

Definition at line 164 of file ConeSurface.cpp.

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

References position, and Acts::Test::transform.

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

Assignment operator

Parameters
otheris the source surface for the assignment

Definition at line 76 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 76 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 Vector3 position,
const Vector3 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 132 of file ConeSurface.cpp.

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

References bounds, Acts::ConeBounds::eAlpha, ActsTests::PropagationDatasets::phi, Acts::VectorHelpers::phi(), 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 178 of file ConeSurface.cpp.

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

References bounds, Acts::detail::VerticesHelper::createSegment(), Acts::detail::FacesHelper::cylindricalFaceMesh(), Acts::ConeBounds::eAveragePhi, Acts::ConeBounds::eHalfPhiSector, Acts::ConeBounds::eMaxZ, Acts::ConeBounds::eMinZ, maxZ, Acts::detail::VerticesHelper::phiSegments(), physmon_track_finding_ttbar::r, Acts::s_onSurfaceTolerance, swap(), Acts::Test::transform, ActsExamples::HepMC3Event::vertices(), and physmon_track_finding_ttbar::z.

+ Here is the call graph for this function:

Acts::RotationMatrix3 Acts::ConeSurface::referenceFrame ( const GeometryContext gctx,
const Vector3 position,
const Vector3 direction 
) const
finalvirtual

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

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

<

Reimplemented from Acts::Surface.

Definition at line 89 of file ConeSurface.cpp.

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

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

Definition at line 84 of file ConeSurface.cpp.

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

References Acts::Test::transform.

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

Return the surface type.

Implements Acts::Surface.

Definition at line 72 of file ConeSurface.cpp.

View newest version in sPHENIX GitHub at line 72 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 239 of file ConeSurface.hpp.

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

Referenced by ConeSurface(), and operator=().

friend Acts::ConeSurface::Surface
private

Definition at line 44 of file ConeSurface.hpp.

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


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