Belle II Software  release-08-01-10
MplTrackRep Class Reference

Monopole track representation. More...

#include <MplTrackRep.h>

Inheritance diagram for MplTrackRep:
Collaboration diagram for MplTrackRep:

Public Member Functions

 MplTrackRep (int pdgCode, float magCharge, char propDir=0)
 
double RKPropagate (M1x7 &state7, M7x7 *jacobian, M1x3 &SA, double S, bool varField=true, bool calcOnlyLastRowOfJ=false) const override
 The actual Runge Kutta propagation. More...
 
double getCharge (const StateOnPlane &state) const override
 Get the (fitted) charge of a state. More...
 
virtual AbsTrackRepclone () const override
 Clone the trackRep.
 
virtual double extrapolateToPlane (StateOnPlane &state, const SharedPlanePtr &plane, bool stopAtBoundary=false, bool calcJacobianNoise=false) const override
 Extrapolates the state to plane, and returns the extrapolation length and, via reference, the extrapolated state. More...
 
virtual double extrapolateToLine (StateOnPlane &state, const TVector3 &linePoint, const TVector3 &lineDirection, bool stopAtBoundary=false, bool calcJacobianNoise=false) const override
 Extrapolates the state to the POCA to a line, and returns the extrapolation length and, via reference, the extrapolated state. More...
 
virtual double extrapolateToLine (StateOnPlane &state, const TVector3 &linePoint, const TVector3 &lineDirection, bool stopAtBoundary=false, bool calcJacobianNoise=false) const=0
 Extrapolates the state to the POCA to a line, and returns the extrapolation length and, via reference, the extrapolated state. More...
 
virtual double extrapolateToLine (StateOnPlane &state, const TVector3 &point1, const TVector3 &point2, TVector3 &poca, TVector3 &dirInPoca, TVector3 &poca_onwire, bool stopAtBoundary=false, bool calcJacobianNoise=false) const
 Resembles the interface of GFAbsTrackRep in old versions of genfit. More...
 
virtual double extrapolateToLine (StateOnPlane &state, const TVector3 &point1, const TVector3 &point2, TVector3 &poca, TVector3 &dirInPoca, TVector3 &poca_onwire, bool stopAtBoundary=false, bool calcJacobianNoise=false) const
 Resembles the interface of GFAbsTrackRep in old versions of genfit. More...
 
virtual double extrapolateToPoint (StateOnPlane &state, const TVector3 &point, bool stopAtBoundary=false, bool calcJacobianNoise=false) const override
 Extrapolates the state to the POCA to a point, and returns the extrapolation length and, via reference, the extrapolated state. More...
 
virtual double extrapolateToPoint (StateOnPlane &state, const TVector3 &point, const TMatrixDSym &G, bool stopAtBoundary=false, bool calcJacobianNoise=false) const override
 Extrapolates the state to the POCA to a point in the metric of G, and returns the extrapolation length and, via reference, the extrapolated state. More...
 
virtual double extrapolateToCylinder (StateOnPlane &state, double radius, const TVector3 &linePoint=TVector3(0., 0., 0.), const TVector3 &lineDirection=TVector3(0., 0., 1.), bool stopAtBoundary=false, bool calcJacobianNoise=false) const override
 Extrapolates the state to the cylinder surface, and returns the extrapolation length and, via reference, the extrapolated state. More...
 
virtual double extrapolateToCone (StateOnPlane &state, double radius, const TVector3 &linePoint=TVector3(0., 0., 0.), const TVector3 &lineDirection=TVector3(0., 0., 1.), bool stopAtBoundary=false, bool calcJacobianNoise=false) const override
 Extrapolates the state to the cone surface, and returns the extrapolation length and, via reference, the extrapolated state. More...
 
virtual double extrapolateToSphere (StateOnPlane &state, double radius, const TVector3 &point=TVector3(0., 0., 0.), bool stopAtBoundary=false, bool calcJacobianNoise=false) const override
 Extrapolates the state to the sphere surface, and returns the extrapolation length and, via reference, the extrapolated state. More...
 
virtual double extrapolateBy (StateOnPlane &state, double step, bool stopAtBoundary=false, bool calcJacobianNoise=false) const override
 Extrapolates the state by step (cm) and returns the extrapolation length and, via reference, the extrapolated state. More...
 
unsigned int getDim () const override
 Get the dimension of the state vector used by the track representation.
 
