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

Abstract Base Class for tracking surfaces. More...

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

+ Inheritance diagram for Acts::Surface:
+ Collaboration diagram for Acts::Surface:

Public Types

enum  SurfaceType {
  Cone = 0, Cylinder = 1, Disc = 2, Perigee = 3,
  Plane = 4, Straw = 5, Curvilinear = 6, Other = 7
}
 

Public Member Functions

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 SurfaceType type () const =0
 Return method for the Surface type to avoid dynamic casts.
 
virtual const Transform3Dtransform (const GeometryContext &gctx) const
 
virtual Vector3D center (const GeometryContext &gctx) const
 
virtual Vector3D normal (const GeometryContext &gctx, const Vector2D &lposition) const =0
 
virtual Vector3D normal (const GeometryContext &gctx, const Vector3D &position) const
 
virtual Vector3D normal (const GeometryContext &gctx) const
 
virtual const SurfaceBoundsbounds () const =0
 
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 Vector3D localToGlobal (const GeometryContext &gctx, const Vector2D &lposition, const Vector3D &momentum) const =0
 
virtual Result< Vector2DglobalToLocal (const GeometryContext &gctx, const Vector3D &position, const Vector3D &momentum) const =0
 
virtual Acts::RotationMatrix3D referenceFrame (const GeometryContext &gctx, const Vector3D &position, const Vector3D &momentum) 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 double pathCorrection (const GeometryContext &gctx, const Vector3D &position, const Vector3D &direction) const =0
 
virtual SurfaceIntersection intersect (const GeometryContext &gctx, const Vector3D &position, const Vector3D &direction, const BoundaryCheck &bcheck) const =0
 
virtual std::ostream & toStream (const GeometryContext &gctx, std::ostream &sl) const
 
virtual std::string name () const =0
 Return properly formatted class name.
 
virtual Polyhedron polyhedronRepresentation (const GeometryContext &gctx, size_t lseg) const =0
 
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
 
virtual
LocalCartesianToBoundLocalMatrix 
localCartesianToBoundLocalDerivative (const GeometryContext &gctx, const Vector3D &position) const =0
 
- 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 Vector3D binningPosition (const GeometryContext &gctx, BinningValue bValue) const =0
 
virtual double binningPositionValue (const GeometryContext &gctx, BinningValue bValue) const
 
void assignGeometryId (const GeometryIdentifier &geometryId)
 

Static Public Member Functions

template<class T , typename... Args>
static std::shared_ptr< TmakeShared (Args &&...args)
 

Protected Member Functions

 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

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
 

Detailed Description

Abstract Base Class for tracking surfaces.

The Surface class builds the core of the Acts Tracking Geometry. All other geometrical objects are either extending the surface or are built from it.

Surfaces are either owned by Detector elements or the Tracking Geometry, in which case they are not copied within the data model objects.

Definition at line 52 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 52 of file Surface.hpp

Member Enumeration Documentation

enum Acts::Surface::SurfaceType

This enumerator simplifies the persistency & calculations, by saving a dynamic_cast, e.g. for persistency

Enumerator:
Cone 
Cylinder 
Disc 
Perigee 
Plane 
Straw 
Curvilinear 
Other 

Definition at line 59 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 59 of file Surface.hpp

Constructor & Destructor Documentation

Acts::Surface::Surface ( const Transform3D transform = Transform3D::Identity())
protected

Constructor with Transform3D as a shared object

Parameters
transformTransform3D positions the surface in 3D global space
Note
also acts as default constructor

Definition at line 17 of file Surface.cpp.

View newest version in sPHENIX GitHub at line 17 of file Surface.cpp

Acts::Surface::Surface ( const Surface other)
protected

Copy constructor

Note
copy construction invalidates the association to detector element and layer
Parameters
otherSource surface for copy.

Definition at line 23 of file Surface.cpp.

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

Acts::Surface::Surface ( const DetectorElementBase detelement)
protected

Constructor fromt DetectorElementBase: Element proxy

Parameters
detelementDetector element which is represented by this surface

Definition at line 20 of file Surface.cpp.

View newest version in sPHENIX GitHub at line 20 of file Surface.cpp

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

Copy constructor with optional shift

