Analysis Software
Documentation for sPHENIX simulation 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 Transform3transform (const GeometryContext &gctx) const
 
virtual Vector3 center (const GeometryContext &gctx) const
 
virtual Vector3 normal (const GeometryContext &gctx, const Vector2 &lposition) const =0
 
virtual Vector3 normal (const GeometryContext &gctx, const Vector3 &position) const
 
virtual Vector3 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 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 Vector3 localToGlobal (const GeometryContext &gctx, const Vector2 &lposition, const Vector3 &direction) const =0
 
virtual Result< Vector2globalToLocal (const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction, double tolerance=s_onSurfaceTolerance) const =0
 
virtual Acts::RotationMatrix3 referenceFrame (const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction) 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 double pathCorrection (const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction) const =0
 
virtual SurfaceMultiIntersection intersect (const GeometryContext &gctx, const Vector3 &position, const Vector3 &direction, const BoundaryCheck &bcheck=false, ActsScalar tolerance=s_onSurfaceTolerance) const =0
 
virtual std::ostream & toStream (const GeometryContext &gctx, std::ostream &sl) const
 
std::string toString (const GeometryContext &gctx) 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 &parameters, const FreeVector &pathDerivative) const
 
virtual AlignmentToPathMatrix alignmentToPathDerivative (const GeometryContext &gctx, const FreeVector &parameters) const
 
virtual ActsMatrix< 2, 3 > localCartesianToBoundLocalDerivative (const GeometryContext &gctx, const Vector3 &position) const =0
 
- 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 Vector3 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)
 

Static Public Attributes

static std::array< std::string,
SurfaceType::Other > 
s_surfaceTypeNames
 Helper strings for screen output.
 

Protected Member Functions

 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

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

AlignmentToBoundMatrix alignmentToBoundDerivativeWithoutCorrection (const GeometryContext &gctx, const FreeVector &parameters) const
 

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 63 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 63 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 70 of file Surface.hpp.

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

Constructor & Destructor Documentation

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

Constructor with Transform3 as a shared object

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

Definition at line 26 of file Surface.cpp.

View newest version in sPHENIX GitHub at line 26 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 32 of file Surface.cpp.

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

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

Constructor from DetectorElementBase: Element proxy

Parameters
detelementDetector element which is represented by this surface

Definition at line 29 of file Surface.cpp.

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

Acts::Surface::Surface ( const GeometryContext gctx,
const Surface other,
const Transform3 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 38 of file Surface.cpp.

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

Acts::Surface::~Surface ( )
virtualdefault

Member Function Documentation

Acts::AlignmentToBoundMatrix Acts::Surface::alignmentToBoundDerivative ( const GeometryContext gctx,
const FreeVector parameters,
const FreeVector pathDerivative 
) 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 change of alignment parameters
parametersis the free parameters
pathDerivativeis the derivative of free parameters w.r.t. path length
Returns
Derivative of bound track parameters w.r.t. local frame alignment parameters

Definition at line 58 of file Surface.cpp.

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

Acts::AlignmentToBoundMatrix Acts::Surface::alignmentToBoundDerivativeWithoutCorrection ( const GeometryContext gctx,
const FreeVector parameters 
) const
private

Calculate the derivative of bound track parameters w.r.t. alignment parameters of its reference surface (i.e. origin in global 3D Cartesian coordinates and its rotation represented with extrinsic Euler angles) without any path correction

Note
This function should be used together with alignment to path derivative to get the full alignment to bound derivatives
Parameters
gctxThe current geometry context object, e.g. alignment
parametersis the free parameters
Returns
Derivative of bound track parameters w.r.t. local frame alignment parameters without path correction

Definition at line 79 of file Surface.cpp.

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

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

+ Here is the call graph for this function:

Acts::AlignmentToPathMatrix Acts::Surface::alignmentToPathDerivative ( const GeometryContext gctx,
const FreeVector parameters 
) const
virtual

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)

Note
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
parametersis the free parameters
Returns
Derivative of path length w.r.t. the alignment parameters

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

Definition at line 118 of file Surface.cpp.

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

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

+ Here is the call graph for this function:

void Acts::Surface::assignDetectorElement ( const DetectorElementBase detelement)

Assign a detector element

Parameters
detelementDetector element which is represented by this surface

Definition at line 354 of file Surface.cpp.

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

References m_transform().

+ Here is the call graph for this function:

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

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 362 of file Surface.cpp.

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

References testing::internal::move().

Referenced by Acts::addLayerProtoMaterial(), Acts::TrackingVolume::assignBoundaryMaterial(), Acts::DD4hepDetectorSurfaceFactory::attachSurfaceMaterial(), Acts::Test::BOOST_AUTO_TEST_CASE(), ActsExamples::Generic::LayerBuilderT< detector_element_t >::centralLayers(), ActsExamples::Generic::LayerBuilderT< detector_element_t >::constructEndcapLayers(), Acts::JsonMaterialDecorator::decorate(), Acts::TrackingVolume::glueTrackingVolume(), Acts::CylinderVolumeHelper::glueTrackingVolumes(), Acts::Test::CylindricalTrackingGeometry::operator()(), Acts::Geant4PhysicalVolumeConverter::surface(), ActsExamples::Telescope::TelescopeDetectorElement::TelescopeDetectorElement(), and Acts::TrackingVolume::updateBoundarySurface().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

const Acts::Layer * Acts::Surface::associatedLayer ( ) const

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

Returns
Layer by plain pointer, can be nullptr

Definition at line 341 of file Surface.cpp.

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

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

+ Here is the caller graph for this function:

void Acts::Surface::associateLayer ( const Layer lay)

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 367 of file Surface.cpp.

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

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:

Acts::BoundToFreeMatrix Acts::Surface::boundToFreeJacobian ( const GeometryContext gctx,
const BoundVector boundParams 
) const
virtual

Calculate the jacobian from local to global which the surface knows best, hence the calculation is done here.

Note
In principle, the input could also be a free parameters vector as it could be transformed to a bound parameters. But the transform might fail in case the parameters is not on surface. To avoid the check inside this function, it takes directly the bound parameters as input (then the check might be done where this function is called).
Parameters
gctxThe current geometry context object, e.g. alignment
boundParamsis the bound parameters vector
Returns
Jacobian from local to global

Reimplemented in Acts::LineSurface.

Definition at line 260 of file Surface.cpp.

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

References Acts::eBoundPhi, Acts::eBoundQOverP, Acts::eBoundTheta, Acts::eBoundTime, Acts::eFreeDir0, Acts::eFreeDir1, Acts::eFreeDir2, Acts::eFreePos0, Acts::eFreeQOverP, Acts::eFreeTime, Acts::VectorHelpers::evaluateTrigonomics(), position, and Acts::detail::transformBoundToFreeParameters().

Referenced by Acts::CovarianceCache::CovarianceCache(), Acts::StraightLineStepper::resetState(), Acts::EigenStepper< extensionlist_t, auctioneer_t >::resetState(), Acts::StraightLineStepper::update(), Acts::EigenStepper< extensionlist_t, auctioneer_t >::update(), and Acts::MultiEigenStepperLoop< extensionlist_t, component_reducer_t, auctioneer_t >::ComponentProxy::update().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

Acts::Vector3 Acts::Surface::center ( const GeometryContext gctx) const
virtual

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 230 of file Surface.cpp.

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

References Acts::Test::transform.

Referenced by Acts::Test::BOOST_AUTO_TEST_CASE(), BOOST_AUTO_TEST_CASE(), Acts::PlaneLayer::buildApproachDescriptor(), Acts::Svg::SurfaceConverter::convert(), Acts::EventDataView3DTest::createDetector(), ActsAlignmentStates::get_projectionXY(), Acts::ImpactPointEstimator< input_track_t, propagator_t, propagator_options_t >::getDistanceAndMomentum(), ActsExamples::SeedingFTFAlgorithm::LayerNumbering(), normal(), Acts::Experimental::detail::CenterReferenceGenerator::references(), Acts::TGeoCylinderDiscSplitter::split(), Acts::PortalJsonConverter::toJsonDetray(), Acts::AtlasStepper::transportCovarianceToBound(), and Acts::MultiEigenStepperLoop< extensionlist_t, component_reducer_t, auctioneer_t >::updateSurfaceStatus().

