GeneralBrokenLines V03-00-00
gblpy3
Public Member Functions | Private Attributes | List of all members
gblpy3.gblfit.GblPoint Class Reference

User supplied point on (initial) trajectory. More...

Inheritance diagram for gblpy3.gblfit.GblPoint:

Public Member Functions

def __init__ (self, aJacobian)
 Create new point. More...
 
def addMeasurement (self, aMeasurement, minPrecision=0.)
 Add a mesurement to a point. More...
 
def hasMeasurement (self)
 Check point for a measurement. More...
 
def getMeasurements (self)
 Retrieve measurements of a point. More...
 
def addScatterer (self, aScatterer)
 Add a (thin or thick) scatterer to a point. More...
 
def getScatDim (self)
 Get scatterer dimension. More...
 
def getScatterer (self)
 Retrieve scatterer of a point. More...
 
def getReducedScatterer (self)
 Retrieve (reduced) scatterer of a point. More...
 
def addLocals (self, derivatives)
 Add local derivatives (to last measurement). More...
 
def addGlobals (self, labels, derivatives)
 Add global derivatives (to last measurement). More...
 
def setLabel (self, aLabel)
 Define label of a point. More...
 
def getLabel (self)
 Retrieve label of a point. More...
 
def setOffset (self, anOffset)
 Define offset of a point and references to previous and next point with offset. More...
 
def getOffset (self)
 Get offset of a point. More...
 
def setType (self, aType)
 Define type of a point. More...
 
def isFirst (self)
 Is first point? More...
 
def isLast (self)
 Is last point? More...
 
def getP2pJacobian (self)
 Retrieve jacobian of a point. More...
 
def addPrevJacobian (self, aJacobian)
 Add jacobian to previous offset. More...
 
def addNextJacobian (self, aJacobian)
 Add jacobian to next offset. More...
 
def getDerivatives (self, index)
 Get derivatives for locally linearized track model (backward or forward propagation). More...
 
def printPoint (self)
 Print point. More...
 

Private Attributes

 __label
 label for referencing point (0,1,..,number(points)-1); int More...
 
 __offset
 >=0: offset number at point, <0: offset number at next point with offset; int More...
 
 __type
 type (-1: first, 0: inner, 1: last) More...
 
 __p2pJacobian
 Point-to-point jacobian from previous point; matrix(float) More...
 
 __jacobians
 jacobians for propagation to previous or next point with offsets; pair(matrix(float)) More...
 
 __measurements
 measurements at point; list(GblMeasurement) More...
 
 __scatDim
 scatterer dimension (valid 2(thin) or 4(thick)) More...
 
 __scatterer
 scatterer at point: (transformation or None,) initial kinks, precision (inverse covariance matrix); list(matrix(float)) More...
 

Detailed Description

User supplied point on (initial) trajectory.

Must have jacobians for propagation to previous or next point with offsets (first, last point, points with scatterer). May have:

  1. Measurement(s) (usually 1D or 2D (or 4D or 5D))
  2. Scatterer (thin, 2D kinks)

Definition at line 182 of file gblfit.py.

Constructor & Destructor Documentation

◆ __init__()

def gblpy3.gblfit.GblPoint.__init__ (   self,
  aJacobian 
)

Create new point.

Parameters
aJacobianjacobian from previous point; matrix(float)

Definition at line 188 of file gblfit.py.

Member Function Documentation

◆ addGlobals()

def gblpy3.gblfit.GblPoint.addGlobals (   self,
  labels,
  derivatives 
)

Add global derivatives (to last measurement).

Parameters
labelsglobal labels; matrix(int)
derivativesglobal derivatives; matrix(float)

Definition at line 341 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__measurements, and gblpy3.gblfit.GblPoint.addGlobals().

Referenced by gblpy3.gblfit.GblPoint.addGlobals().

◆ addLocals()

def gblpy3.gblfit.GblPoint.addLocals (   self,
  derivatives 
)

Add local derivatives (to last measurement).

Parameters
derivativeslocal derivatives; matrix(float)

Definition at line 332 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__measurements, and gblpy3.gblfit.GblPoint.addLocals().

Referenced by gblpy3.gblfit.GblPoint.addLocals().

◆ addMeasurement()

def gblpy3.gblfit.GblPoint.addMeasurement (   self,
  aMeasurement,
  minPrecision = 0. 
)

Add a mesurement to a point.