Note
copy construction invalidates the association to detector element and layer
Parameters
gctxThe current geometry context object, e.g. alignment
otherSource surface for copy
shiftAdditional transform applied as: shift * transform

Definition at line 29 of file Surface.cpp.

View newest version in sPHENIX GitHub at line 29 of file Surface.cpp

Acts::Surface::~Surface ( )
virtualdefault

Member Function Documentation

Acts::AlignmentToBoundMatrix Acts::Surface::alignmentToBoundDerivative ( const GeometryContext gctx,
const FreeVector derivatives,
const Vector3D position,
const Vector3D direction 
) const

The derivative of bound track parameters w.r.t. alignment parameters of its reference 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
derivativesPath length derivatives of the free, nominal parameters to help evaluate change of free track parameters caused by change of alignment parameters
positionThe position of the paramters in global
directionThe direction of the track
Returns
Derivative of bound track parameters w.r.t. local frame alignment parameters

Definition at line 50 of file Surface.cpp.

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

References Acts::eAlignmentCenter0, Acts::eAlignmentRotation0, Acts::eAlignmentSize, Acts::eBoundLoc0, Acts::eBoundPhi, Acts::eFreePos0, Acts::eFreeSize, Acts::eX, Acts::eY, Acts::eZ, Acts::detail::rotationToLocalAxesDerivative(), and Acts::Test::transform.

+ Here is the call graph for this function:

Acts::AlignmentRowVector Acts::Surface::alignmentToPathDerivative ( const GeometryContext gctx,
const RotationMatrix3D rotToLocalZAxis,
const Vector3D position,
const Vector3D direction 
) const
virtual

Calculate the derivative of path length 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)

Re-implementation is needed for surface whose intersection with track is not its local xy plane, e.g. LineSurface, CylinderSurface and ConeSurface

Parameters
gctxThe current geometry context object, e.g. alignment
rotToLocalZAxisThe derivative of local frame z axis vector w.r.t. its rotation
positionThe position of the paramters in global
directionThe direction of the track
Returns
Derivative of path length w.r.t. the alignment parameters

Reimplemented in Acts::LineSurface.

Definition at line 109 of file Surface.cpp.

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

References Acts::eAlignmentCenter0, Acts::eAlignmentRotation0, and Acts::Test::transform.

void Surface::assignSurfaceMaterial ( std::shared_ptr< const ISurfaceMaterial material)
inline

Assign the surface material description

The material is usually derived in a complicated way and loaded from a framework given source. As various surfaces may share the same source this is provided by a shared pointer

Parameters
materialMaterial description associated to this surface

Definition at line 134 of file Surface.ipp.

View newest version in sPHENIX GitHub at line 134 of file Surface.ipp

Referenced by Acts::addLayerProtoMaterial(), Acts::TrackingVolume::assignBoundaryMaterial(), Acts::Test::BOOST_AUTO_TEST_CASE(), ActsExamples::Generic::LayerBuilderT< detector_element_t >::centralLayers(), ActsExamples::Generic::LayerBuilderT< detector_element_t >::constructEndcapLayers(), ActsExamples::Generic::GenericDetectorElement::GenericDetectorElement(), Acts::TrackingVolume::glueTrackingVolume(), Acts::CylinderVolumeHelper::glueTrackingVolumes(), Acts::Test::CylindricalTrackingGeometry::operator()(), and Acts::TrackingVolume::updateBoundarySurface().

+ Here is the caller graph for this function:

const DetectorElementBase * Surface::associatedDetectorElement ( ) const
inline

Return method for the associated Detector Element

Returns
plain pointer to the DetectorElement, can be nullptr

Definition at line 117 of file Surface.ipp.

View newest version in sPHENIX GitHub at line 117 of file Surface.ipp

Referenced by Acts::Test::BOOST_AUTO_TEST_CASE(), Acts::KalmanFitter< propagator_t, updater_t, smoother_t, outlier_finder_t, calibrator_t >::Actor< source_link_t, parameters_t >::filter(), Acts::CombinatorialKalmanFilter< propagator_t, updater_t, smoother_t, source_link_selector_t, branch_stopper_t, calibrator_t >::Actor< source_link_t, parameters_t >::filter(), Acts::SurfaceSelector::operator()(), and Acts::EventDataView3DTest::testMultiTrajectory().