virtual TVector3 getPos (const StateOnPlane &state) const override
 Get the cartesian position of a state.
 
virtual TVector3 getMom (const StateOnPlane &state) const override
 Get the cartesian momentum vector of a state.
 
virtual void getPosMom (const StateOnPlane &state, TVector3 &pos, TVector3 &mom) const override
 Get cartesian position and momentum vector of a state.
 
virtual double getMomMag (const StateOnPlane &state) const override
 get the magnitude of the momentum in GeV.
 
virtual double getMomVar (const MeasuredStateOnPlane &state) const override
 get the variance of the absolute value of the momentum .
 
virtual TMatrixDSym get6DCov (const MeasuredStateOnPlane &state) const override
 Get the 6D covariance.
 
virtual void getPosMomCov (const MeasuredStateOnPlane &state, TVector3 &pos, TVector3 &mom, TMatrixDSym &cov) const override
 Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance.
 
virtual double getQop (const StateOnPlane &state) const override
 Get charge over momentum.
 
double getSpu (const StateOnPlane &state) const
 
double getTime (const StateOnPlane &state) const override
 Get the time corresponding to the StateOnPlane. Extrapolation.
 
virtual void getForwardJacobianAndNoise (TMatrixD &jacobian, TMatrixDSym &noise, TVectorD &deltaState) const override
 Get the jacobian and noise matrix of the last extrapolation.
 
virtual void getBackwardJacobianAndNoise (TMatrixD &jacobian, TMatrixDSym &noise, TVectorD &deltaState) const override
 Get the jacobian and noise matrix of the last extrapolation if it would have been done in opposite direction.
 
std::vector< genfit::MatStepgetSteps () const override
 Get stepsizes and material properties of crossed materials of the last extrapolation.
 
virtual double getRadiationLenght () const override
 Get the accumulated X/X0 (path / radiation length) of the material crossed in the last extrapolation.
 
virtual void setPosMom (StateOnPlane &state, const TVector3 &pos, const TVector3 &mom) const override
 Set position and momentum of state.
 
virtual void setPosMom (StateOnPlane &state, const TVectorD &state6) const override
 Set position and momentum of state.
 
virtual void setPosMomErr (MeasuredStateOnPlane &state, const TVector3 &pos, const TVector3 &mom, const TVector3 &posErr, const TVector3 &momErr) const override
 Set position and momentum and error of state.
 
virtual void setPosMomCov (MeasuredStateOnPlane &state, const TVector3 &pos, const TVector3 &mom, const TMatrixDSym &cov6x6) const override
 Set position, momentum and covariance of state.
 
virtual void setPosMomCov (MeasuredStateOnPlane &state, const TVectorD &state6, const TMatrixDSym &cov6x6) const override
 Set position, momentum and covariance of state.
 
virtual void setChargeSign (StateOnPlane &state, double charge) const override
 Set the sign of the charge according to charge.
 
virtual void setQop (StateOnPlane &state, double qop) const override
 Set charge/momentum.
 
void setSpu (StateOnPlane &state, double spu) const
 
void setTime (StateOnPlane &state, double time) const override
 Set time at which the state was defined.
 
virtual bool isSameType (const AbsTrackRep *other) override
 check if other is of same type (e.g. RKTrackRep).
 
virtual bool isSame (const AbsTrackRep *other) override
 check if other is of same type (e.g. RKTrackRep) and has same pdg code.
 
double extrapolateToMeasurement (StateOnPlane &state, const AbsMeasurement *measurement, bool stopAtBoundary=false, bool calcJacobianNoise=false) const
 extrapolate to an AbsMeasurement
 
TVector3 getDir (const StateOnPlane &state) const
 Get the direction vector of a state.
 
void getPosDir (const StateOnPlane &state, TVector3 &pos, TVector3 &dir) const
 Get cartesian position and direction vector of a state.
 
virtual TVectorD get6DState (const StateOnPlane &state) const
 Get the 6D state vector (x, y, z, p_x, p_y, p_z).
 
virtual void get6DStateCov (const MeasuredStateOnPlane &state, TVectorD &stateVec, TMatrixDSym &cov) const
 Translates MeasuredStateOnPlane into 6D state vector (x, y, z, p_x, p_y, p_z) and 6x6 covariance.
 
int getPDG () const
 Get the pdg code.
 
double getPDGCharge () const
 Get the charge of the particle of the pdg code.
 
double getMass (const StateOnPlane &state) const
 Get tha particle mass in GeV/c^2.
 