Add measurement with arbitrary precision (inverse covariance) matrix. Will be diagonalized.

Parameters
aMeasurementmeasurement (projection (or None), residuals, precision (diagonal of or full matrix)); list(matrix(float))
minPrecisionMinimal precision to accept measurement

Definition at line 220 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__measurements.

◆ addNextJacobian()

def gblpy3.gblfit.GblPoint.addNextJacobian (   self,
  aJacobian 
)

Add jacobian to next offset.

Parameters
aJacobianjacobian for propagation to next point with offsets; matrix(float)

Definition at line 420 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__jacobians.

◆ addPrevJacobian()

def gblpy3.gblfit.GblPoint.addPrevJacobian (   self,
  aJacobian 
)

Add jacobian to previous offset.

Parameters
aJacobianjacobian for propagation to previous point with offsets; matrix(float)

Definition at line 413 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__jacobians.

◆ addScatterer()

def gblpy3.gblfit.GblPoint.addScatterer (   self,
  aScatterer 
)

Add a (thin or thick) scatterer to a point.

Add scatterer with arbitrary precision (inverse covariance) matrix. Will be diagonalized

Thin scatterer: Changes local track direction.

The (2x2) precision matrix for the local slopes is defined by the angular scattering error theta_0 and the scalar products c_1, c_2 of the offset directions in the local frame with the track direction:

       (1 - c_1*c_1 - c_2*c_2)   |  1 - c_1*c_1     - c_1*c_2  |
  P =  ~~~~~~~~~~~~~~~~~~~~~~~ * |                             |
           theta_0*theta_0       |    - c_1*c_2   1 - c_2*c_2  |

For a diagonal precision matrix at least one of the scalar products c_1, c_2 must be zero.

One offset is defined at the point. The comparison with the previous and next offsets allows to define a kink.

Thick scatterer: Local track direction and position change (correlated).

Combination of several scatterers extrapolated to point (yielding dense (4x4) covariance matrix). E.g. sum of thin scatterers (covariance V_i, jacobian J_i for slopes and offsets):

  P = [sum(J_i*V_i*J_i^t)]^-1

For each thin scatterer the covarinace matrix V_i is defined by the angular scattering error theta_0,i and the scalar products c_1,i, c_2,i of the offset directions in the local frame with the track direction:

                                    | 1 - c_2,i^2  c_1,i*c_2,i     0     0 |
                                    |                                      |
                theta_0,i^2         | c_1,i*c_2,i  1 - c_1,i^2     0     0 |
  V_i = ------------------------- * |                                      |
        (1 - c_1,i^2 - c_2,i^2)^2   |      0            0          0     0 |
                                    |                                      |
                                    |      0            0          0     0 |

Two offsets are defined at the point to describe a step (in the trajectory). The comparison with the previous and next offsets allows to define a kink.

Parameters
aScattererscatterer (kinks(+steps), precision (diagonal of or full matrix, thin:2x2, thick:4x4)); list(matrix(float))

Definition at line 284 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__scatDim, and gblpy3.gblfit.GblPoint.__scatterer.

◆ getDerivatives()

def gblpy3.gblfit.GblPoint.getDerivatives (   self,
  index 
)

Get derivatives for locally linearized track model (backward or forward propagation).

Parameters
index0 (previous) or 1 (next point with offsets); int
Returns
derivatives; list(matrix(float))

Definition at line 428 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__jacobians.

◆ getLabel()

def gblpy3.gblfit.GblPoint.getLabel (   self)

Retrieve label of a point.

Returns
label; int

Definition at line 356 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__label, and gblpy3.gblfit.GblData.__label.

◆ getMeasurements()

def gblpy3.gblfit.GblPoint.getMeasurements (   self)

Retrieve measurements of a point.

Returns
measurements; list(GblMeasurements)

Definition at line 234 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__measurements.

◆ getOffset()

def gblpy3.gblfit.GblPoint.getOffset (   self)

Get offset of a point.

Returns
offset number (at point (>=0) or at next point with offset (<0)); int

Definition at line 370 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__offset.

◆ getP2pJacobian()

def gblpy3.gblfit.GblPoint.getP2pJacobian (   self)

Retrieve jacobian of a point.

Returns
point-to-point jacobian (from previous point); matrix(float)

Definition at line 406 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__p2pJacobian.

◆ getReducedScatterer()

def gblpy3.gblfit.GblPoint.getReducedScatterer (   self)