+ Here is the caller graph for this function:

const Layer * Surface::associatedLayer ( ) const
inline

Return method for the associated Layer in which the surface is embedded

Returns
Layer by plain pointer, can be nullptr

Definition at line 121 of file Surface.ipp.

View newest version in sPHENIX GitHub at line 121 of file Surface.ipp

Referenced by Acts::Test::BOOST_AUTO_TEST_CASE().

+ Here is the caller graph for this function:

void Surface::associateLayer ( const Layer lay)
inline

Set Associated Layer Many surfaces can be associated to a Layer, but it might not be known yet during construction of the layer, this can be set afterwards

Parameters
laythe assignment Layer by reference

Definition at line 139 of file Surface.ipp.

View newest version in sPHENIX GitHub at line 139 of file Surface.ipp

Referenced by Acts::LayerCreator::associateSurfacesToLayer(), Acts::Test::BOOST_AUTO_TEST_CASE(), Acts::PlaneLayer::buildApproachDescriptor(), Acts::CylinderLayer::CylinderLayer(), Acts::Test::CubicTrackingGeometry::operator()(), Acts::PlaneLayer::PlaneLayer(), and Acts::GenericApproachDescriptor::registerLayer().

+ Here is the caller graph for this function:

virtual const SurfaceBounds& Acts::Surface::bounds ( ) const
pure virtual
Vector3D Surface::center ( const GeometryContext gctx) const
inlinevirtual

Return method for the surface center by reference

Note
the center is always recalculated in order to not keep a cache
Parameters
gctxThe current geometry context object, e.g. alignment
Returns
center position by value

Definition at line 9 of file Surface.ipp.

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

References Acts::Test::transform.

Referenced by Acts::Test::BOOST_AUTO_TEST_CASE(), Acts::PlaneLayer::buildApproachDescriptor(), Acts::AtlasStepper< bfield_t >::covarianceTransport(), Acts::ImpactPointEstimator< input_track_t, propagator_t, propagator_options_t >::getDistanceAndMomentum(), normal(), and Acts::EventDataView3DTest::testMultiTrajectory().

+ Here is the caller graph for this function:

BoundRowVector Surface::derivativeFactors ( const GeometryContext gctx,
const Vector3D position,
const Vector3D direction,
const RotationMatrix3D rft,
const BoundToFreeMatrix jacobian 
) const
inlinevirtual

Calculate the form factors for the derivatives the calculation is identical for all surfaces where the reference frame does not depend on the direction

Parameters
gctxThe current geometry context object, e.g. alignment
positionis the position of the paramters in global
directionis the direction of the track
rftis the transposed reference frame (avoids recalculation)
jacobianis the transport jacobian
Returns
a five-dim vector

Reimplemented in Acts::LineSurface.

Definition at line 106 of file Surface.ipp.

View newest version in sPHENIX GitHub at line 106 of file Surface.ipp

References Acts::eBoundSize.

std::shared_ptr< Acts::Surface > Acts::Surface::getSharedPtr ( )

Retrieve a std::shared_ptr for this surface (non-const version)

Note
Will error if this was not created through the makeShared factory since it needs access to the original reference. In C++14 this is undefined behavior (but most likely implemented as a bad_weak_ptr exception), in C++17 it is defined as that exception.
Only call this if you need shared ownership of this object.
Returns
The shared pointer

Definition at line 132 of file Surface.cpp.

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

Referenced by Acts::Test::PropagatorState::Stepper::boundState(), Acts::detail::boundState(), Acts::AtlasStepper< bfield_t >::boundState(), Acts::DD4hepLayerBuilder::createSensitiveSurface(), Acts::KalmanFitter< propagator_t, updater_t, smoother_t, outlier_finder_t, calibrator_t >::Actor< source_link_t, parameters_t >::filter(), ActsExamples::SimSourceLink::operator*(), ActsExamples::CsvPlanarClusterReader::read(), and Acts::EventDataView3DTest::testMultiTrajectory().

+ Here is the caller graph for this function:

std::shared_ptr< const Acts::Surface > Acts::Surface::getSharedPtr ( ) const

Retrieve a std::shared_ptr for this surface (const version)

Note
Will error if this was not created through the makeShared factory since it needs access to the original reference. In C++14 this is undefined behavior, but most likely implemented as a bad_weak_ptr exception, in C++17 it is defined as that exception.
Only call this if you need shared ownership of this object.
Returns
The shared pointer

Definition at line 136 of file Surface.cpp.

View newest version in sPHENIX GitHub at line 136 of file Surface.cpp

virtual Result<Vector2D> Acts::Surface::globalToLocal ( const GeometryContext gctx,
const Vector3D position,
const Vector3D momentum 
) const
pure virtual

Global to local transformation Generalized global to local transformation for the surface types. Since some surface types need the global momentum/direction to resolve sign ambiguity this is also provided

Parameters
gctxThe current geometry context object, e.g. alignment
positionglobal 3D position - considered to be on surface but not inside bounds (check is done)
momentumglobal 3D momentum representation (optionally ignored)
Returns
a Result<Vector2D> which can be !ok() if the operation fails

Implemented in Acts::LineSurface, Acts::CylinderSurface, Acts::ConeSurface, Acts::PlaneSurface, and Acts::SurfaceStub.

Referenced by ActsExamples::HitSmearing::execute(), ActsFatras::detail::Interactor< generator_t, physics_list_t, hit_surface_selector_t >::operator()(), ActsExamples::CsvPlanarClusterReader::read(), and Acts::detail::transformFreeToBoundParameters().

+ Here is the caller graph for this function:

void Surface::initJacobianToGlobal ( const GeometryContext gctx,
BoundToFreeMatrix jacobian,
const Vector3D position,
const Vector3D direction,
const BoundVector pars 
) const
inlinevirtual

Initialize the jacobian from local to global the surface knows best, hence the calculation is done here. The jacobian is assumed to be initialised, so only the relevant entries are filled

Parameters
gctxThe current geometry context object, e.g. alignment
jacobianis the jacobian to be initialized
positionis the global position of the parameters
directionis the direction at of the parameters
parsis the parameter vector

Reimplemented in Acts::LineSurface.

Definition at line 39 of file Surface.ipp.

View newest version in sPHENIX GitHub at line 39 of file Surface.ipp

References Acts::eBoundPhi, Acts::eBoundQOverP, Acts::eBoundTheta, Acts::eBoundTime, x, y, and z.

Referenced by Acts::StraightLineStepper::resetState(), and Acts::EigenStepper< bfield_t, extensionlist_t, auctioneer_t >::resetState().

+ Here is the caller graph for this function:

RotationMatrix3D Surface::initJacobianToLocal ( const GeometryContext gctx,
FreeToBoundMatrix jacobian,
const Vector3D position,
const Vector3D direction 
) const
inlinevirtual

Initialize the jacobian from global to local the surface knows best, hence the calculation is done here. The jacobian is assumed to be initialised, so only the relevant entries are filled

Parameters
jacobianis the jacobian to be initialized
positionis the global position of the parameters
directionis the direction at of the parameters
gctxThe current geometry context object, e.g. alignment
Returns
the transposed reference frame (avoids recalculation)

Definition at line 75 of file Surface.ipp.

View newest version in sPHENIX GitHub at line 75 of file Surface.ipp

References Acts::eBoundPhi, Acts::eBoundQOverP, Acts::eBoundTheta, Acts::eBoundTime, x, y, and z.

bool Surface::insideBounds ( const Vector2D lposition,
const BoundaryCheck bcheck = true 
) const
inlinevirtual

The insideBounds method for local positions

Parameters
lpositionThe local position to check
bcheckBoundaryCheck directive for this onSurface check
Returns
boolean indication if operation was successful

Definition at line 28 of file Surface.ipp.

View newest version in sPHENIX GitHub at line 28 of file Surface.ipp

Referenced by Acts::Test::BOOST_AUTO_TEST_CASE().

+ Here is the caller graph for this function:

virtual SurfaceIntersection Acts::Surface::intersect ( const GeometryContext gctx,
const Vector3D position,
const Vector3D direction,
const BoundaryCheck bcheck 
) const
pure virtual

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
Returns
SurfaceIntersection object (contains intersection & surface)

