20 #include "RKTrackRep.h"
30 #include <TDecompLU.h>
35 #define MINSTEP 0.001 // minimum step [cm] for Runge Kutta and iteration to POCA
39 const bool useInvertFast =
false;
47 lastStartState_(
this),
62 lastStartState_(
this),
83 bool calcJacobianNoise)
const {
86 debugOut <<
"RKTrackRep::extrapolateToPlane()\n";
92 debugOut <<
"state is already defined at plane. Do nothing! \n";
100 M1x7 state7 = {{0, 0, 0, 0, 0, 0, 0}};
103 TMatrixDSym* covPtr(
nullptr);
104 bool fillExtrapSteps(
false);
105 if (dynamic_cast<MeasuredStateOnPlane*>(&state) !=
nullptr) {
107 fillExtrapSteps =
true;
109 else if (calcJacobianNoise)
110 fillExtrapSteps =
true;
113 bool isAtBoundary(
false);
114 double flightTime( 0. );
115 double coveredDistance(
Extrap(*(state.
getPlane()), *plane,
getCharge(state),
getMass(state), isAtBoundary, state7, flightTime, fillExtrapSteps, covPtr,
false, stopAtBoundary) );
117 if (stopAtBoundary && isAtBoundary) {
119 TVector3(state7[3], state7[4], state7[5]))));
131 return coveredDistance;
136 const TVector3& linePoint,
137 const TVector3& lineDirection,
139 bool calcJacobianNoise)
const {
142 debugOut <<
"RKTrackRep::extrapolateToLine()\n";
147 static const unsigned int maxIt(1000);
153 bool fillExtrapSteps(
false);
154 if (dynamic_cast<MeasuredStateOnPlane*>(&state) !=
nullptr) {
155 fillExtrapSteps =
true;
157 else if (calcJacobianNoise)
158 fillExtrapSteps =
true;
160 double step(0.), lastStep(0.), maxStep(1.E99), angle(0), distToPoca(0), tracklength(0);
163 double flightTime = 0;
164 TVector3 dir(state7[3], state7[4], state7[5]);
165 TVector3 lastDir(0,0,0);
166 TVector3 poca, poca_onwire;
167 bool isAtBoundary(
false);
171 unsigned int iterations(0);
174 if(++iterations == maxIt) {
175 Exception exc(
"RKTrackRep::extrapolateToLine ==> extrapolation to line failed, maximum number of iterations reached",__LINE__,__FILE__);
183 step = this->
Extrap(startPlane, *plane, charge, mass, isAtBoundary, state7, flightTime,
false,
nullptr,
true, stopAtBoundary, maxStep);
186 dir.SetXYZ(state7[3], state7[4], state7[5]);
187 poca.SetXYZ(state7[0], state7[1], state7[2]);
188 poca_onwire =
pocaOnLine(linePoint, lineDirection, poca);
191 if (stopAtBoundary && isAtBoundary) {
192 plane->setON(dir, poca);
196 angle = fabs(dir.Angle((poca_onwire-poca))-TMath::PiOver2());
197 distToPoca = (poca_onwire-poca).Mag();
198 if (angle*distToPoca < 0.1*
MINSTEP)
break;
202 if (lastStep*
step < 0){
204 maxStep = 0.5*fabs(lastStep);
208 plane->setU(dir.Cross(lineDirection));
211 if (fillExtrapSteps) {
226 debugOut <<
"RKTrackRep::extrapolateToLine(): Reached POCA after " << iterations+1 <<
" iterations. Distance: " << (poca_onwire-poca).Mag() <<
" cm. Angle deviation: " << dir.Angle((poca_onwire-poca))-TMath::PiOver2() <<
" rad \n";
236 const TVector3&
point,
237 const TMatrixDSym* G,
239 bool calcJacobianNoise)
const {
242 debugOut <<
"RKTrackRep::extrapolateToPoint()\n";
247 static const unsigned int maxIt(1000);
253 bool fillExtrapSteps(
false);
254 if (dynamic_cast<MeasuredStateOnPlane*>(&state) !=
nullptr) {
255 fillExtrapSteps =
true;
257 else if (calcJacobianNoise)
258 fillExtrapSteps =
true;
260 double step(0.), lastStep(0.), maxStep(1.E99), angle(0), distToPoca(0), tracklength(0);
261 TVector3 dir(state7[3], state7[4], state7[5]);
263 if(G->GetNrows() != 3) {
264 Exception exc(
"RKTrackRep::extrapolateToLine ==> G is not 3x3",__LINE__,__FILE__);
268 dir = TMatrix(*G) * dir;
270 TVector3 lastDir(0,0,0);
273 bool isAtBoundary(
false);
277 unsigned int iterations(0);
280 double flightTime = 0;
283 if(++iterations == maxIt) {
284 Exception exc(
"RKTrackRep::extrapolateToPoint ==> extrapolation to point failed, maximum number of iterations reached",__LINE__,__FILE__);
292 step = this->
Extrap(startPlane, *plane, charge, mass, isAtBoundary, state7, flightTime,
false,
nullptr,
true, stopAtBoundary, maxStep);
295 dir.SetXYZ(state7[3], state7[4], state7[5]);
297 dir = TMatrix(*G) * dir;
299 poca.SetXYZ(state7[0], state7[1], state7[2]);
302 if (stopAtBoundary && isAtBoundary) {
303 plane->setON(dir, poca);
307 angle = fabs(dir.Angle((point-poca))-TMath::PiOver2());
308 distToPoca = (point-poca).Mag();
309 if (angle*distToPoca < 0.1*
MINSTEP)
break;
313 if (lastStep*
step < 0){
319 maxStep = 0.5*fabs(lastStep);
323 plane->setNormal(dir);
326 if (fillExtrapSteps) {
342 debugOut <<
"RKTrackRep::extrapolateToPoint(): Reached POCA after " << iterations+1 <<
" iterations. Distance: " << (point-poca).Mag() <<
" cm. Angle deviation: " << dir.Angle((point-poca))-TMath::PiOver2() <<
" rad \n";
353 const TVector3& linePoint,
354 const TVector3& lineDirection,
356 bool calcJacobianNoise)
const {
359 debugOut <<
"RKTrackRep::extrapolateToCylinder()\n";
364 static const unsigned int maxIt(1000);
370 bool fillExtrapSteps(
false);
371 if (dynamic_cast<MeasuredStateOnPlane*>(&state) !=
nullptr) {
372 fillExtrapSteps =
true;
374 else if (calcJacobianNoise)
375 fillExtrapSteps =
true;
377 double tracklength(0.), maxStep(1.E99);
379 TVector3 dest,
pos, dir;
381 bool isAtBoundary(
false);
385 unsigned int iterations(0);
388 double flightTime = 0;
391 if(++iterations == maxIt) {
392 Exception exc(
"RKTrackRep::extrapolateToCylinder ==> maximum number of iterations reached",__LINE__,__FILE__);
397 pos.SetXYZ(state7[0], state7[1], state7[2]);
398 dir.SetXYZ(state7[3], state7[4], state7[5]);
401 TVector3 AO = (pos - linePoint);
402 TVector3 AOxAB = (AO.Cross(lineDirection));
403 TVector3 VxAB = (dir.Cross(lineDirection));
404 float ab2 = (lineDirection * lineDirection);
405 float a = (VxAB * VxAB);
406 float b = 2 * (VxAB * AOxAB);
407 float c = (AOxAB * AOxAB) - (radius*radius * ab2);
408 double arg = b*b - 4.*a*
c;
410 Exception exc(
"RKTrackRep::extrapolateToCylinder ==> cannot solve",__LINE__,__FILE__);
414 double term = sqrt(arg);
417 k1 = (-b + term)/(2.*a);
418 k2 = 2.*c/(-b + term);
421 k1 = 2.*c/(-b - term);
422 k2 = (-b - term)/(2.*a);
427 if (fabs(k2)<fabs(k))
431 debugOut <<
"RKTrackRep::extrapolateToCylinder(); k = " << k <<
"\n";
434 dest = pos + k * dir;
437 plane->setUV((dest-linePoint).Cross(lineDirection), lineDirection);
439 tracklength += this->
Extrap(startPlane, *plane, charge, mass, isAtBoundary, state7, flightTime,
false,
nullptr,
true, stopAtBoundary, maxStep);
442 if (stopAtBoundary && isAtBoundary) {
443 pos.SetXYZ(state7[0], state7[1], state7[2]);
444 dir.SetXYZ(state7[3], state7[4], state7[5]);
446 plane->setUV((pos-linePoint).Cross(lineDirection), lineDirection);
456 if (fillExtrapSteps) {
478 const TVector3& conePoint,
479 const TVector3& coneDirection,
481 bool calcJacobianNoise)
const {
484 debugOut <<
"RKTrackRep::extrapolateToCone()\n";
489 static const unsigned int maxIt(1000);
495 bool fillExtrapSteps(
false);
496 if (dynamic_cast<MeasuredStateOnPlane*>(&state) !=
nullptr) {
497 fillExtrapSteps =
true;
499 else if (calcJacobianNoise)
500 fillExtrapSteps =
true;
502 double tracklength(0.), maxStep(1.E99);
504 TVector3 dest,
pos, dir;
506 bool isAtBoundary(
false);
510 unsigned int iterations(0);
513 double flightTime = 0;
516 if(++iterations == maxIt) {
517 Exception exc(
"RKTrackRep::extrapolateToCone ==> maximum number of iterations reached",__LINE__,__FILE__);
522 pos.SetXYZ(state7[0], state7[1], state7[2]);
523 dir.SetXYZ(state7[3], state7[4], state7[5]);
530 TVector3 cDirection = coneDirection.Unit();
531 TVector3
Delta = (pos - conePoint);
532 double DirDelta = cDirection *
Delta;
534 double UDir = dir * cDirection;
535 double UDelta = dir *
Delta;
536 double U2 = dir * dir;
537 double cosAngle2 =
cos(openingAngle)*
cos(openingAngle);
538 double a = UDir*UDir - cosAngle2*U2;
539 double b = UDir*DirDelta - cosAngle2*UDelta;
540 double c = DirDelta*DirDelta - cosAngle2*
Delta2;
542 double arg = b*b - a*
c;
544 Exception exc(
"RKTrackRep::extrapolateToCone ==> cannot solve",__LINE__,__FILE__);
551 double term = sqrt(arg);
553 k1 = (-b + term) / a;
554 k2 = (-b - term) / a;
558 if(fabs(k2) < fabs(k)) {
563 debugOut <<
"RKTrackRep::extrapolateToCone(); k = " << k <<
"\n";
566 dest = pos + k * dir;
571 plane->setUV((dest-conePoint).Cross(coneDirection), dest-conePoint);
573 tracklength += this->
Extrap(startPlane, *plane, charge, mass, isAtBoundary, state7, flightTime,
false,
nullptr,
true, stopAtBoundary, maxStep);
576 if (stopAtBoundary && isAtBoundary) {
577 pos.SetXYZ(state7[0], state7[1], state7[2]);
578 dir.SetXYZ(state7[3], state7[4], state7[5]);
580 plane->setUV((pos-conePoint).Cross(coneDirection), pos-conePoint);
590 if (fillExtrapSteps) {
612 const TVector3&
point,
614 bool calcJacobianNoise)
const {
617 debugOut <<
"RKTrackRep::extrapolateToSphere()\n";
622 static const unsigned int maxIt(1000);
628 bool fillExtrapSteps(
false);
629 if (dynamic_cast<MeasuredStateOnPlane*>(&state) !=
nullptr) {
630 fillExtrapSteps =
true;
632 else if (calcJacobianNoise)
633 fillExtrapSteps =
true;
635 double tracklength(0.), maxStep(1.E99);
637 TVector3 dest,
pos, dir;
639 bool isAtBoundary(
false);
643 unsigned int iterations(0);
646 double flightTime = 0;
649 if(++iterations == maxIt) {
650 Exception exc(
"RKTrackRep::extrapolateToSphere ==> maximum number of iterations reached",__LINE__,__FILE__);
655 pos.SetXYZ(state7[0], state7[1], state7[2]);
656 dir.SetXYZ(state7[3], state7[4], state7[5]);
659 TVector3 AO = (pos -
point);
660 double dirAO = dir * AO;
661 double arg = dirAO*dirAO - AO*AO + radius*
radius;
663 Exception exc(
"RKTrackRep::extrapolateToSphere ==> cannot solve",__LINE__,__FILE__);
667 double term = sqrt(arg);
674 if (fabs(k2)<fabs(k))
678 debugOut <<
"RKTrackRep::extrapolateToSphere(); k = " << k <<
"\n";
681 dest = pos + k * dir;
683 plane->setON(dest, dest-point);
685 tracklength += this->
Extrap(startPlane, *plane, charge, mass, isAtBoundary, state7, flightTime,
false,
nullptr,
true, stopAtBoundary, maxStep);
688 if (stopAtBoundary && isAtBoundary) {
689 pos.SetXYZ(state7[0], state7[1], state7[2]);
690 dir.SetXYZ(state7[3], state7[4], state7[5]);
691 plane->setON(pos, pos-point);
701 if (fillExtrapSteps) {
724 bool calcJacobianNoise)
const {
727 debugOut <<
"RKTrackRep::extrapolateBy()\n";
732 static const unsigned int maxIt(1000);
738 bool fillExtrapSteps(
false);
739 if (dynamic_cast<MeasuredStateOnPlane*>(&state) !=
nullptr) {
740 fillExtrapSteps =
true;
742 else if (calcJacobianNoise)
743 fillExtrapSteps =
true;
745 double tracklength(0.);
747 TVector3 dest,
pos, dir;
749 bool isAtBoundary(
false);
753 unsigned int iterations(0);
756 double flightTime = 0;
759 if(++iterations == maxIt) {
760 Exception exc(
"RKTrackRep::extrapolateBy ==> maximum number of iterations reached",__LINE__,__FILE__);
765 pos.SetXYZ(state7[0], state7[1], state7[2]);
766 dir.SetXYZ(state7[3], state7[4], state7[5]);
768 dest = pos + 1.5*(step-tracklength) * dir;
770 plane->setON(dest, dir);
772 tracklength += this->
Extrap(startPlane, *plane, charge, mass, isAtBoundary, state7, flightTime,
false,
nullptr,
true, stopAtBoundary, (step-tracklength));
775 if (stopAtBoundary && isAtBoundary) {
776 pos.SetXYZ(state7[0], state7[1], state7[2]);
777 dir.SetXYZ(state7[3], state7[4], state7[5]);
778 plane->setON(pos, dir);
782 if (fabs(tracklength-step) <
MINSTEP) {
784 debugOut <<
"RKTrackRep::extrapolateBy(): reached after " << iterations <<
" iterations. \n";
786 pos.SetXYZ(state7[0], state7[1], state7[2]);
787 dir.SetXYZ(state7[3], state7[4], state7[5]);
788 plane->setON(pos, dir);
796 if (fillExtrapSteps) {
820 return TVector3(state7[0], state7[1], state7[2]);
825 M1x7 state7 = {{0, 0, 0, 0, 0, 0, 0}};
828 TVector3
mom(state7[3], state7[4], state7[5]);
835 M1x7 state7 = {{0, 0, 0, 0, 0, 0, 0}};
838 pos.SetXYZ(state7[0], state7[1], state7[2]);
839 mom.SetXYZ(state7[3], state7[4], state7[5]);
861 if (dynamic_cast<const MeasurementOnPlane*>(&state) !=
nullptr) {
862 Exception exc(
"RKTrackRep::getCharge - cannot get charge from MeasurementOnPlane",__LINE__,__FILE__);
870 if (state.
getState()(0) * pdgCharge < 0)
887 if (dynamic_cast<const MeasurementOnPlane*>(&state) !=
nullptr) {
888 Exception exc(
"RKTrackRep::getMomVar - cannot get momVar from MeasurementOnPlane",__LINE__,__FILE__);
907 if (dynamic_cast<const MeasurementOnPlane*>(&state) !=
nullptr) {
908 Exception exc(
"RKTrackRep::getSpu - cannot get spu from MeasurementOnPlane",__LINE__,__FILE__);
914 if (auxInfo.GetNrows() == 2
915 || auxInfo.GetNrows() == 1)
924 if (auxInfo.GetNrows() == 2)
932 const M1x7& destState7,
const DetPlane& destPlane)
const {
935 debugOut <<
"RKTrackRep::calcForwardJacobianAndNoise " << std::endl;
939 Exception exc(
"RKTrackRep::calcForwardJacobianAndNoise ==> cache is empty. Extrapolation must run with a MeasuredStateOnPlane.",__LINE__,__FILE__);
944 TMatrixD jac(TMatrixD::kTransposed, TMatrixD(7, 7,
ExtrapSteps_.back().jac7_.begin()));
945 TMatrixDSym noise(7,
ExtrapSteps_.back().noise7_.begin());
948 jac *= TMatrixD(TMatrixD::kTransposed, TMatrixD(7, 7,
ExtrapSteps_[i].jac7_.begin()));
952 M1x3 pTilde = {{startState7[3], startState7[4], startState7[5]}};
953 const TVector3& normal = startPlane.
getNormal();
954 double pTildeW = pTilde[0] * normal.X() + pTilde[1] * normal.Y() + pTilde[2] * normal.Z();
955 double spu = pTildeW > 0 ? 1 : -1;
956 for (
unsigned int i=0; i<3; ++i) {
957 pTilde[i] *= spu/pTildeW;
979 jacobian.ResizeTo(5,5);
986 deltaState.ResizeTo(5);
990 deltaState *= jacobian;
996 debugOut <<
"delta state : "; deltaState.Print();
1004 debugOut <<
"RKTrackRep::getBackwardJacobianAndNoise " << std::endl;
1008 Exception exc(
"RKTrackRep::getBackwardJacobianAndNoise ==> cache is empty. Extrapolation must run with a MeasuredStateOnPlane.",__LINE__,__FILE__);
1012 jacobian.ResizeTo(5,5);
1014 if (!useInvertFast) {
1015 bool status = TDecompLU::InvertLU(jacobian, 0.0);
1017 Exception e(
"cannot invert matrix, status = 0", __LINE__,__FILE__);
1023 jacobian.InvertFast(&det);
1025 Exception e(
"cannot invert matrix, status = 0", __LINE__,__FILE__);
1031 noise.ResizeTo(5,5);
1033 noise.Similarity(jacobian);
1036 deltaState.ResizeTo(5);
1046 Exception exc(
"RKTrackRep::getSteps ==> cache is empty.",__LINE__,__FILE__);
1050 std::vector<MatStep> retVal;
1053 for (
unsigned int i = 0; i<
RKSteps_.size(); ++i) {
1054 retVal.push_back(
RKSteps_[i].matStep_);
1066 Exception exc(
"RKTrackRep::getRadiationLenght ==> cache is empty.",__LINE__,__FILE__);
1072 for (
unsigned int i = 0; i<
RKSteps_.size(); ++i) {
1073 radLen +=
RKSteps_.at(i).matStep_.stepSize_ /
RKSteps_.at(i).matStep_.material_.radiationLength;
1083 if (state.
getRep() !=
this){
1084 Exception exc(
"RKTrackRep::setPosMom ==> state is defined wrt. another TrackRep",__LINE__,__FILE__);
1088 if (dynamic_cast<MeasurementOnPlane*>(&state) !=
nullptr) {
1089 Exception exc(
"RKTrackRep::setPosMom - cannot set pos/mom of a MeasurementOnPlane",__LINE__,__FILE__);
1094 if (mom.Mag2() == 0) {
1095 Exception exc(
"RKTrackRep::setPosMom - momentum is 0",__LINE__,__FILE__);
1102 if (auxInfo.GetNrows() != 2) {
1103 bool alreadySet = auxInfo.GetNrows() == 1;
1104 auxInfo.ResizeTo(2);
1113 state7[0] = pos.X();
1114 state7[1] = pos.Y();
1115 state7[2] = pos.Z();
1117 state7[3] = mom.X();
1118 state7[4] = mom.Y();
1119 state7[5] = mom.Z();
1122 double norm = 1. / sqrt(state7[3]*state7[3] + state7[4]*state7[4] + state7[5]*state7[5]);
1123 for (
unsigned int i=3; i<6; ++i)
1137 TVectorD& state5(state.
getState());
1152 if (state6.GetNrows()!=6){
1153 Exception exc(
"RKTrackRep::setPosMom ==> state has to be 6d (x, y, z, px, py, pz)",__LINE__,__FILE__);
1156 setPosMom(state, TVector3(state6(0), state6(1), state6(2)), TVector3(state6(3), state6(4), state6(5)));
1166 const TVector3& U(state.
getPlane()->getU());
1167 const TVector3& V(state.
getPlane()->getV());
1168 TVector3 W(state.
getPlane()->getNormal());
1170 double pw = mom * W;
1171 double pu = mom * U;
1172 double pv = mom * V;
1176 cov(0,0) = pow(
getCharge(state), 2) / pow(mom.Mag(), 6) *
1177 (mom.X()*mom.X() * momErr.X()*momErr.X()+
1178 mom.Y()*mom.Y() * momErr.Y()*momErr.Y()+
1179 mom.Z()*mom.Z() * momErr.Z()*momErr.Z());
1181 cov(1,1) = pow((U.X()/pw - W.X()*pu/(pw*pw)),2.) * momErr.X()*momErr.X() +
1182 pow((U.Y()/pw - W.Y()*pu/(pw*pw)),2.) * momErr.Y()*momErr.Y() +
1183 pow((U.Z()/pw - W.Z()*pu/(pw*pw)),2.) * momErr.Z()*momErr.Z();
1185 cov(2,2) = pow((V.X()/pw - W.X()*pv/(pw*pw)),2.) * momErr.X()*momErr.X() +
1186 pow((V.Y()/pw - W.Y()*pv/(pw*pw)),2.) * momErr.Y()*momErr.Y() +
1187 pow((V.Z()/pw - W.Z()*pv/(pw*pw)),2.) * momErr.Z()*momErr.Z();
1189 cov(3,3) = posErr.X()*posErr.X() * U.X()*U.X() +
1190 posErr.Y()*posErr.Y() * U.Y()*U.Y() +
1191 posErr.Z()*posErr.Z() * U.Z()*U.Z();
1193 cov(4,4) = posErr.X()*posErr.X() * V.X()*V.X() +
1194 posErr.Y()*posErr.Y() * V.Y()*V.Y() +
1195 posErr.Z()*posErr.Z() * V.Z()*V.Z();
1204 if (cov6x6.GetNcols()!=6 || cov6x6.GetNrows()!=6){
1205 Exception exc(
"RKTrackRep::setPosMomCov ==> cov has to be 6x6 (x, y, z, px, py, pz)",__LINE__,__FILE__);
1214 const M6x6& cov6x6_( *((
M6x6*) cov6x6.GetMatrixArray()) );
1222 if (state6.GetNrows()!=6){
1223 Exception exc(
"RKTrackRep::setPosMomCov ==> state has to be 6d (x, y, z, px, py, pz)",__LINE__,__FILE__);
1227 if (cov6x6.GetNcols()!=6 || cov6x6.GetNrows()!=6){
1228 Exception exc(
"RKTrackRep::setPosMomCov ==> cov has to be 6x6 (x, y, z, px, py, pz)",__LINE__,__FILE__);
1232 TVector3
pos(state6(0), state6(1), state6(2));
1233 TVector3
mom(state6(3), state6(4), state6(5));
1239 const M6x6& cov6x6_( *((
M6x6*) cov6x6.GetMatrixArray()) );
1248 if (dynamic_cast<MeasurementOnPlane*>(&state) !=
nullptr) {
1249 Exception exc(
"RKTrackRep::setChargeSign - cannot set charge of a MeasurementOnPlane",__LINE__,__FILE__);
1254 if (state.
getState()(0) * charge < 0) {
1277 bool calcOnlyLastRowOfJ)
const {
1292 static const double EC ( 0.000149896229 );
1293 static const double P3 ( 1./3. );
1294 static const double DLT ( .0002 );
1298 double S3(0), S4(0), PS2(0);
1299 M1x3 H0 = {{0.,0.,0.}}, H1 = {{0.,0.,0.}}, H2 = {{0.,0.,0.}};
1300 M1x3 r = {{0.,0.,0.}};
1302 double A0(0), A1(0),
A2(0), A3(0), A4(0), A5(0), A6(0);
1303 double B0(0), B1(0), B2(0), B3(0), B4(0), B5(0), B6(0);
1304 double C0(0), C1(0), C2(0), C3(0), C4(0), C5(0), C6(0);
1311 PS2 = state7[6]*EC * S;
1314 r[0] = R[0]; r[1] = R[1]; r[2]=R[2];
1316 H0[0] *= PS2; H0[1] *= PS2; H0[2] *= PS2;
1317 A0 = A[1]*H0[2]-A[2]*H0[1]; B0 = A[2]*H0[0]-A[0]*H0[2]; C0 = A[0]*H0[1]-A[1]*H0[0];
1318 A2 = A[0]+A0 ; B2 = A[1]+B0 ; C2 = A[2]+C0 ;
1319 A1 =
A2+A[0] ; B1 = B2+A[1] ; C1 = C2+A[2] ;
1323 r[0] += A1*S4; r[1] += B1*S4; r[2] += C1*S4;
1325 H1[0] *= PS2; H1[1] *= PS2; H1[2] *= PS2;
1328 A3 = B2*H1[2]-C2*H1[1]+A[0]; B3 = C2*H1[0]-
A2*H1[2]+A[1]; C3 =
A2*H1[1]-B2*H1[0]+A[2];
1329 A4 = B3*H1[2]-C3*H1[1]+A[0]; B4 = C3*H1[0]-A3*H1[2]+A[1]; C4 = A3*H1[1]-B3*H1[0]+A[2];
1330 A5 = A4-A[0]+A4 ; B5 = B4-A[1]+B4 ; C5 = C4-A[2]+C4 ;
1334 r[0]=R[0]+S*A4; r[1]=R[1]+S*B4; r[2]=R[2]+S*C4;
1336 H2[0] *= PS2; H2[1] *= PS2; H2[2] *= PS2;
1339 A6 = B5*H2[2]-C5*H2[1]; B6 = C5*H2[0]-A5*H2[2]; C6 = A5*H2[1]-B5*H2[0];
1345 if(jacobianT !=
nullptr){
1355 M7x7& J = *jacobianT;
1357 double dA0(0), dA2(0), dA3(0), dA4(0), dA5(0), dA6(0);
1358 double dB0(0), dB2(0), dB3(0), dB4(0), dB5(0), dB6(0);
1359 double dC0(0), dC2(0), dC3(0), dC4(0), dC5(0), dC6(0);
1363 if (!calcOnlyLastRowOfJ) {
1367 J(0, 0) = 1; J(1, 1) = 1; J(2, 2) = 1;
1373 for(
int i=start; i<6; ++i) {
1376 dA0 = H0[2]*J(i, 4)-H0[1]*J(i, 5);
1377 dB0 = H0[0]*J(i, 5)-H0[2]*J(i, 3);
1378 dC0 = H0[1]*J(i, 3)-H0[0]*J(i, 4);
1385 dA3 = J(i, 3)+dB2*H1[2]-dC2*H1[1];
1386 dB3 = J(i, 4)+dC2*H1[0]-dA2*H1[2];
1387 dC3 = J(i, 5)+dA2*H1[1]-dB2*H1[0];
1389 dA4 = J(i, 3)+dB3*H1[2]-dC3*H1[1];
1390 dB4 = J(i, 4)+dC3*H1[0]-dA3*H1[2];
1391 dC4 = J(i, 5)+dA3*H1[1]-dB3*H1[0];
1394 dA5 = dA4+dA4-J(i, 3);
1395 dB5 = dB4+dB4-J(i, 4);
1396 dC5 = dC4+dC4-J(i, 5);
1398 dA6 = dB5*H2[2]-dC5*H2[1];
1399 dB6 = dC5*H2[0]-dA5*H2[2];
1400 dC6 = dA5*H2[1]-dB5*H2[0];
1403 J(i, 0) += (dA2+dA3+dA4)*S3; J(i, 3) = ((dA0+2.*dA3)+(dA5+dA6))*P3;
1404 J(i, 1) += (dB2+dB3+dB4)*S3; J(i, 4) = ((dB0+2.*dB3)+(dB5+dB6))*P3;
1405 J(i, 2) += (dC2+dC3+dC4)*S3; J(i, 5) = ((dC0+2.*dC3)+(dC5+dC6))*P3;
1410 J(6, 3) *= state7[6]; J(6, 4) *= state7[6]; J(6, 5) *= state7[6];
1413 dA0 = H0[2]*J(6, 4)-H0[1]*J(6, 5) + A0;
1414 dB0 = H0[0]*J(6, 5)-H0[2]*J(6, 3) + B0;
1415 dC0 = H0[1]*J(6, 3)-H0[0]*J(6, 4) + C0;
1422 dA3 = J(6, 3)+dB2*H1[2]-dC2*H1[1] + (A3-A[0]);
1423 dB3 = J(6, 4)+dC2*H1[0]-dA2*H1[2] + (B3-A[1]);
1424 dC3 = J(6, 5)+dA2*H1[1]-dB2*H1[0] + (C3-A[2]);
1426 dA4 = J(6, 3)+dB3*H1[2]-dC3*H1[1] + (A4-A[0]);
1427 dB4 = J(6, 4)+dC3*H1[0]-dA3*H1[2] + (B4-A[1]);
1428 dC4 = J(6, 5)+dA3*H1[1]-dB3*H1[0] + (C4-A[2]);
1431 dA5 = dA4+dA4-J(6, 3);
1432 dB5 = dB4+dB4-J(6, 4);
1433 dC5 = dC4+dC4-J(6, 5);
1435 dA6 = dB5*H2[2]-dC5*H2[1] + A6;
1436 dB6 = dC5*H2[0]-dA5*H2[2] + B6;
1437 dC6 = dA5*H2[1]-dB5*H2[0] + C6;
1440 J(6, 0) += (dA2+dA3+dA4)*S3/state7[6]; J(6, 3) = ((dA0+2.*dA3)+(dA5+dA6))*P3/state7[6];
1441 J(6, 1) += (dB2+dB3+dB4)*S3/state7[6]; J(6, 4) = ((dB0+2.*dB3)+(dB5+dB6))*P3/state7[6];
1442 J(6, 2) += (dC2+dC3+dC4)*S3/state7[6]; J(6, 5) = ((dC0+2.*dC3)+(dC5+dC6))*P3/state7[6];
1449 R[0] += (
A2+A3+A4)*S3; A[0] += (SA[0]=((A0+2.*A3)+(A5+A6))*P3-A[0]);
1450 R[1] += (B2+B3+B4)*S3; A[1] += (SA[1]=((B0+2.*B3)+(B5+B6))*P3-A[1]);
1451 R[2] += (C2+C3+C4)*S3; A[2] += (SA[2]=((C0+2.*C3)+(C5+C6))*P3-A[2]);
1454 double CBA ( 1./sqrt(A[0]*A[0]+A[1]*A[1]+A[2]*A[2]) );
1455 A[0] *= CBA; A[1] *= CBA; A[2] *= CBA;
1459 double EST ( fabs((A1+A6)-(A3+A4)) +
1460 fabs((B1+B6)-(B3+B4)) +
1461 fabs((C1+C6)-(C3+C4)) );
1463 debugOut <<
" RKTrackRep::RKPropagate. Step = "<< S <<
"; quality EST = " << EST <<
" \n";
1473 return pow(DLT/EST, 1./5.);
1481 for (
unsigned int i=0; i<7; ++i)
1499 if (dynamic_cast<const MeasurementOnPlane*>(&state) !=
nullptr) {
1500 Exception exc(
"RKTrackRep::getState7 - cannot get pos or mom from a MeasurementOnPlane",__LINE__,__FILE__);
1505 const TVector3& U(state.
getPlane()->getU());
1506 const TVector3& V(state.
getPlane()->getV());
1507 const TVector3&
O(state.
getPlane()->getO());
1508 const TVector3& W(state.
getPlane()->getNormal());
1510 assert(state.
getState().GetNrows() == 5);
1511 const double* state5 = state.
getState().GetMatrixArray();
1513 double spu =
getSpu(state);
1515 state7[0] =
O.X() + state5[3]*U.X() + state5[4]*V.X();
1516 state7[1] =
O.Y() + state5[3]*U.Y() + state5[4]*V.Y();
1517 state7[2] =
O.Z() + state5[3]*U.Z() + state5[4]*V.Z();
1519 state7[3] = spu * (W.X() + state5[1]*U.X() + state5[2]*V.X());
1520 state7[4] = spu * (W.Y() + state5[1]*U.Y() + state5[2]*V.Y());
1521 state7[5] = spu * (W.Z() + state5[1]*U.Z() + state5[2]*V.Z());
1524 double norm = 1. / sqrt(state7[3]*state7[3] + state7[4]*state7[4] + state7[5]*state7[5]);
1525 for (
unsigned int i=3; i<6; ++i) state7[i] *= norm;
1527 state7[6] = state5[0];
1537 const TVector3&
O(state.
getPlane()->getO());
1538 const TVector3& U(state.
getPlane()->getU());
1539 const TVector3& V(state.
getPlane()->getV());
1540 const TVector3& W(state.
getPlane()->getNormal());
1543 double AtW( state7[3]*W.X() + state7[4]*W.Y() + state7[5]*W.Z() );
1550 double* state5 = state.
getState().GetMatrixArray();
1552 state5[0] = state7[6];
1553 state5[1] = (state7[3]*U.X() + state7[4]*U.Y() + state7[5]*U.Z()) / AtW;
1554 state5[2] = (state7[3]*V.X() + state7[4]*V.Y() + state7[5]*V.Z()) / AtW;
1555 state5[3] = ((state7[0]-
O.X())*U.X() +
1556 (state7[1]-
O.Y())*U.Y() +
1557 (state7[2]-
O.Z())*U.Z());
1558 state5[4] = ((state7[0]-
O.X())*V.X() +
1559 (state7[1]-
O.Y())*V.Y() +
1560 (state7[2]-
O.Z())*V.Z());
1575 std::fill(J_pM.
begin(), J_pM.
end(), 0);
1577 const double pTildeMag = sqrt(pTilde[0]*pTilde[0] + pTilde[1]*pTilde[1] + pTilde[2]*pTilde[2]);
1578 const double pTildeMag2 = pTildeMag*pTildeMag;
1580 const double utpTildeOverpTildeMag2 = (U.X()*pTilde[0] + U.Y()*pTilde[1] + U.Z()*pTilde[2]) / pTildeMag2;
1581 const double vtpTildeOverpTildeMag2 = (V.X()*pTilde[0] + V.Y()*pTilde[1] + V.Z()*pTilde[2]) / pTildeMag2;
1596 double fact = spu / pTildeMag;
1597 J_pM[10] = fact * ( U.X() - pTilde[0]*utpTildeOverpTildeMag2 );
1598 J_pM[11] = fact * ( U.Y() - pTilde[1]*utpTildeOverpTildeMag2 );
1599 J_pM[12] = fact * ( U.Z() - pTilde[2]*utpTildeOverpTildeMag2 );
1601 J_pM[17] = fact * ( V.X() - pTilde[0]*vtpTildeOverpTildeMag2 );
1602 J_pM[18] = fact * ( V.Y() - pTilde[1]*vtpTildeOverpTildeMag2 );
1603 J_pM[19] = fact * ( V.Z() - pTilde[2]*vtpTildeOverpTildeMag2 );
1612 M6x6& out6x6)
const {
1615 const TVector3& U(state.
getPlane()->getU());
1616 const TVector3& V(state.
getPlane()->getV());
1617 const TVector3& W(state.
getPlane()->getNormal());
1619 const TVectorD& state5(state.
getState());
1620 double spu(
getSpu(state));
1623 pTilde[0] = spu * (W.X() + state5(1)*U.X() + state5(2)*V.X());
1624 pTilde[1] = spu * (W.Y() + state5(1)*U.Y() + state5(2)*V.Y());
1625 pTilde[2] = spu * (W.Z() + state5(1)*U.Z() + state5(2)*V.Z());
1627 const double pTildeMag = sqrt(pTilde[0]*pTilde[0] + pTilde[1]*pTilde[1] + pTilde[2]*pTilde[2]);
1628 const double pTildeMag2 = pTildeMag*pTildeMag;
1630 const double utpTildeOverpTildeMag2 = (U.X()*pTilde[0] + U.Y()*pTilde[1] + U.Z()*pTilde[2]) / pTildeMag2;
1631 const double vtpTildeOverpTildeMag2 = (V.X()*pTilde[0] + V.Y()*pTilde[1] + V.Z()*pTilde[2]) / pTildeMag2;
1635 const double qop = state5(0);
1639 std::fill(J_pM_5x6.
begin(), J_pM_5x6.
end(), 0);
1642 double fact = -1. * p / (pTildeMag * qop);
1643 J_pM_5x6[3] = fact * pTilde[0];
1644 J_pM_5x6[4] = fact * pTilde[1];
1645 J_pM_5x6[5] = fact * pTilde[2];
1647 fact = p * spu / pTildeMag;
1648 J_pM_5x6[9] = fact * ( U.X() - pTilde[0]*utpTildeOverpTildeMag2 );
1649 J_pM_5x6[10] = fact * ( U.Y() - pTilde[1]*utpTildeOverpTildeMag2 );
1650 J_pM_5x6[11] = fact * ( U.Z() - pTilde[2]*utpTildeOverpTildeMag2 );
1652 J_pM_5x6[15] = fact * ( V.X() - pTilde[0]*vtpTildeOverpTildeMag2 );
1653 J_pM_5x6[16] = fact * ( V.Y() - pTilde[1]*vtpTildeOverpTildeMag2 );
1654 J_pM_5x6[17] = fact * ( V.Z() - pTilde[2]*vtpTildeOverpTildeMag2 );
1656 J_pM_5x6[18] = U.X();
1657 J_pM_5x6[19] = U.Y();
1658 J_pM_5x6[20] = U.Z();
1660 J_pM_5x6[24] = V.X();
1661 J_pM_5x6[25] = V.Y();
1662 J_pM_5x6[26] = V.Z();
1682 std::fill(J_Mp.
begin(), J_Mp.
end(), 0);
1684 const double AtU = A[0]*U.X() + A[1]*U.Y() + A[2]*U.Z();
1685 const double AtV = A[0]*V.X() + A[1]*V.Y() + A[2]*V.Z();
1686 const double AtW = A[0]*W.X() + A[1]*W.Y() + A[2]*W.Z();
1691 double fact = 1./(AtW*AtW);
1692 J_Mp[16] = fact * (U.X()*AtW - W.X()*AtU);
1693 J_Mp[21] = fact * (U.Y()*AtW - W.Y()*AtU);
1694 J_Mp[26] = fact * (U.Z()*AtW - W.Z()*AtU);
1696 J_Mp[17] = fact * (V.X()*AtW - W.X()*AtV);
1697 J_Mp[22] = fact * (V.Y()*AtW - W.Y()*AtV);
1698 J_Mp[27] = fact * (V.Z()*AtW - W.Z()*AtV);
1722 const TVector3& U(state.
getPlane()->getU());
1723 const TVector3& V(state.
getPlane()->getV());
1724 const TVector3& W(state.
getPlane()->getNormal());
1726 const double AtU = state7[3]*U.X() + state7[4]*U.Y() + state7[5]*U.Z();
1727 const double AtV = state7[3]*V.X() + state7[4]*V.Y() + state7[5]*V.Z();
1728 const double AtW = state7[3]*W.X() + state7[4]*W.Y() + state7[5]*W.Z();
1732 const double qop = state7[6];
1736 std::fill(J_Mp_6x5.
begin(), J_Mp_6x5.
end(), 0);
1739 J_Mp_6x5[3] = U.X();
1740 J_Mp_6x5[8] = U.Y();
1741 J_Mp_6x5[13] = U.Z();
1743 J_Mp_6x5[4] = V.X();
1744 J_Mp_6x5[9] = V.Y();
1745 J_Mp_6x5[14] = V.Z();
1747 double fact = (-1.) * qop /
p;
1748 J_Mp_6x5[15] = fact * state7[3];
1749 J_Mp_6x5[20] = fact * state7[4];
1750 J_Mp_6x5[25] = fact * state7[5];
1752 fact = 1./(p*AtW*AtW);
1753 J_Mp_6x5[16] = fact * (U.X()*AtW - W.X()*AtU);
1754 J_Mp_6x5[21] = fact * (U.Y()*AtW - W.Y()*AtU);
1755 J_Mp_6x5[26] = fact * (U.Z()*AtW - W.Z()*AtU);
1757 J_Mp_6x5[17] = fact * (V.X()*AtW - W.X()*AtV);
1758 J_Mp_6x5[22] = fact * (V.Y()*AtW - W.Y()*AtV);
1759 J_Mp_6x5[27] = fact * (V.Z()*AtW - W.Z()*AtV);
1805 M1x7* J_MMT_unprojected_lastRow,
1806 double& coveredDistance,
1809 M7x7& noiseProjection,
1812 bool calcOnlyLastRowOfJ)
const {
1815 static const double Wmax ( 3000. );
1816 static const double AngleMax ( 6.3 );
1817 static const double Pmin ( 4.E-3 );
1818 static const unsigned int maxNumIt ( 1000 );
1822 M1x3 SA = {{0.,0.,0.}};
1824 double momentum ( fabs(charge/state7[6]) );
1825 double relMomLoss ( 0 );
1826 double deltaAngle ( 0. );
1827 double An(0), S(0), Sl(0), CBA(0);
1831 debugOut <<
"RKTrackRep::RKutta \n";
1832 debugOut <<
"position: "; TVector3(R[0], R[1], R[2]).Print();
1833 debugOut <<
"direction: "; TVector3(A[0], A[1], A[2]).Print();
1834 debugOut <<
"momentum: " << momentum <<
" GeV\n";
1838 checkJacProj =
false;
1841 if(momentum < Pmin){
1842 std::ostringstream sstream;
1843 sstream <<
"RKTrackRep::RKutta ==> momentum too low: " << momentum*1000. <<
" MeV";
1844 Exception exc(sstream.str(),__LINE__,__FILE__);
1849 unsigned int counter(0);
1852 S =
estimateStep(state7, SU, plane, charge, relMomLoss, limits);
1857 while (fabs(S) >=
MINSTEP || counter == 0) {
1859 if(++counter > maxNumIt){
1860 Exception exc(
"RKTrackRep::RKutta ==> maximum number of iterations exceeded",__LINE__,__FILE__);
1866 debugOut <<
"------ RKutta main loop nr. " << counter-1 <<
" ------\n";
1869 M1x3 ABefore = {{ A[0], A[1], A[2] }};
1870 RKPropagate(state7, jacobianT, SA, S,
true, calcOnlyLastRowOfJ);
1873 coveredDistance += S;
1876 double beta = 1/hypot(1, mass*state7[6]/charge);
1877 flightTime += S / beta / 29.9792458;
1881 std::ostringstream sstream;
1882 sstream <<
"RKTrackRep::RKutta ==> Total extrapolation length is longer than length limit : " << Way <<
" cm !";
1883 Exception exc(sstream.str(),__LINE__,__FILE__);
1888 if (onlyOneStep)
return(
true);
1893 debugOut<<
" momLossExceeded -> return(true); \n";
1901 debugOut<<
" at boundary -> return(true); \n";
1913 S =
estimateStep(state7, SU, plane, charge, relMomLoss, limits);
1918 debugOut<<
" (at Plane && fabs(S) < MINSTEP) -> break and do linear extrapolation \n";
1925 debugOut<<
" (momLossExceeded && fabs(S) < MINSTEP) -> return(true), no linear extrapolation; \n";
1933 double arg = ABefore[0]*A[0] + ABefore[1]*A[1] + ABefore[2]*A[2];
1934 arg = arg > 1 ? 1 : arg;
1935 arg = arg < -1 ? -1 : arg;
1936 deltaAngle += acos(arg);
1937 if (fabs(deltaAngle) > AngleMax){
1938 std::ostringstream sstream;
1939 sstream <<
"RKTrackRep::RKutta ==> Do not get to an active plane! Already extrapolated " << deltaAngle * 180 / TMath::Pi() <<
"°.";
1940 Exception exc(sstream.str(),__LINE__,__FILE__);
1947 if (S *
RKSteps_.at(counter-1).matStep_.stepSize_ < 0 &&
1948 RKSteps_.at(counter-1).matStep_.stepSize_*
RKSteps_.at(counter-2).matStep_.stepSize_ < 0 &&
1949 RKSteps_.at(counter-2).matStep_.stepSize_*
RKSteps_.at(counter-3).matStep_.stepSize_ < 0){
1950 Exception exc(
"RKTrackRep::RKutta ==> Do not get closer to plane!",__LINE__,__FILE__);
1964 if (fabs(Sl) > 0.001*
MINSTEP){
1966 debugOut <<
" RKutta - linear extrapolation to surface\n";
1971 SA[0]*=Sl; SA[1]*=Sl; SA[2]*=Sl;
1979 CBA = 1./sqrt(A[0]*A[0]+A[1]*A[1]+A[2]*A[2]);
1980 A[0] *= CBA; A[1] *= CBA; A[2] *= CBA;
1982 R[0] += S*(A[0]-0.5*S*SA[0]);
1983 R[1] += S*(A[1]-0.5*S*SA[1]);
1984 R[2] += S*(A[2]-0.5*S*SA[2]);
1987 coveredDistance += S;
1990 double beta = 1/hypot(1, mass*state7[6]/charge);
1991 flightTime += S / beta / 29.9792458;
1994 debugOut <<
" RKutta - last stepsize too small -> can't do linear extrapolation! \n";
2000 if (jacobianT !=
nullptr) {
2011 if (checkJacProj &&
RKSteps_.size()>0){
2012 Exception exc(
"RKTrackRep::Extrap ==> covariance is projected onto destination plane again",__LINE__,__FILE__);
2019 debugOut <<
" Project Jacobian of extrapolation onto destination plane\n";
2021 An = A[0]*SU[0] + A[1]*SU[1] + A[2]*SU[2];
2022 An = (fabs(An) > 1.E-7 ? 1./An : 0);
2025 if (calcOnlyLastRowOfJ)
2028 double* jacPtr = jacobianT->
begin();
2030 for(
unsigned int j=42; j<49; j+=7) {
2031 (*J_MMT_unprojected_lastRow)[j-42] = jacPtr[j];
2035 norm = (jacPtr[i]*SU[0] + jacPtr[i+1]*SU[1] + jacPtr[i+2]*SU[2]) * An;
2036 jacPtr[i] -= norm*A [0]; jacPtr[i+1] -= norm*A [1]; jacPtr[i+2] -= norm*A [2];
2037 jacPtr[i+3] -= norm*SA[0]; jacPtr[i+4] -= norm*SA[1]; jacPtr[i+5] -= norm*SA[2];
2039 checkJacProj =
true;
2047 if (!calcOnlyLastRowOfJ) {
2048 for (
int iRow = 0; iRow < 3; ++iRow) {
2049 for (
int iCol = 0; iCol < 3; ++iCol) {
2050 noiseProjection[iRow*7 + iCol] = (iRow == iCol) - An * SU[iCol] * A[iRow];
2051 noiseProjection[(iRow + 3)*7 + iCol] = - An * SU[iCol] * SA[iRow];
2105 debugOut <<
" RKTrackRep::estimateStep \n";
2106 debugOut <<
" position: "; TVector3(state7[0], state7[1], state7[2]).Print();
2107 debugOut <<
" direction: "; TVector3(state7[3], state7[4], state7[5]).Print();
2111 double Dist ( SU[3] - (state7[0]*SU[0] +
2114 double An ( state7[3]*SU[0] +
2119 if (fabs(An) > 1.E-10)
2122 SLDist = Dist*1.E10;
2123 if (An<0) SLDist *= -1.;
2130 debugOut <<
" Distance to plane: " << Dist <<
"\n";
2131 debugOut <<
" SL distance to plane: " << SLDist <<
"\n";
2133 debugOut <<
" Direction is pointing towards surface.\n";
2135 debugOut <<
" Direction is pointing away from surface.\n";
2144 std::pair<double, double> distVsStep (9.E99, 9.E99);
2146 static const unsigned int maxNumIt = 10;
2147 unsigned int counter(0);
2149 while (fabs(fieldCurvLimit) >
MINSTEP) {
2151 if(++counter > maxNumIt){
2155 fieldCurvLimit *= 0.5;
2159 M1x7 state7_temp = state7;
2160 M1x3 SA = {{0, 0, 0}};
2162 double q (
RKPropagate(state7_temp,
nullptr, SA, fieldCurvLimit,
true) );
2164 debugOut <<
" maxStepArg = " << fieldCurvLimit <<
"; q = " << q <<
" \n";
2169 Dist = SU[3] - (state7_temp[0] * SU[0] +
2170 state7_temp[1] * SU[1] +
2171 state7_temp[2] * SU[2]);
2173 An = state7_temp[3] * SU[0] +
2174 state7_temp[4] * SU[1] +
2175 state7_temp[5] * SU[2];
2177 if (fabs(Dist/An) < fabs(distVsStep.first)) {
2178 distVsStep.first = Dist/An;
2179 distVsStep.second = fieldCurvLimit;
2186 fieldCurvLimit *= 2;
2190 fieldCurvLimit *= q * 0.95;
2192 if (fabs(q-1) < 0.25 ||
2196 if (fabs(fieldCurvLimit) <
MINSTEP)
2202 if (fabs(distVsStep.first) < 8.E99) {
2203 stepToPlane = distVsStep.first + distVsStep.second;
2214 debugOut <<
" auto select direction";
2222 debugOut <<
" straight line approximation is fine.\n";
2226 if( plane.
isInActive(state7[0], state7[1], state7[2], state7[3], state7[4], state7[5]) ) {
2228 debugOut <<
" direction is pointing to active part of surface. \n";
2236 debugOut <<
" we are near the plane, but not pointing to the active area. make a big step! \n";
2246 debugOut <<
" invert Step according to propDir_ and make a big step. \n";
2253 static const RKStep defaultRKStep;
2254 RKSteps_.push_back( defaultRKStep );
2255 std::vector<RKStep>::iterator lastStep =
RKSteps_.end() - 1;
2256 lastStep->state7_ = state7;
2260 M1x7 state7_temp = {{ state7[0], state7[1], state7[2], state7[3], state7[4], state7[5], state7[6] }};
2267 lastStep->matStep_.material_,
2272 lastStep->matStep_.material_ = (lastStep - 1)->matStep_.material_;
2283 lastStep->matStep_.stepSize_ = finalStep;
2284 lastStep->limits_ = limits;
2287 debugOut <<
" --> Step to be used: " << finalStep <<
"\n";
2297 TVector3 retVal(lineDirection);
2299 double t = 1./(retVal.Mag2()) * ((point*retVal) - (linePoint*retVal));
2301 retVal += linePoint;
2314 bool fillExtrapSteps,
2317 bool stopAtBoundary,
2318 double maxStep)
const
2321 static const unsigned int maxNumIt(500);
2322 unsigned int numIt(0);
2324 double coveredDistance(0.);
2327 const TVector3 W(destPlane.
getNormal());
2328 M1x4 SU = {{W.X(), W.Y(), W.Z(), destPlane.
distance(0., 0., 0.)}};
2331 if (W*destPlane.
getO() < 0) {
2338 M1x7 startState7 = state7;
2343 debugOut <<
"\n============ RKTrackRep::Extrap loop nr. " << numIt <<
" ============\n";
2345 debugOut <<
"fillExtrapSteps " << fillExtrapSteps <<
"\n";
2348 if(++numIt > maxNumIt){
2349 Exception exc(
"RKTrackRep::Extrap ==> maximum number of iterations exceeded",__LINE__,__FILE__);
2355 for(
int i = 0; i < 7*7; ++i)
J_MMT_[i] = 0;
2356 for(
int i=0; i<7; ++i)
J_MMT_[8*i] = 1.;
2358 M7x7* noise =
nullptr;
2359 isAtBoundary =
false;
2362 bool checkJacProj =
false;
2366 M1x7 J_MMT_unprojected_lastRow = {{0, 0, 0, 0, 0, 0, 1}};
2368 if( !
RKutta(SU, destPlane, charge, mass, state7, &
J_MMT_, &J_MMT_unprojected_lastRow,
2370 limits_, onlyOneStep, !fillExtrapSteps) ) {
2371 Exception exc(
"RKTrackRep::Extrap ==> Runge Kutta propagation failed",__LINE__,__FILE__);
2378 isAtBoundary =
true;
2384 debugOut <<
"stepSize = " <<
it->matStep_.stepSize_ <<
"\t";
2385 it->matStep_.material_.Print();
2393 if(fillExtrapSteps) {
2404 fabs(charge/state7[6]),
2411 debugOut <<
"momLoss: " << momLoss <<
" GeV; relative: " << momLoss/fabs(charge/state7[6])
2412 <<
"; coveredDistance = " << coveredDistance <<
"\n";
2413 if (
debugLvl_ > 1 && noise !=
nullptr) {
2420 if(fabs(state7[6])>1.E-10) {
2423 debugOut <<
"correct state7 with dx/dqop, dy/dqop ...\n";
2426 dqop = charge/(fabs(charge/state7[6])-momLoss) - state7[6];
2433 if (checkJacProj && fabs(coveredDistance) >
MINSTEP) {
2434 M1x3 state7_correction_unprojected = {{0, 0, 0}};
2435 for (
unsigned int i=0; i<3; ++i) {
2436 state7_correction_unprojected[i] = 0.5 * dqop * J_MMT_unprojected_lastRow[i];
2441 M1x3 state7_correction_projected = {{0, 0, 0}};
2442 for (
unsigned int i=0; i<3; ++i) {
2443 state7_correction_projected[i] = 0.5 * dqop *
J_MMT_[6*7 + i];
2449 M1x3 delta_state = {{0, 0, 0}};
2450 for (
unsigned int i=0; i<3; ++i) {
2451 delta_state[i] = state7_correction_unprojected[i] - state7_correction_projected[i];
2454 double Dist( sqrt(delta_state[0]*delta_state[0]
2455 + delta_state[1]*delta_state[1]
2456 + delta_state[2]*delta_state[2] ) );
2459 if (delta_state[0]*state7[3] + delta_state[1]*state7[4] + delta_state[2]*state7[5] > 0)
2462 double correctionFactor( 1. + Dist / coveredDistance );
2463 flightTime *= correctionFactor;
2464 momLoss *= correctionFactor;
2465 coveredDistance = coveredDistance + Dist;
2468 debugOut <<
"correctionFactor-1 = " << correctionFactor-1. <<
"; Dist = " << Dist <<
"\n";
2469 debugOut <<
"corrected momLoss: " << momLoss <<
" GeV; relative: " << momLoss/fabs(charge/state7[6])
2470 <<
"; corrected coveredDistance = " << coveredDistance <<
"\n";
2475 double qop( charge/(fabs(charge/state7[6])-momLoss) );
2476 dqop = qop - state7[6];
2479 for (
unsigned int i=0; i<6; ++i) {
2480 state7[i] += 0.5 * dqop *
J_MMT_[6*7 + i];
2483 double norm( 1. / sqrt(state7[3]*state7[3] + state7[4]*state7[4] + state7[5]*state7[5]) );
2484 for (
unsigned int i=3; i<6; ++i)
2492 if (fillExtrapSteps) {
2495 std::vector<ExtrapStep>::iterator lastStep =
ExtrapSteps_.end() - 1;
2498 lastStep->jac7_ =
J_MMT_;
2500 if( checkJacProj ==
true ){
2505 debugOut <<
"7D noise projected onto plane: \n";
2526 if (stopAtBoundary and isAtBoundary) {
2528 debugOut <<
"stopAtBoundary -> break; \n ";
2535 debugOut <<
"onlyOneStep -> break; \n ";
2543 debugOut <<
"arrived at destPlane with a distance of " << destPlane.
distance(state7[0], state7[1], state7[2]) <<
" cm left. ";
2544 if (destPlane.
isInActive(state7[0], state7[1], state7[2], state7[3], state7[4], state7[5]))
2545 debugOut <<
"In active area of destPlane. \n";
2547 debugOut <<
"NOT in active area of plane. \n";
2549 debugOut <<
" position: "; TVector3(state7[0], state7[1], state7[2]).Print();
2550 debugOut <<
" direction: "; TVector3(state7[3], state7[4], state7[5]).Print();
2557 if (fillExtrapSteps) {
2561 if (cov !=
nullptr) {
2567 if (cov !=
nullptr) {
2568 debugOut <<
"final covariance matrix after Extrap: "; cov->Print();
2573 return coveredDistance;
2579 if (state.
getRep() !=
this){
2580 Exception exc(
"RKTrackRep::checkCache ==> state is defined wrt. another TrackRep",__LINE__,__FILE__);
2585 if (dynamic_cast<const MeasurementOnPlane*>(&state) !=
nullptr) {
2586 Exception exc(
"RKTrackRep::checkCache - cannot extrapolate MeasurementOnPlane",__LINE__,__FILE__);
2607 double firstStep(0);
2608 for (
unsigned int i=0; i<
RKSteps_.size(); ++i) {
2610 firstStep =
RKSteps_.at(0).matStep_.stepSize_;
2613 if (
RKSteps_.at(i).matStep_.stepSize_ * firstStep < 0) {
2614 if (
RKSteps_.at(i-1).matStep_.material_ ==
RKSteps_.at(i).matStep_.material_) {
2622 debugOut <<
"RKTrackRep::checkCache: use cached material and step values.\n";
2628 debugOut <<
"RKTrackRep::checkCache: can NOT use cached material and step values.\n";
2630 if (plane !=
nullptr) {
2632 debugOut <<
"state.getPlane() != lastStartState_.getPlane()\n";
2636 debugOut <<
"state.getState() != lastStartState_.getState()\n";
2655 return fabs(1/state7[6]);
2660 if (dynamic_cast<const RKTrackRep*>(other) ==
nullptr)
2675 void RKTrackRep::Streamer(TBuffer &R__b)
2680 typedef ::genfit::RKTrackRep thisClass;
2682 if (R__b.IsReading()) {
2683 ::genfit::AbsTrackRep::Streamer(R__b);
2684 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
if (R__v) { }
2685 R__b.CheckByteCount(R__s, R__c, thisClass::IsA());
2689 ::genfit::AbsTrackRep::Streamer(R__b);
2690 R__c = R__b.WriteVersion(thisClass::IsA(), kTRUE);
2691 R__b.SetByteCount(R__c, kTRUE);