char getPropDir () const
 Get propagation direction. (-1, 0, 1) -> (backward, auto, forward).
 
void calcJacobianNumerically (const genfit::StateOnPlane &origState, const genfit::SharedPlanePtr destPlane, TMatrixD &jacobian) const
 Calculate Jacobian of transportation numerically. More...
 
bool switchPDGSign ()
 try to multiply pdg code with -1. (Switch from particle to anti-particle and vice versa).
 
void setPropDir (int dir)
 Set propagation direction. (-1, 0, 1) -> (backward, auto, forward).
 
void switchPropDir ()
 Switch propagation direction. Has no effect if propDir_ is set to 0.
 
virtual void setDebugLvl (unsigned int lvl=1)
 
virtual void Print (const Option_t *="") const
 

Protected Attributes

StateOnPlane lastStartState_
 
StateOnPlane lastEndState_
 state where the last extrapolation has started
 
int pdgCode_
 Particle code.
 
char propDir_
 propagation direction (-1, 0, 1) -> (backward, auto, forward)
 
unsigned int debugLvl_
 

Private Member Functions

void initArrays () const
 
virtual double extrapToPoint (StateOnPlane &state, const TVector3 &point, const TMatrixDSym *G=nullptr, bool stopAtBoundary=false, bool calcJacobianNoise=false) const
 
void getState7 (const StateOnPlane &state, M1x7 &state7) const
 
void getState5 (StateOnPlane &state, const M1x7 &state7) const
 
void calcJ_pM_5x7 (M5x7 &J_pM, const TVector3 &U, const TVector3 &V, const M1x3 &pTilde, double spu) const
 
void transformPM6 (const MeasuredStateOnPlane &state, M6x6 &out6x6) const
 
void calcJ_Mp_7x5 (M7x5 &J_Mp, const TVector3 &U, const TVector3 &V, const TVector3 &W, const M1x3 &A) const
 
void calcForwardJacobianAndNoise (const M1x7 &startState7, const DetPlane &startPlane, const M1x7 &destState7, const DetPlane &destPlane) const
 
void transformM6P (const M6x6 &in6x6, const M1x7 &state7, MeasuredStateOnPlane &state) const
 
bool RKutta (const M1x4 &SU, const DetPlane &plane, double charge, double mass, M1x7 &state7, M7x7 *jacobianT, M1x7 *J_MMT_unprojected_lastRow, double &coveredDistance, double &flightTime, bool &checkJacProj, M7x7 &noiseProjection, StepLimits &limits, bool onlyOneStep=false, bool calcOnlyLastRowOfJ=false) const
 Propagates the particle through the magnetic field. More...
 
double estimateStep (const M1x7 &state7, const M1x4 &SU, const DetPlane &plane, const double &charge, double &relMomLoss, StepLimits &limits) const
 
TVector3 pocaOnLine (const TVector3 &linePoint, const TVector3 &lineDirection, const TVector3 &point) const
 
double Extrap (const DetPlane &startPlane, const DetPlane &destPlane, double charge, double mass, bool &isAtBoundary, M1x7 &state7, double &flightTime, bool fillExtrapSteps, TMatrixDSym *cov=nullptr, bool onlyOneStep=false, bool stopAtBoundary=false, double maxStep=1.E99) const
 Handles propagation and material effects. More...
 
void checkCache (const StateOnPlane &state, const SharedPlanePtr *plane) const
 
double momMag (const M1x7 &state7) const
 

Private Attributes

const double m_magCharge
 
const double m_mass
 
std::vector< RKStepRKSteps_
 state where the last extrapolation has ended
 
int RKStepsFXStart_
 RungeKutta steps made in the last extrapolation.
 
int RKStepsFXStop_
 
std::vector< ExtrapStepExtrapSteps_
 
TMatrixD fJacobian_
 steps made in Extrap during last extrapolation
 
TMatrixDSym fNoise_
 
bool useCache_
 
unsigned int cachePos_
 use cached RKSteps_ for extrapolation
 
StepLimits limits_
 
M7x7 noiseArray_
 
M7x7 noiseProjection_
 noise matrix of the last extrapolation
 
M7x7 J_MMT_
 

Detailed Description

Monopole track representation.

It is a minimal modification of RKTrackRep with a different equations of motion for magnetic charges.

In the current implementation the states on plane are 5-d: u, v, u', v', q/p except that q in this case is magnetic, and the monopole has no electic charge.