Retrieve (reduced) scatterer of a point.

Reduce to positional scatterer.

Returns
scatterer (transformation, steps, precision); list(matrix(float))

Definition at line 315 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__scatterer.

◆ getScatDim()

def gblpy3.gblfit.GblPoint.getScatDim (   self)

Get scatterer dimension.

Returns
dimension; int

Definition at line 299 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__scatDim.

◆ getScatterer()

def gblpy3.gblfit.GblPoint.getScatterer (   self)

Retrieve scatterer of a point.

Returns
scatterer (transformation, kinks(+steps), precision); list(matrix(float))

Definition at line 306 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__scatterer.

◆ hasMeasurement()

def gblpy3.gblfit.GblPoint.hasMeasurement (   self)

Check point for a measurement.

Returns
flag; bool

Definition at line 227 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__measurements.

◆ isFirst()

def gblpy3.gblfit.GblPoint.isFirst (   self)

Is first point?

Returns
flag; bool

Definition at line 384 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__type, and gblpy3.gblfit.GblData.__type.

◆ isLast()

def gblpy3.gblfit.GblPoint.isLast (   self)

Is last point?

Returns
flag; bool

Definition at line 391 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__type, and gblpy3.gblfit.GblData.__type.

◆ printPoint()

def gblpy3.gblfit.GblPoint.printPoint (   self)

◆ setLabel()

def gblpy3.gblfit.GblPoint.setLabel (   self,
  aLabel 
)

Define label of a point.

Parameters
aLabellabel; int

Definition at line 349 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__label, and gblpy3.gblfit.GblData.__label.

◆ setOffset()

def gblpy3.gblfit.GblPoint.setOffset (   self,
  anOffset 
)

Define offset of a point and references to previous and next point with offset.

Parameters
anOffsetoffset number (at point (>=0) or at next point with offset (<0)); int

Definition at line 363 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__offset.

◆ setType()

def gblpy3.gblfit.GblPoint.setType (   self,
  aType 
)

Define type of a point.

Parameters
aTypetype (-1: first, 0: inner, 1: last); int

Definition at line 377 of file gblfit.py.

References gblpy3.gblfit.GblPoint.__type, and gblpy3.gblfit.GblData.__type.

Member Data Documentation

◆ __jacobians

gblpy3.gblfit.GblPoint.__jacobians
private

jacobians for propagation to previous or next point with offsets; pair(matrix(float))

Definition at line 198 of file gblfit.py.

Referenced by gblpy3.gblfit.GblPoint.addNextJacobian(), gblpy3.gblfit.GblPoint.addPrevJacobian(), and gblpy3.gblfit.GblPoint.getDerivatives().

◆ __label

gblpy3.gblfit.GblPoint.__label
private

◆ __measurements

gblpy3.gblfit.GblPoint.__measurements
private

◆ __offset

gblpy3.gblfit.GblPoint.__offset
private

>=0: offset number at point, <0: offset number at next point with offset; int

Definition at line 192 of file gblfit.py.

Referenced by gblpy3.gblfit.GblPoint.getOffset(), gblpy3.gblfit.GblPoint.printPoint(), and gblpy3.gblfit.GblPoint.setOffset().

◆ __p2pJacobian

gblpy3.gblfit.GblPoint.__p2pJacobian
private

Point-to-point jacobian from previous point; matrix(float)

Definition at line 196 of file gblfit.py.

Referenced by gblpy3.gblfit.GblPoint.getP2pJacobian().

◆ __scatDim

gblpy3.gblfit.GblPoint.__scatDim
private

scatterer dimension (valid 2(thin) or 4(thick))

Definition at line 202 of file gblfit.py.

Referenced by gblpy3.gblfit.GblPoint.addScatterer(), gblpy3.gblfit.GblPoint.getScatDim(), and gblpy3.gblfit.GblPoint.printPoint().

◆ __scatterer

gblpy3.gblfit.GblPoint.__scatterer
private

scatterer at point: (transformation or None,) initial kinks, precision (inverse covariance matrix); list(matrix(float))

Definition at line 204 of file gblfit.py.

Referenced by gblpy3.gblfit.GblPoint.addScatterer(), gblpy3.gblfit.GblPoint.getReducedScatterer(), and gblpy3.gblfit.GblPoint.getScatterer().

◆ __type

gblpy3.gblfit.GblPoint.__type
private

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