+ Here is the caller graph for this function:

Acts::FreeToBoundMatrix Acts::Surface::freeToBoundJacobian ( const GeometryContext gctx,
const FreeVector parameters 
) const
virtual

Calculate the jacobian from global to local which the surface knows best, hence the calculation is done here.

Note
It assumes the input free parameters is on surface, hence no onSurface check is done inside this function.
Parameters
gctxThe current geometry context object, e.g. alignment
parametersis the free parameters
Returns
Jacobian from global to local

Definition at line 290 of file Surface.cpp.

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

References Acts::eBoundLoc0, Acts::eBoundPhi, Acts::eBoundQOverP, Acts::eBoundTheta, Acts::eBoundTime, Acts::eFreeDir0, Acts::eFreeDir1, Acts::eFreeDir2, Acts::eFreePos0, Acts::eFreeQOverP, Acts::eFreeTime, Acts::VectorHelpers::evaluateTrigonomics(), and position.

Referenced by Acts::detail::boundToBoundTransportJacobian(), and Acts::detail::freeToBoundTransportJacobian().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

Acts::FreeToPathMatrix Acts::Surface::freeToPathDerivative ( const GeometryContext gctx,
const FreeVector parameters 
) const
virtual

Calculate the derivative of path length at the geometry constraint or point-of-closest-approach w.r.t. free parameters. 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
parametersis the free parameters
Returns
Derivative of path length w.r.t. free parameters

Reimplemented in Acts::LineSurface.

Definition at line 318 of file Surface.cpp.

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

References dz, Acts::eFreeDir0, Acts::eFreePos0, and position.

Referenced by Acts::detail::boundToBoundTransportJacobian(), and Acts::detail::freeToBoundTransportJacobian().

+ Here is the caller graph for this function:

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 145 of file Surface.cpp.

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

Referenced by Acts::detail::GsfActor< bethe_heitler_approx_t, traj_t >::addCombinedState(), Acts::Test::BOOST_AUTO_TEST_CASE(), BOOST_AUTO_TEST_CASE(), Acts::detail::boundState(), Acts::Test::PropagatorState::Stepper::boundState(), Acts::AtlasStepper::boundState(), Acts::EDM4hepUtil::detail::convertTrackParametersToEdm4hep(), ActsExamples::TracksToTrajectories::execute(), ActsExamples::TrackParamsEstimationAlgorithm::execute(), Acts::Test::fillTrackState(), Acts::CombinatorialKalmanFilter< propagator_t, traj_t >::findTracks(), PHCosmicsTrkFitter::getTrackFitResult(), PHActsTrkFitter::getTrackFitResult(), Acts::detail::kalmanHandleMeasurement(), Acts::detail::kalmanHandleNoMeasurement(), ActsExamples::CsvPlanarClusterReader::read(), Acts::KalmanFitter< propagator_t, traj_t >::Actor< parameters_t >::reversedFilter(), Acts::MultiEigenStepperLoop< extensionlist_t, component_reducer_t, auctioneer_t >::State::State(), and PHActsGSF::updateTrack().

+ 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 149 of file Surface.cpp.

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

virtual Result<Vector2> Acts::Surface::globalToLocal ( const GeometryContext gctx,
const Vector3 position,
const Vector3 direction,
double  tolerance = s_onSurfaceTolerance 
) 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)
directionglobal 3D momentum direction
toleranceoptional tolerance within which a point is considered valid on surface
Returns
a Result<Vector2> which can be !ok() if the operation fails

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

Referenced by ActsExamples::averageSimHits(), ActsFatras::Particle::boundParameters(), Acts::estimateTrackParamsFromSeed(), ActsFatras::detail::SimulationActor< generator_t, decay_t, interactions_t, hit_surface_selector_t >::operator()(), ActsExamples::CsvPlanarClusterReader::read(), Acts::detail::transformFreeToBoundParameters(), ActsExamples::RootPlanarClusterWriter::writeT(), and ActsExamples::RootTrajectorySummaryWriter::writeT().