Definition at line 33 of file MplTrackRep.h.

Member Function Documentation

◆ calcJacobianNumerically()

void calcJacobianNumerically ( const genfit::StateOnPlane origState,
const genfit::SharedPlanePtr  destPlane,
TMatrixD &  jacobian 
) const
inherited

Calculate Jacobian of transportation numerically.

Slow but accurate. Can be used to validate (semi)analytic calculations.

Definition at line 105 of file AbsTrackRep.cc.

◆ Extrap()

double Extrap ( const DetPlane startPlane,
const DetPlane destPlane,
double  charge,
double  mass,
bool &  isAtBoundary,
M1x7 state7,
double &  flightTime,
bool  fillExtrapSteps,
TMatrixDSym *  cov = nullptr,
bool  onlyOneStep = false,
bool  stopAtBoundary = false,
double  maxStep = 1.E99 
) const
privateinherited

Handles propagation and material effects.

extrapolateToPlane(), extrapolateToPoint() and extrapolateToLine() etc. call this function. Extrap() needs a plane as an argument, hence extrapolateToPoint() and extrapolateToLine() create virtual detector planes. In this function, RKutta() is called and the resulting points and point paths are filtered so that the direction doesn't change and tiny steps are filtered out. After the propagation the material effects are called via the MaterialEffects singleton. Extrap() will loop until the plane is reached, unless the propagation fails or the maximum number of iterations is exceeded.

fNoMaterial &&

Definition at line 2311 of file RKTrackRep.cc.