Implemented in Acts::LineSurface, Acts::CylinderSurface, Acts::PlaneSurface, Acts::ConeSurface, and Acts::SurfaceStub.

Referenced by Acts::Layer::compatibleSurfaces(), Acts::TrackingVolume::compatibleSurfacesFromHierarchy(), Acts::SurfaceReached::operator()(), Acts::Layer::surfaceOnApproach(), Acts::detail::updateSingleSurfaceStatus(), and ActsExamples::RootMaterialTrackWriter::writeT().

+ Here is the caller graph for this function:

bool Acts::Surface::isOnSurface ( const GeometryContext gctx,
const Vector3D position,
const Vector3D momentum,
const BoundaryCheck bcheck = true 
) const

The geometric onSurface method

Geometrical check whether position is on Surface

Parameters
gctxThe current geometry context object, e.g. alignment
positionglobal position to be evaludated
momentumglobal momentum (required for line-type surfaces)
bcheckBoundaryCheck directive for this onSurface check
Returns
boolean indication if operation was successful

Definition at line 38 of file Surface.cpp.

View newest version in sPHENIX GitHub at line 38 of file Surface.cpp

Referenced by Acts::Test::BOOST_AUTO_TEST_CASE(), and Acts::BoundarySurfaceT< volume_t >::onBoundary().

+ Here is the caller graph for this function:

virtual LocalCartesianToBoundLocalMatrix Acts::Surface::localCartesianToBoundLocalDerivative ( const GeometryContext gctx,
const Vector3D position 
) const
pure virtual

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

Implemented in Acts::LineSurface, Acts::CylinderSurface, Acts::ConeSurface, Acts::PlaneSurface, and Acts::SurfaceStub.

virtual Vector3D Acts::Surface::localToGlobal ( const GeometryContext gctx,
const Vector2D lposition,
const Vector3D momentum 
) const
pure virtual

Local to global transformation Generalized local to global transformation for the surface types. Since some surface types need the global momentum/direction to resolve sign ambiguity this is also provided

Parameters
gctxThe current geometry context object, e.g. alignment
lpositionlocal 2D position in specialized surface frame
momentumglobal 3D momentum representation (optionally ignored)
Returns
The global position by value

Implemented in Acts::LineSurface, Acts::CylinderSurface, Acts::ConeSurface, Acts::PlaneSurface, and Acts::SurfaceStub.

Referenced by Acts::SurfaceArrayCreator::makeGlobalVertices(), PHActsSiliconSeeding::makeSpacePoint(), and Acts::detail::transformBoundToFreeParameters().

+ Here is the caller graph for this function:

template<class T , typename... Args>
static std::shared_ptr<T> Acts::Surface::makeShared ( Args &&...  args)
inlinestatic

Factory for producing memory managed instances of Surface. Will forward all parameters and will attempt to find a suitable constructor.

Definition at line 108 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 108 of file Surface.hpp

References charm_jet_coverage::args, and T.

virtual std::string Acts::Surface::name ( ) const
pure virtual

Return properly formatted class name.

Implemented in Acts::LineSurface, Acts::ConeSurface, Acts::CylinderSurface, Acts::PlaneSurface, Acts::SurfaceStub, Acts::StrawSurface, and Acts::PerigeeSurface.

Referenced by eicpy.verify.PythiaHistograms::__init__(), eicpy.verify.DjangohHistograms::__init__(), Acts::JsonGeometryConverter::DefaultBin(), and Acts::detail::printBoundParameters().

+ Here is the caller graph for this function:

virtual Vector3D Acts::Surface::normal ( const GeometryContext gctx,
const Vector2D lposition 
) const
pure virtual

Return method for the normal vector of the surface The normal vector can only be generally defined at a given local position It requires a local position to be given (in general)

Parameters
gctxThe current geometry context object, e.g. alignment
lpositionis the local position where the normal vector is constructed
Returns
normal vector by value

Implemented in Acts::CylinderSurface, Acts::ConeSurface, Acts::PlaneSurface, Acts::LineSurface, and Acts::SurfaceStub.