+ Here is the caller graph for this function:

bool Acts::Surface::insideBounds ( const Vector2 lposition,
const BoundaryCheck bcheck = true 
) const
virtual

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 249 of file Surface.cpp.

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

References bounds.

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

+ Here is the caller graph for this function:

virtual SurfaceMultiIntersection Acts::Surface::intersect ( const GeometryContext gctx,
const Vector3 position,
const Vector3 direction,
const BoundaryCheck bcheck = false,
ActsScalar  tolerance = s_onSurfaceTolerance 
) 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
tolerancethe tolerance used for the intersection
Returns
SurfaceMultiIntersection 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::CombinatorialKalmanFilter< propagator_t, traj_t >::Actor< source_link_accessor_t, parameters_t >::finalize(), meanFromFree(), Acts::MultiStepperSurfaceReached::operator()(), Acts::detail::CorrectedFreeToBoundTransformer::operator()(), Acts::SurfaceReached::operator()(), Acts::Layer::surfaceOnApproach(), Acts::Experimental::updateCandidates(), Acts::detail::updateSingleSurfaceStatus(), Acts::MultiEigenStepperLoop< extensionlist_t, component_reducer_t, auctioneer_t >::updateStepSize(), ActsExamples::RootTrajectorySummaryWriter::writeT(), and ActsExamples::RootMaterialTrackWriter::writeT().

+ Here is the caller graph for this function:

bool Acts::Surface::isOnSurface ( const GeometryContext gctx,
const Vector3 position,
const Vector3 direction,
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
directionglobal momentum direction (required for line-type surfaces)
bcheckBoundaryCheck directive for this onSurface check
Returns
boolean indication if operation was successful

Definition at line 46 of file Surface.cpp.

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

References bounds.

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

+ Here is the caller graph for this function:

virtual ActsMatrix<2, 3> Acts::Surface::localCartesianToBoundLocalDerivative ( const GeometryContext gctx,
const Vector3 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 parameters 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 Vector3 Acts::Surface::localToGlobal ( const GeometryContext gctx,
const Vector2 lposition,
const Vector3 direction 
) 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
directionglobal 3D momentum direction
Returns
The global position by value

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

Referenced by ActsExamples::HoughTransformSeeder::addMeasurements(), Acts::EDM4hepUtil::detail::convertTrackParametersToEdm4hep(), Acts::SpacePointUtility::globalCoords(), ActsExamples::SeedingFTFAlgorithm::LayerNumbering(), Acts::SurfaceArrayCreator::makeGlobalVertices(), ActsExamples::TGeoITkModuleSplitter::splitBarrelModule(), Acts::Test::stripEnds(), Acts::detail::transformBoundToFreeParameters(), ActsExamples::CsvPlanarClusterWriter::writeT(), and ActsExamples::RootPlanarClusterWriter::writeT().

+ 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 122 of file Surface.hpp.

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

References check_smearing_config::args, and Acts::UnitConstants::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 Acts::detail::printBoundParameters().

+ Here is the caller graph for this function:

virtual Vector3 Acts::Surface::normal ( const GeometryContext gctx,
const Vector2 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::LineSurface, Acts::PlaneSurface, and Acts::SurfaceStub.

Referenced by Acts::PlaneLayer::buildApproachDescriptor(), ActsAlignmentStates::get_projectionXY(), normal(), ActsFatras::detail::SimulationActor< generator_t, decay_t, interactions_t, hit_surface_selector_t >::operator()(), Acts::detail::CorrectedFreeToBoundTransformer::operator()(), Acts::PlaneSurface::pathCorrection(), and Acts::PortalJsonConverter::toJsonDetray().

+ Here is the caller graph for this function:

Acts::Vector3 Acts::Surface::normal ( const GeometryContext gctx,
const Vector3 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.

Definition at line 236 of file Surface.cpp.

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

virtual Vector3 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 224 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 224 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 226 of file Surface.cpp.

View newest version in sPHENIX GitHub at line 226 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 153 of file Surface.cpp.

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

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

Referenced by Acts::PlaneSurface::operator=(), Acts::LineSurface::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 165 of file Surface.cpp.

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

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

+ Here is the call graph for this function:

virtual double Acts::Surface::pathCorrection ( const GeometryContext gctx,
const Vector3 position,
const Vector3 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::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::StrawSurface, Acts::SurfaceStub, Acts::PerigeeSurface, and Acts::Test::LineSurfaceStub.

Referenced by Acts::Svg::SurfaceConverter::convert(), Acts::GeometryView3D::drawSurface(), Acts::SurfaceBinningMatcher::operator()(), and Acts::Experimental::detail::PolyhedronReferenceGenerator< nSEGS, aBARY >::references().

+ Here is the caller graph for this function:

Acts::RotationMatrix3 Acts::Surface::referenceFrame ( const GeometryContext gctx,
const Vector3 position,
const Vector3 direction 
) const
virtual

Return method 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)
directionglobal 3D momentum direction (optionally ignored)
Returns
RotationMatrix3 which defines the three axes of the measurement frame

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

Definition at line 254 of file Surface.cpp.

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

References Acts::Test::transform.

Referenced by Acts::Test::BOOST_AUTO_TEST_CASE(), Acts::SpacePointUtility::globalCoords(), Acts::SpacePointUtility::rhoZCovariance(), Acts::AtlasStepper::transportCovarianceToBound(), Acts::AtlasStepper::update(), and ActsExamples::RootMeasurementWriter::writeT().

+ Here is the caller graph for this function:

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

Return method for the shared pointer to the associated Material

Returns
SurfaceMaterial as shared_pointer, can be nullptr

Definition at line 350 of file Surface.cpp.

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

Referenced by Acts::SurfaceMaterialMapper::checkAndInsert(), 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 196 of file Surface.cpp.

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

References bounds, perf_headwind::name, and Acts::Test::transform.

std::string Acts::Surface::toString ( const GeometryContext gctx) const

Output into a std::string

Parameters
gctxThe current geometry context object, e.g. alignment

Definition at line 220 of file Surface.cpp.

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

const Acts::Transform3 & Acts::Surface::transform ( const GeometryContext gctx) const
virtual

Return method for the surface Transform3 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 241 of file Surface.cpp.

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

References m_transform().

Referenced by Acts::adjustBinUtility(), Acts::Experimental::detail::PortalHelper::attachDetectorVolumesUpdator(), Acts::Test::BOOST_AUTO_TEST_CASE(), Acts::PlaneLayer::buildApproachDescriptor(), Acts::Svg::PortalConverter::convert(), Acts::PodioUtil::convertSurfaceToPodio(), Acts::LayerArrayCreator::createNavigationSurface(), ActsAlignmentStates::get_projectionXY(), Acts::ImpactPointEstimator< input_track_t, propagator_t, propagator_options_t >::getVertexCompatibility(), ActsExamples::TGeoITkModuleSplitter::splitBarrelModule(), ActsExamples::TGeoITkModuleSplitter::splitDiscModule(), DD4hepTestsHelper::surfaceToXML(), Acts::SurfaceJsonConverter::toJson(), Acts::PortalJsonConverter::toJsonDetray(), Acts::SurfaceJsonConverter::toJsonDetray(), and ActsFatras::PlanarSurfaceDrift::toReadout().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

Member Data Documentation

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

Pointer to the a DetectorElementBase.

Definition at line 498 of file Surface.hpp.

View newest version in sPHENIX GitHub at line 498 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 502 of file Surface.hpp.

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

Referenced by operator=().

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

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

Definition at line 506 of file Surface.hpp.

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

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

Possibility to attach a material description.

Definition at line 509 of file Surface.hpp.

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

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

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

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

Definition at line 495 of file Surface.hpp.

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

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

std::array< std::string, Acts::Surface::SurfaceType::Other > Acts::Surface::s_surfaceTypeNames
static
Initial value:
= {
"Cone", "Cylinder", "Disc", "Perigee", "Plane", "Straw", "Curvilinear"}

Helper strings for screen output.

Definition at line 82 of file Surface.hpp.

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

Referenced by Acts::Experimental::LayerStructureBuilder::construct().


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