Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Acts::CylinderSurface Class Reference

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

+ Inheritance diagram for Acts::CylinderSurface:
+ Collaboration diagram for Acts::CylinderSurface:

Public Member Functions

 ~CylinderSurface () override=default
 
 CylinderSurface ()=delete
 
CylinderSurfaceoperator= (const CylinderSurface &other)
 
Vector3 binningPosition (const GeometryContext &gctx, BinningValue bValue) const final
 
RotationMatrix3 referenceFrame (const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction) const final
 
SurfaceType type () const override
 Return the surface type.
 
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 CylinderBoundsbounds () const final
 This method returns the CylinderBounds 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, ActsScalar tolerance=s_onSurfaceTolerance) const final
 
double pathCorrection (const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction) const final
 
std::string name () const override
 Return method for properly formatted output string.
 
Polyhedron polyhedronRepresentation (const GeometryContext &gctx, size_t lseg) const override
 
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

 CylinderSurface (const Transform3 &transform, double radius, double halfz, double halfphi=M_PI, double avphi=0., double bevelMinZ=0., double bevelMaxZ=0.)
 
 CylinderSurface (const Transform3 &transform, std::shared_ptr< const CylinderBounds > cbounds)
 
 CylinderSurface (std::shared_ptr< const CylinderBounds > cbounds, const DetectorElementBase &detelement)
 
 CylinderSurface (const CylinderSurface &other)
 
 CylinderSurface (const GeometryContext &gctx, const CylinderSurface &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
CylinderBounds
m_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 Transform3 &transform, 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 CylinderSurface in the TrackingGeometry. It inherits from Surface.

The cylinder surface has a special role in the TrackingGeometry, since it builds the surfaces of all TrackingVolumes at container level for a cylindrical tracking geometry.

CylinderSurface.png

Definition at line 43 of file CylinderSurface.hpp.

View newest version in sPHENIX GitHub at line 43 of file CylinderSurface.hpp

Constructor & Destructor Documentation

Acts::CylinderSurface::CylinderSurface ( const Transform3 transform,
double  radius,
double  halfz,
double  halfphi = M_PI,
double  avphi = 0.,
double  bevelMinZ = 0.,
double  bevelMaxZ = 0. 
)
protected

Constructor from Transform3 and CylinderBounds

Parameters
transformThe transform to position the surface
radiusThe radius of the cylinder
halfzThe half length in z
halfphiThe half opening angle
avphiThe phi value from which the opening angle spans (both sides)
bevelMinZ(optional) The bevel on the negative z side
bevelMaxZ(optional) The bevel on the positive z sid The bevel on the positive z side

Definition at line 40 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 40 of file CylinderSurface.cpp

Acts::CylinderSurface::CylinderSurface ( const Transform3 transform,
std::shared_ptr< const CylinderBounds cbounds 
)
protected

Constructor from Transform3 and CylinderBounds arguments

Parameters
transformThe transform to position the surface
cboundsis a shared pointer to a cylindeer bounds object, it must exist (assert test)

Definition at line 56 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 56 of file CylinderSurface.cpp

References m_bounds, and throw_assert.

Acts::CylinderSurface::CylinderSurface ( std::shared_ptr< const CylinderBounds cbounds,
const DetectorElementBase detelement 
)
protected

Constructor from DetectorElementBase: Element proxy

Parameters
cboundsare the provided cylinder bounds (shared)
detelementis the linked detector element to this surface

surfaces representing a detector element must have bounds

Definition at line 48 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 48 of file CylinderSurface.cpp

References m_bounds, and throw_assert.

Acts::CylinderSurface::CylinderSurface ( const CylinderSurface other)
protected

Copy constructor

Parameters
otheris the source cylinder for the copy

Definition at line 32 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 32 of file CylinderSurface.cpp

Acts::CylinderSurface::CylinderSurface ( const GeometryContext gctx,
const CylinderSurface 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 35 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 35 of file CylinderSurface.cpp

Acts::CylinderSurface::~CylinderSurface ( )
overridedefault
Acts::CylinderSurface::CylinderSurface ( )
delete

Member Function Documentation

Acts::AlignmentToPathMatrix Acts::CylinderSurface::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 287 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 287 of file CylinderSurface.cpp

References 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::CylinderSurface::binningPosition ( const GeometryContext gctx,
BinningValue  bValue 
) const
finalvirtual

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

Parameters
gctxThe current geometry context object, e.g. alignment
bValueis the type of global binning to be done
Returns
is the global position to be used for binning

Implements Acts::GeometryObject.

Definition at line 72 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 72 of file CylinderSurface.cpp

References Acts::binR, Acts::binRPhi, bounds, Acts::CylinderBounds::eAveragePhi, Acts::CylinderBounds::eR, ActsTests::PropagationDatasets::phi, and Acts::IntegrationTest::R.

const Acts::CylinderBounds & Acts::CylinderSurface::bounds ( ) const
finalvirtual

This method returns the CylinderBounds by reference.

Implements Acts::Surface.

Definition at line 171 of file CylinderSurface.cpp.

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

Referenced by Acts::CylinderVolumeBuilder::analyzeContent().

+ Here is the caller graph for this function:

Acts::Result< Acts::Vector2 > Acts::CylinderSurface::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 119 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 119 of file CylinderSurface.cpp

References bounds, Acts::CylinderBounds::eR, Acts::VectorHelpers::perp(), ActsTests::PropagationDatasets::phi, Acts::s_onSurfaceTolerance, Acts::Test::tolerance, and Acts::Test::transform.

+ Here is the call graph for this function:

Acts::SurfaceMultiIntersection Acts::CylinderSurface::intersect ( const GeometryContext gctx,
const Vector3 position,
const Vector3 direction,
const BoundaryCheck bcheck = false,
ActsScalar  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
SurfaceIntersection object (contains intersection & surface)

Implements Acts::Surface.

Definition at line 219 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 219 of file CylinderSurface.cpp

References Acts::Test::boundaryCheck, bounds, Acts::eBoundLoc1, Acts::CylinderBounds::eHalfLengthZ, Acts::Intersection< 3 >::invalid(), Acts::Intersection< DIM >::pathLength(), utils::status, Acts::BoundaryCheck::tolerance(), Acts::Test::transform, and Acts::BoundaryCheck::type().

+ Here is the call graph for this function:

Acts::detail::RealQuadraticEquation Acts::CylinderSurface::intersectionSolver ( const Transform3 transform,
const Vector3 position,
const Vector3 direction 
) const
private

Implementation of the intersection solver

mathematical motivation:

The cylinder is given by :

  • cylinder center: ccenter (C)
  • the direction of the cylinder axis: cdirection (DZ)
  • the radius r The line is given by :
  • a reference position : lposition (L0)
  • the line direction: ldirection (DL) the parametric form for the line is then : L(t) = L0 + t * DL

Any point P on infinite cylinder if : ((P - C) x DZ)^2 = r^2 * DZ^2 We know that DZ is a unit vector: DZ^2 == 1 When expanded with the line equation, this is : ((L0 - C) x DZ + t * (DL x DZ))^2 = r^2 * DZ^2 which is a quadratic equation in the form (X + t * Y)^2 = d, where : X = (L0 - C) x DZ Y = DL x DZ d = r^2 * (DZ)^2 Now, expand the equation : t^2 * (Y . Y) + t * (2 * (X . Y)) + (X . X) - d = 0 => second order equation in the form : a*t^2 + b*t + c = 0 where a = (Y . Y) b = 2 * (X . Y) c = (X . X) - d finally solve the second order equation : a*t^2 + b*t + c = 0 reinsertion into the line equation.

Returns
the quadratic equation

Definition at line 197 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 197 of file CylinderSurface.cpp

References KFPMath::a, KFPMath::b, bounds, Acts::PhysicalConstants::c, Acts::CylinderBounds::eR, and Acts::IntegrationTest::R.

Acts::ActsMatrix< 2, 3 > Acts::CylinderSurface::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 337 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 337 of file CylinderSurface.cpp

References bounds, Acts::CylinderBounds::eR, 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::CylinderSurface::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 109 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 109 of file CylinderSurface.cpp

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

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

Return method for properly formatted output string.

Implements Acts::Surface.

Definition at line 141 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 141 of file CylinderSurface.cpp

Acts::Vector3 Acts::CylinderSurface::normal ( const GeometryContext gctx,
const Vector2 lposition 
) const
finalvirtual

Return method for surface normal information

Note
for a Cylinder a local position is always required for the normal vector
Parameters
gctxThe current geometry context object, e.g. alignment
lpositionis the local position for which the normal vector is requested
Returns
normal vector at the local position by value

Implements Acts::Surface.

Definition at line 145 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 145 of file CylinderSurface.cpp

References Acts::eBoundLoc0, Acts::CylinderBounds::eR, ActsTests::PropagationDatasets::phi, and Acts::Test::transform.

Acts::Vector3 Acts::CylinderSurface::normal ( const GeometryContext gctx,
const Vector3 position 
) const
finalvirtual

Return method for surface normal information

Note
for a Cylinder a local position is always required for the normal vector
Parameters
gctxThe current geometry context object, e.g. alignment
positionis the global position for which the normal vector is requested
Returns
normal vector at the global position by value

Reimplemented from Acts::Surface.

Definition at line 152 of file CylinderSurface.cpp.

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

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

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

Assignment operator

Parameters
otheris the source cylinder for the copy

Definition at line 62 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 62 of file CylinderSurface.cpp

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

+ Here is the call graph for this function:

double Acts::CylinderSurface::pathCorrection ( const GeometryContext gctx,
const Vector3 position,
const Vector3 direction 
) const
finalvirtual

Path correction due to incident of the track

Parameters
gctxThe current geometry context object, e.g. alignment
positionis the global position as a starting point
directionis the global momentum direction at the starting point
Returns
is the correction factor due to incident

Implements Acts::Surface.

Definition at line 163 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 163 of file CylinderSurface.cpp

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

Return a Polyhedron for a cylinder

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
Returns
A list of vertices and a face/facett description of it

Implements Acts::Surface.

Definition at line 175 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 175 of file CylinderSurface.cpp

References bounds, Acts::detail::FacesHelper::cylindricalFaceMesh(), Acts::Test::transform, and ActsExamples::HepMC3Event::vertices().

+ Here is the call graph for this function:

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

Return the measurement frame - this is needed for alignment, in particular The measurement frame of a cylinder is the tangential plane at a given position

Parameters
gctxThe current geometry context object, e.g. alignment
positionis the position where the measurement frame is defined
directionis the momentum direction vector (ignored)
Returns
rotation matrix that defines the measurement frame

Reimplemented from Acts::Surface.

Definition at line 86 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 86 of file CylinderSurface.cpp

Acts::Vector3 Acts::CylinderSurface::rotSymmetryAxis ( const GeometryContext gctx) const
virtual

Return method for the rotational symmetry axis

Parameters
gctxThe current geometry context object, e.g. alignment
Returns
the z-Axis of transform

Definition at line 191 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 191 of file CylinderSurface.cpp

References Acts::Test::transform.

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

Return the surface type.

Implements Acts::Surface.

Definition at line 105 of file CylinderSurface.cpp.

View newest version in sPHENIX GitHub at line 105 of file CylinderSurface.cpp

References Acts::Surface::Cylinder.

Member Data Documentation

std::shared_ptr<const CylinderBounds> Acts::CylinderSurface::m_bounds
protected

bounds (shared)

Definition at line 248 of file CylinderSurface.hpp.

View newest version in sPHENIX GitHub at line 248 of file CylinderSurface.hpp

Referenced by Acts::CylinderLayer::CylinderLayer(), CylinderSurface(), and operator=().

friend Acts::CylinderSurface::Surface
private

Definition at line 45 of file CylinderSurface.hpp.

View newest version in sPHENIX GitHub at line 45 of file CylinderSurface.hpp


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