Referenced by Acts::PlaneLayer::buildApproachDescriptor(), normal(), and ActsFatras::detail::Interactor< generator_t, physics_list_t, hit_surface_selector_t >::operator()().

+ Here is the caller graph for this function:

virtual Vector3D Acts::Surface::normal ( const GeometryContext gctx,
const Vector3D position 
) const
virtual

Return method for the normal vector of the surface The normal vector can only be generally defined at a given local position It requires a local position to be given (in general)

Parameters
positionis the global position where the normal vector is constructed
gctxThe current geometry context object, e.g. alignment
Returns
normal vector by value

Reimplemented in Acts::CylinderSurface, Acts::ConeSurface, and Acts::SurfaceStub.

virtual Vector3D Acts::Surface::normal ( const GeometryContext gctx) const
inlinevirtual

Return method for the normal vector of the surface

It will return a normal vector at the center() position

Parameters
gctxThe current geometry context object, e.g. alignment
Returns
normal vector by value

Reimplemented in Acts::SurfaceStub.

Definition at line 210 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 210 of file Surface.hpp

References center(), and normal().

+ Here is the call graph for this function:

bool Acts::Surface::operator!= ( const Surface sf) const
virtual

Comparison (non-equality) operator

Parameters
sfSource surface for the comparison

Definition at line 207 of file Surface.cpp.

View newest version in sPHENIX GitHub at line 207 of file Surface.cpp

References Acts::operator==().

+ Here is the call graph for this function:

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

Assignment operator

Note
copy construction invalidates the association to detector element and layer
Parameters
otherSource surface for the assignment

Definition at line 140 of file Surface.cpp.

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

References m_associatedDetElement, m_associatedLayer, m_surfaceMaterial, m_transform(), m_transform, and Acts::GeometryObject::operator=().

Referenced by Acts::LineSurface::operator=(), Acts::PlaneSurface::operator=(), Acts::ConeSurface::operator=(), and Acts::CylinderSurface::operator=().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

bool Acts::Surface::operator== ( const Surface other) const
virtual

Comparison (equality) operator The strategy for comparison is (a) first pointer comparison (b) then type comparison (c) then bounds comparison (d) then transform comparison

Parameters
othersource surface for the comparison

Definition at line 152 of file Surface.cpp.

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

References bounds(), Acts::UnitConstants::e, m_associatedDetElement, m_surfaceMaterial, m_transform(), m_transform, and type().

+ Here is the call graph for this function:

virtual double Acts::Surface::pathCorrection ( const GeometryContext gctx,
const Vector3D position,
const Vector3D direction 
) const
pure virtual

Calucation of the path correction for incident

Parameters
gctxThe current geometry context object, e.g. alignment
positionglobal 3D position - considered to be on surface but not inside bounds (check is done)
directionglobal 3D momentum direction
Returns
Path correction with respect to the nominal incident.

Implemented in Acts::LineSurface, Acts::CylinderSurface, Acts::ConeSurface, Acts::PlaneSurface, and Acts::SurfaceStub.

Referenced by Acts::Layer::compatibleSurfaces(), and Acts::detail::PointwiseMaterialInteraction::evaluateMaterialSlab().

+ Here is the caller graph for this function:

virtual Polyhedron Acts::Surface::polyhedronRepresentation ( const GeometryContext gctx,
size_t  lseg 
) const
pure virtual

Return a Polyhedron for this object

Parameters
gctxThe current geometry context object, e.g. alignment
lsegNumber of segments along curved lines, if the lseg is set to one, only the corners and the extrema are given, otherwise it represents the number of segments for a full 2*M_PI circle and is scaled to the relevant sector
Note
An internal surface transform can invalidate the extrema in the transformed space
Returns
A list of vertices and a face/facett description of it

Implemented in Acts::CylinderSurface, Acts::ConeSurface, Acts::PlaneSurface, Acts::SurfaceStub, Acts::StrawSurface, Acts::PerigeeSurface, and Acts::Test::LineSurfaceStub.

Referenced by Acts::GeometryView3D::drawSurface(), and Acts::SurfaceBinningMatcher::operator()().

+ Here is the caller graph for this function:

RotationMatrix3D Surface::referenceFrame ( const GeometryContext gctx,
const Vector3D position,
const Vector3D momentum 
) const
inlinevirtual

Return mehtod for the reference frame This is the frame in which the covariance matrix is defined (specialized by all surfaces)

Parameters
gctxThe current geometry context object, e.g. alignment
positionglobal 3D position - considered to be on surface but not inside bounds (check is done)
momentumglobal 3D momentum representation (optionally ignored)
Returns
RotationMatrix3D which defines the three axes of the measurement frame

Reimplemented in Acts::LineSurface, Acts::ConeSurface, and Acts::CylinderSurface.

Definition at line 33 of file Surface.ipp.

View newest version in sPHENIX GitHub at line 33 of file Surface.ipp

References Acts::Test::transform.

Referenced by Acts::Test::BOOST_AUTO_TEST_CASE(), Acts::AtlasStepper< bfield_t >::covarianceTransport(), and Acts::AtlasStepper< bfield_t >::resetState().

+ Here is the caller graph for this function:

const std::shared_ptr< const ISurfaceMaterial > & Surface::surfaceMaterialSharedPtr ( ) const
inline

Return method for the shared pointer to the associated Material

Returns
SurfaceMaterial as shared_pointer, can be nullptr

Definition at line 130 of file Surface.ipp.

View newest version in sPHENIX GitHub at line 130 of file Surface.ipp

Referenced by ActsExamples::RootMaterialWriter::collectMaterial(), and Acts::TrackingVolume::glueTrackingVolume().

+ Here is the caller graph for this function:

std::ostream & Acts::Surface::toStream ( const GeometryContext gctx,
std::ostream &  sl 
) const
virtual

Output Method for std::ostream, to be overloaded by child classes

Parameters
gctxThe current geometry context object, e.g. alignment
slis the ostream to be dumped into

Reimplemented in Acts::PerigeeSurface.

Definition at line 183 of file Surface.cpp.

View newest version in sPHENIX GitHub at line 183 of file Surface.cpp

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

Referenced by PHActsTrkFitter::getSourceLinks().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

const Transform3D & Surface::transform ( const GeometryContext gctx) const
inlinevirtual

Return method for the surface Transform3D by reference In case a detector element is associated the surface transform is just forwarded to the detector element in order to keep the (mis-)alignment cache cetrally handled

Parameters
gctxThe current geometry context object, e.g. alignment
Returns
the contextual transform

Definition at line 20 of file Surface.ipp.

View newest version in sPHENIX GitHub at line 20 of file Surface.ipp

References m_transform().

Referenced by Acts::JsonGeometryConverter::addSurfaceToJson(), Acts::adjustBinUtility(), Acts::Test::BOOST_AUTO_TEST_CASE(), Acts::PlaneLayer::buildApproachDescriptor(), Acts::LayerArrayCreator::createNavigationSurface(), and Acts::ImpactPointEstimator< input_track_t, propagator_t, propagator_options_t >::get3dVertexCompatibility().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

virtual SurfaceType Acts::Surface::type ( ) const
pure virtual

Member Data Documentation

const DetectorElementBase* Acts::Surface::m_associatedDetElement {nullptr}
protected

Pointer to the a DetectorElementBase.

Definition at line 481 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 481 of file Surface.hpp

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

const Layer* Acts::Surface::m_associatedLayer {nullptr}
protected

The associated layer Layer - layer in which the Surface is be embedded, nullptr if not associated

Definition at line 485 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 485 of file Surface.hpp

Referenced by operator=().

const TrackingVolume* Acts::Surface::m_associatedTrackingVolume {nullptr}
protected

The assoicated TrackingVolume - tracking volume in case the surface is a boundary surface, nullptr if not associated

Definition at line 489 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 489 of file Surface.hpp

std::shared_ptr<const ISurfaceMaterial> Acts::Surface::m_surfaceMaterial
protected

Possibility to attach a material descrption.

Definition at line 492 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 492 of file Surface.hpp

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

Transform3D Acts::Surface::m_transform = Transform3D::Identity()
protected

Transform3D definition that positions (translation, rotation) the surface in global space

Definition at line 478 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 478 of file Surface.hpp

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


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