2323 {
2324 
2325  static const unsigned int maxNumIt(500);
2326  unsigned int numIt(0);
2327 
2328  double coveredDistance(0.);
2329  // cppcheck-suppress unreadVariable
2330  double dqop(0.);
2331 
2332  const TVector3 W(destPlane.getNormal());
2333  M1x4 SU = {{W.X(), W.Y(), W.Z(), destPlane.distance(0., 0., 0.)}};
2334 
2335  // make SU vector point away from origin
2336  if (W*destPlane.getO() < 0) {
2337  SU[0] *= -1;
2338  SU[1] *= -1;
2339  SU[2] *= -1;
2340  }
2341 
2342 
2343  M1x7 startState7 = state7;
2344 
2345  while(true){
2346 
2347  if (debugLvl_ > 0) {
2348  debugOut << "\n============ RKTrackRep::Extrap loop nr. " << numIt << " ============\n";
2349  debugOut << "Start plane: "; startPlane.Print();
2350  debugOut << "fillExtrapSteps " << fillExtrapSteps << "\n";
2351  }
2352 
2353  if(++numIt > maxNumIt){
2354  Exception exc("RKTrackRep::Extrap ==> maximum number of iterations exceeded",__LINE__,__FILE__);
2355  exc.setFatal();
2356  throw exc;
2357  }
2358 
2359  // initialize jacobianT with unit matrix
2360  for(int i = 0; i < 7*7; ++i) J_MMT_[i] = 0;
2361  for(int i=0; i<7; ++i) J_MMT_[8*i] = 1.;
2362 
2363  M7x7* noise = nullptr;
2364  isAtBoundary = false;
2365 
2366  // propagation
2367  bool checkJacProj = false;
2368  limits_.reset();
2369  limits_.setLimit(stp_sMaxArg, maxStep-fabs(coveredDistance));
2370 
2371  M1x7 J_MMT_unprojected_lastRow = {{0, 0, 0, 0, 0, 0, 1}};
2372 
2373  if( ! RKutta(SU, destPlane, charge, mass, state7, &J_MMT_, &J_MMT_unprojected_lastRow,
2374  coveredDistance, flightTime, checkJacProj, noiseProjection_,
2375  limits_, onlyOneStep, !fillExtrapSteps) ) {
2376  Exception exc("RKTrackRep::Extrap ==> Runge Kutta propagation failed",__LINE__,__FILE__);
2377  exc.setFatal();
2378  throw exc;
2379  }
2380 
2381  bool atPlane(limits_.getLowestLimit().first == stp_plane);
2382  if (limits_.getLowestLimit().first == stp_boundary)
2383  isAtBoundary = true;
2384 
2385 
2386  if (debugLvl_ > 0) {
2387  debugOut<<"RKSteps \n";
2388  for (std::vector<RKStep>::iterator it = RKSteps_.begin(); it != RKSteps_.end(); ++it){
2389  debugOut << "stepSize = " << it->matStep_.stepSize_ << "\t";
2390  it->matStep_.material_.Print();
2391  }
2392  debugOut<<"\n";
2393  }
2394 
2395 
2396 
2397  // call MatFX
2398  if(fillExtrapSteps) {
2399  noise = &noiseArray_;
2400  for(int i = 0; i < 7*7; ++i) noiseArray_[i] = 0; // set noiseArray_ to 0
2401  }
2402 
2403  unsigned int nPoints(RKStepsFXStop_ - RKStepsFXStart_);
2404  if ( nPoints>0){
2405  // momLoss has a sign - negative loss means momentum gain
2406  double momLoss = MaterialEffects::getInstance()->effects(RKSteps_,
2408  RKStepsFXStop_,
2409  fabs(charge/state7[6]), // momentum
2410  pdgCode_,
2411  noise);
2412 
2413  RKStepsFXStart_ = RKStepsFXStop_;
2414 
2415  if (debugLvl_ > 0) {
2416  debugOut << "momLoss: " << momLoss << " GeV; relative: " << momLoss/fabs(charge/state7[6])
2417  << "; coveredDistance = " << coveredDistance << "\n";
2418  if (debugLvl_ > 1 && noise != nullptr) {
2419  debugOut << "7D noise: \n";
2420  RKTools::printDim(noise->begin(), 7, 7);
2421  }
2422  }
2423 
2424  // do momLoss only for defined 1/momentum .ne.0
2425  if(fabs(state7[6])>1.E-10) {
2426 
2427  if (debugLvl_ > 0) {
2428  debugOut << "correct state7 with dx/dqop, dy/dqop ...\n";
2429  }
2430 
2431  dqop = charge/(fabs(charge/state7[6])-momLoss) - state7[6];
2432 
2433  // Correct coveredDistance and flightTime and momLoss if checkJacProj == true
2434  // The idea is to calculate the state correction (based on the mometum loss) twice:
2435  // Once with the unprojected Jacobian (which preserves coveredDistance),
2436  // and once with the projected Jacobian (which is constrained to the plane and does NOT preserve coveredDistance).
2437  // The difference of these two corrections can then be used to calculate a correction factor.
2438  if (checkJacProj && fabs(coveredDistance) > MINSTEP) {
2439  M1x3 state7_correction_unprojected = {{0, 0, 0}};
2440  for (unsigned int i=0; i<3; ++i) {
2441  state7_correction_unprojected[i] = 0.5 * dqop * J_MMT_unprojected_lastRow[i];
2442  //debugOut << "J_MMT_unprojected_lastRow[i] " << J_MMT_unprojected_lastRow[i] << "\n";
2443  //debugOut << "state7_correction_unprojected[i] " << state7_correction_unprojected[i] << "\n";
2444  }
2445 
2446  M1x3 state7_correction_projected = {{0, 0, 0}};
2447  for (unsigned int i=0; i<3; ++i) {
2448  state7_correction_projected[i] = 0.5 * dqop * J_MMT_[6*7 + i];
2449  //debugOut << "J_MMT_[6*7 + i] " << J_MMT_[6*7 + i] << "\n";
2450  //debugOut << "state7_correction_projected[i] " << state7_correction_projected[i] << "\n";
2451  }
2452 
2453  // delta distance
2454  M1x3 delta_state = {{0, 0, 0}};
2455  for (unsigned int i=0; i<3; ++i) {
2456  delta_state[i] = state7_correction_unprojected[i] - state7_correction_projected[i];
2457  }
2458 
2459  double Dist( sqrt(delta_state[0]*delta_state[0]
2460  + delta_state[1]*delta_state[1]
2461  + delta_state[2]*delta_state[2] ) );
2462 
2463  // sign: delta * a
2464  if (delta_state[0]*state7[3] + delta_state[1]*state7[4] + delta_state[2]*state7[5] > 0)
2465  Dist *= -1.;
2466 
2467  double correctionFactor( 1. + Dist / coveredDistance );
2468  flightTime *= correctionFactor;
2469  momLoss *= correctionFactor;
2470  coveredDistance = coveredDistance + Dist;
2471 
2472  if (debugLvl_ > 0) {
2473  debugOut << "correctionFactor-1 = " << correctionFactor-1. << "; Dist = " << Dist << "\n";
2474  debugOut << "corrected momLoss: " << momLoss << " GeV; relative: " << momLoss/fabs(charge/state7[6])
2475  << "; corrected coveredDistance = " << coveredDistance << "\n";
2476  }
2477  }
2478 
2479  // correct state7 with dx/dqop, dy/dqop ... Greatly improves extrapolation accuracy!
2480  double qop( charge/(fabs(charge/state7[6])-momLoss) );
2481  dqop = qop - state7[6];
2482  state7[6] = qop;
2483 
2484  for (unsigned int i=0; i<6; ++i) {
2485  state7[i] += 0.5 * dqop * J_MMT_[6*7 + i];
2486  }
2487  // normalize direction, just to make sure
2488  double norm( 1. / sqrt(state7[3]*state7[3] + state7[4]*state7[4] + state7[5]*state7[5]) );
2489  for (unsigned int i=3; i<6; ++i)
2490  state7[i] *= norm;
2491 
2492  }
2493  } // finished MatFX
2494 
2495 
2496  // fill ExtrapSteps_
2497  if (fillExtrapSteps) {
2498  static const ExtrapStep defaultExtrapStep;
2499  ExtrapSteps_.push_back(defaultExtrapStep);
2500  std::vector<ExtrapStep>::iterator lastStep = ExtrapSteps_.end() - 1;
2501 
2502  // Store Jacobian of this step for final calculation.
2503  lastStep->jac7_ = J_MMT_;
2504 
2505  if( checkJacProj == true ){
2506  //project the noise onto the destPlane
2507  RKTools::Np_N_NpT(noiseProjection_, noiseArray_);
2508 
2509  if (debugLvl_ > 1) {
2510  debugOut << "7D noise projected onto plane: \n";
2511  RKTools::printDim(noiseArray_.begin(), 7, 7);
2512  }
2513  }
2514 
2515  // Store this step's noise for final calculation.
2516  lastStep->noise7_ = noiseArray_;
2517 
2518  if (debugLvl_ > 2) {
2519  debugOut<<"ExtrapSteps \n";
2520  for (std::vector<ExtrapStep>::iterator it = ExtrapSteps_.begin(); it != ExtrapSteps_.end(); ++it){
2521  debugOut << "7D Jacobian: "; RKTools::printDim((it->jac7_.begin()), 5,5);
2522  debugOut << "7D noise: "; RKTools::printDim((it->noise7_.begin()), 5,5);
2523  }
2524  debugOut<<"\n";
2525  }
2526  }
2527 
2528 
2529 
2530  // check if at boundary
2531  if (stopAtBoundary and isAtBoundary) {
2532  if (debugLvl_ > 0) {
2533  debugOut << "stopAtBoundary -> break; \n ";
2534  }
2535  break;
2536  }
2537 
2538  if (onlyOneStep) {
2539  if (debugLvl_ > 0) {
2540  debugOut << "onlyOneStep -> break; \n ";
2541  }
2542  break;
2543  }
2544 
2545  //break if we arrived at destPlane
2546  if(atPlane) {
2547  if (debugLvl_ > 0) {
2548  debugOut << "arrived at destPlane with a distance of " << destPlane.distance(state7[0], state7[1], state7[2]) << " cm left. ";
2549  if (destPlane.isInActive(state7[0], state7[1], state7[2], state7[3], state7[4], state7[5]))
2550  debugOut << "In active area of destPlane. \n";
2551  else
2552  debugOut << "NOT in active area of plane. \n";
2553 
2554  debugOut << " position: "; TVector3(state7[0], state7[1], state7[2]).Print();
2555  debugOut << " direction: "; TVector3(state7[3], state7[4], state7[5]).Print();
2556  }
2557  break;
2558  }
2559 
2560  }
2561 
2562  if (fillExtrapSteps) {
2563  // propagate cov and add noise
2564  calcForwardJacobianAndNoise(startState7, startPlane, state7, destPlane);
2565 
2566  if (cov != nullptr) {
2567  cov->Similarity(fJacobian_);
2568  *cov += fNoise_;
2569  }
2570 
2571  if (debugLvl_ > 0) {
2572  if (cov != nullptr) {
2573  debugOut << "final covariance matrix after Extrap: "; cov->Print();
2574  }
2575  }
2576  }
2577 
2578  return coveredDistance;
2579 }
int pdgCode_
Particle code.
Definition: AbsTrackRep.h:364
M7x7 noiseProjection_
noise matrix of the last extrapolation
Definition: RKTrackRep.h:309
std::vector< RKStep > RKSteps_
state where the last extrapolation has ended
Definition: RKTrackRep.h:294
int RKStepsFXStart_
RungeKutta steps made in the last extrapolation.
Definition: RKTrackRep.h:295
bool RKutta(const M1x4 &SU, const DetPlane &plane, double charge, double mass, M1x7 &state7, M7x7 *jacobianT, M1x7 *J_MMT_unprojected_lastRow, double &coveredDistance, double &flightTime, bool &checkJacProj, M7x7 &noiseProjection, StepLimits &limits, bool onlyOneStep=false, bool calcOnlyLastRowOfJ=false) const
Propagates the particle through the magnetic field.
Definition: RKTrackRep.cc:1801
TMatrixD fJacobian_
steps made in Extrap during last extrapolation
Definition: RKTrackRep.h:299
std::pair< StepLimitType, double > getLowestLimit(double margin=1.E-3) const
Get the lowest limit.
Definition: StepLimits.cc:44
void setLimit(StepLimitType type, double value)
absolute of value will be taken! If limit is already lower, it will be set to value anyway.
Definition: StepLimits.h:91
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
double charge(int pdgCode)
Returns electric charge of a particle with given pdg code.
Definition: EvtPDLUtil.cc:44
std::ostream debugOut
Default stream for debug output.

◆ extrapolateBy()

double extrapolateBy ( StateOnPlane state,
double  step,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
overridevirtualinherited

Extrapolates the state by step (cm) and returns the extrapolation length and, via reference, the extrapolated state.

If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.

If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.

Implements AbsTrackRep.

Definition at line 723 of file RKTrackRep.cc.

◆ extrapolateToCone()

double extrapolateToCone ( StateOnPlane state,
double  radius,
const TVector3 &  linePoint = TVector3(0., 0., 0.),
const TVector3 &  lineDirection = TVector3(0., 0., 1.),
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
overridevirtualinherited

Extrapolates the state to the cone surface, and returns the extrapolation length and, via reference, the extrapolated state.

If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.

If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.

Implements AbsTrackRep.

Definition at line 478 of file RKTrackRep.cc.

◆ extrapolateToCylinder()

double extrapolateToCylinder ( StateOnPlane state,
double  radius,
const TVector3 &  linePoint = TVector3(0., 0., 0.),
const TVector3 &  lineDirection = TVector3(0., 0., 1.),
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
overridevirtualinherited

Extrapolates the state to the cylinder surface, and returns the extrapolation length and, via reference, the extrapolated state.

If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.

If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.

Implements AbsTrackRep.

Definition at line 353 of file RKTrackRep.cc.

◆ extrapolateToLine() [1/4]

double extrapolateToLine ( StateOnPlane state,
const TVector3 &  linePoint,
const TVector3 &  lineDirection,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
overridevirtualinherited

Extrapolates the state to the POCA to a line, and returns the extrapolation length and, via reference, the extrapolated state.

If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.

If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.

Implements AbsTrackRep.

Definition at line 135 of file RKTrackRep.cc.

◆ extrapolateToLine() [2/4]

virtual double extrapolateToLine
inherited

Extrapolates the state to the POCA to a line, and returns the extrapolation length and, via reference, the extrapolated state.

If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.

If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.

◆ extrapolateToLine() [3/4]

virtual double extrapolateToLine ( StateOnPlane state,
const TVector3 &  point1,
const TVector3 &  point2,
TVector3 &  poca,
TVector3 &  dirInPoca,
TVector3 &  poca_onwire,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
inlinevirtualinherited

Resembles the interface of GFAbsTrackRep in old versions of genfit.

This interface to extrapolateToLine is intended to resemble the interface of GFAbsTrackRep in old versions of genfit and is implemented by default via the preceding function.

If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.

If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.

Definition at line 120 of file AbsTrackRep.h.

◆ extrapolateToLine() [4/4]

virtual double extrapolateToLine
inlineinherited

Resembles the interface of GFAbsTrackRep in old versions of genfit.

This interface to extrapolateToLine is intended to resemble the interface of GFAbsTrackRep in old versions of genfit and is implemented by default via the preceding function.

If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.

If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.

Definition at line 120 of file AbsTrackRep.h.

◆ extrapolateToPlane()

double extrapolateToPlane ( StateOnPlane state,
const SharedPlanePtr plane,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
overridevirtualinherited

Extrapolates the state to plane, and returns the extrapolation length and, via reference, the extrapolated state.

If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.

If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.

Implements AbsTrackRep.

Definition at line 80 of file RKTrackRep.cc.

◆ extrapolateToPoint() [1/2]

virtual double extrapolateToPoint ( StateOnPlane state,
const TVector3 &  point,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
inlineoverridevirtualinherited

Extrapolates the state to the POCA to a point, and returns the extrapolation length and, via reference, the extrapolated state.

If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.

If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.

Implements AbsTrackRep.

Definition at line 100 of file RKTrackRep.h.

103  {
104  return extrapToPoint(state, point, nullptr, stopAtBoundary, calcJacobianNoise);
105  }

◆ extrapolateToPoint() [2/2]

virtual double extrapolateToPoint ( StateOnPlane state,
const TVector3 &  point,
const TMatrixDSym &  G,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
inlineoverridevirtualinherited

Extrapolates the state to the POCA to a point in the metric of G, and returns the extrapolation length and, via reference, the extrapolated state.

If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.

If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.

Implements AbsTrackRep.

Definition at line 107 of file RKTrackRep.h.

◆ extrapolateToSphere()

double extrapolateToSphere ( StateOnPlane state,
double  radius,
const TVector3 &  point = TVector3(0., 0., 0.),
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
overridevirtualinherited

Extrapolates the state to the sphere surface, and returns the extrapolation length and, via reference, the extrapolated state.

If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.

If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.

Implements AbsTrackRep.

Definition at line 612 of file RKTrackRep.cc.

◆ getCharge()

double getCharge ( const StateOnPlane state) const
overridevirtual

Get the (fitted) charge of a state.

This is not always equal the pdg charge (e.g. if the charge sign was flipped during the fit).

Reimplemented from RKTrackRep.

Definition at line 44 of file MplTrackRep.cc.

44  {
45 
46  if (dynamic_cast<const MeasurementOnPlane*>(&state) != nullptr) {
47  Exception exc("RKTrackRep::getCharge - cannot get charge from MeasurementOnPlane",__LINE__,__FILE__);
48  exc.setFatal();
49  throw exc;
50  }
51 
52  double pdgCharge( m_magCharge * (this->getPDGCharge() > 0 ? 1.0 : -1.0));
53 
54  // return pdgCharge with sign of q/p
55  if (state.getState()(0) * pdgCharge < 0)
56  return -pdgCharge;
57  else
58  return pdgCharge;
59 }
double getPDGCharge() const
Get the charge of the particle of the pdg code.
Definition: AbsTrackRep.cc:93
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Definition: Exception.h:48
Measured coordinates on a plane.

◆ RKPropagate()

double RKPropagate ( M1x7 state7,
M7x7 jacobian,
M1x3 SA,
double  S,
bool  varField = true,
bool  calcOnlyLastRowOfJ = false 
) const
overridevirtual

The actual Runge Kutta propagation.

propagate state7 with step S. Fills SA (Start directions derivatives dA/S). This is a single Runge-Kutta step. If jacobian is nullptr, only the state is propagated, otherwise also the 7x7 jacobian is calculated. If varField is false, the magnetic field will only be evaluated at the starting position. The return value is an estimation on how good the extrapolation is, and it is usually fine if it is > 1. It gives a suggestion how you must scale S so that the quality will be sufficient.

Reimplemented from RKTrackRep.

Definition at line 61 of file MplTrackRep.cc.

◆ RKutta()

bool RKutta ( const M1x4 SU,
const DetPlane plane,
double  charge,
double  mass,
M1x7 state7,
M7x7 jacobianT,
M1x7 J_MMT_unprojected_lastRow,
double &  coveredDistance,
double &  flightTime,
bool &  checkJacProj,
M7x7 noiseProjection,
StepLimits limits,
bool  onlyOneStep = false,
bool  calcOnlyLastRowOfJ = false 
) const
privateinherited

Propagates the particle through the magnetic field.

If the propagation is successful and the plane is reached, the function returns true. Propagated state and the jacobian of the extrapolation are written to state7 and jacobianT. The jacobian is only calculated if jacobianT != nullptr. In the main loop of the Runge Kutta algorithm, the estimateStep() is called and may reduce the estimated stepsize so that a maximum momentum loss will not be exceeded, and stop at material boundaries. If this is the case, RKutta() will only propagate the reduced distance and then return. This is to ensure that material effects, which are calculated after the propagation, are taken into account properly.

Definition at line 1801 of file RKTrackRep.cc.


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