59 #include "GblFitter.h"
60 #include "../include/GblFitStatus.h"
61 #include "GblFitterInfo.h"
64 #include "ICalibrationParametersDerivatives.h"
72 #include <FieldManager.h>
75 #include <Math/SMatrix.h>
77 #include <TVectorDfwd.h>
90 GblFitter::~GblFitter() {
91 if (m_segmentController) {
92 delete m_segmentController;
93 m_segmentController =
nullptr;
99 if (m_segmentController) {
100 delete m_segmentController;
101 m_segmentController =
nullptr;
103 m_segmentController = controler;
108 cleanGblInfo(trk, rep);
115 bool fitQoverP =
true;
118 if (!(Bfield > 1.e-16))
125 double lostWeight = 0.;
130 trk->setFitStatus(gblfs, rep);
131 gblfs->setCurvature(fitQoverP);
139 constructGblInfo(trk, rep)
144 gblfs->setIsFittedWithReferenceTrack(
true);
145 gblfs->setNumIterations(0);
146 if (m_externalIterations < 1)
151 unsigned int nFailed = 0;
154 std::vector<std::string> gblIterations;
155 gblIterations.push_back(m_gblInternalIterations);
159 for (
unsigned int iIter = 0; iIter < m_externalIterations; iIter++) {
161 int nscat = 0, nmeas = 0, ndummy = 0;
162 std::vector<gbl::GblPoint> points = collectGblPoints(trk, rep);
163 for(
unsigned int ip = 0;ip<points.size(); ip++) {
165 if (p.hasScatterer())
167 if (p.hasMeasurement())
169 if(!p.hasMeasurement()&&!p.hasScatterer())
174 fitRes = traj.
fit(Chi2, Ndf, lostWeight, (iIter == m_externalIterations - 1) ? m_gblInternalIterations :
"");
177 updateGblInfo(traj, trk, rep);
181 if (m_recalcJacobians > iIter) {
184 for (
unsigned int ip = 0; ip < trk->getNumPoints(); ip++) {
185 if (trk->getPoint(ip)->hasFitterInfo(rep) && (currFitterInfo =
dynamic_cast<GblFitterInfo*
>(trk->getPoint(ip)->
getFitterInfo(rep)))) {
188 prevFitterInfo = currFitterInfo;
193 gblfs->setIsFitted(
true);
194 gblfs->setIsFitConvergedPartially(fitRes == 0);
195 nFailed = trk->getNumPointsWithMeasurement() - nmeas;
196 gblfs->setNFailedPoints(nFailed);
197 gblfs->setIsFitConvergedFully(fitRes == 0 && nFailed == 0);
198 gblfs->setNumIterations(iIter + 1);
199 gblfs->setChi2(Chi2);
204 int npoints_meas = trk->getNumPointsWithMeasurement();
205 int npoints_all = trk->getNumPoints();
207 cout <<
"-------------------------------------------------------" << endl;
208 cout <<
" GBL processed genfit::Track " << endl;
209 cout <<
"-------------------------------------------------------" << endl;
210 cout <<
" # Track Points : " << npoints_all << endl;
211 cout <<
" # Meas. Points : " << npoints_meas << endl;
214 cout <<
" (" << ndummy <<
" dummy) ";
216 cout <<
" # GBL points meas : " << nmeas << endl;
217 cout <<
" # GBL points scat : " << nscat << endl;
218 cout <<
"-------------- GBL Fit Results ----------- Iteration " << iIter+1 <<
" " << ((iIter < gblIterations.size()) ? gblIterations[iIter] :
"") << endl;
219 cout <<
" Fit q/p parameter : " << (gblfs->hasCurvature() ? (
"True") : (
"False")) << endl;
220 cout <<
" Valid trajectory : " << ((traj.
isValid()) ? (
"True") : (
"False")) << endl;
221 cout <<
" Fit result : " << fitRes <<
" (0 for success)" << endl;
222 cout <<
" GBL track NDF : " << Ndf <<
" (-1 for failure)" << endl;
223 cout <<
" GBL track Chi2 : " << Chi2 << endl;
224 cout <<
" GBL track P-value : " << TMath::Prob(Chi2, Ndf) << endl;
225 cout <<
"-------------------------------------------------------" << endl;
235 for (
int ip = trk->getNumPoints() - 1; ip >=0; ip--) {
236 trk->getPoint(ip)->setScatterer(
nullptr);
237 trk->getPoint(ip)->deleteFitterInfo(rep);
239 if (!trk->getPoint(ip)->hasRawMeasurements())
240 trk->deletePoint(ip);
246 int npoints_meas = trk->getNumPointsWithMeasurement();
249 rep->
setTime(reference, trk->getTimeSeed());
250 rep->
setPosMom(reference, trk->getStateSeed());
253 reference.extrapolateToPlane(firstPlane);
255 double arcLenPos = 0;
258 for (
int ipoint_meas = 0; ipoint_meas < npoints_meas - 1; ipoint_meas++) {
260 TrackPoint* point_meas = trk->getPointWithMeasurement(ipoint_meas);
266 point_meas->setSortingParameter(arcLenPos);
267 arcLenPos += reference.extrapolateToPlane(nextPlane);
270 trk->getPointWithMeasurement(npoints_meas - 1)->setSortingParameter(arcLenPos);
276 std::vector<gbl::GblPoint> thePoints;
280 for (
unsigned int ip = 0; ip < trk->getNumPoints(); ip++) {
296 for (
unsigned int ip = 0; ip < trk->getNumPoints(); ip++) {
315 void GblFitter::getScattererFromMatList(
double& length,
316 double& theta,
double& s,
double& ds,
317 const double p,
const double mass,
const double charge,
318 const std::vector<genfit::MatStep>& steps)
const {
319 theta = 0.; s = 0.; ds = 0.; length = 0;
320 if (steps.empty())
return;
335 for (
unsigned int i = 0; i < steps.size(); i++) {
336 const MatStep step = steps.at(i);
338 double rho = 1. / step.material_.radiationLength;
339 len += fabs(step.stepSize_);
341 xmax = xmin + fabs(step.stepSize_);
345 sumxx += rho * (xmax - xmin);
347 sumx2x2 += rho * (xmax * xmax - xmin * xmin) / 2.;
349 sumx3x3 += rho * (xmax * xmax * xmax - xmin * xmin * xmin) / 3.;
352 if (sumxx < 1.0e-10)
return;
356 double beta = p /
sqrt(p * p + mass * mass);
357 theta = (0.0136 / p / beta) * fabs(charge) *
sqrt(sumxx) * (1. + 0.038 * log(sumxx));
363 double N = 1. / sumxx;
368 double ds_2 = N * (sumx3x3 - 2. * sumx2x2 * s + sumxx * s * s);
379 int npoints_meas = trk->getNumPointsWithMeasurement();
383 TMatrixD jacPointToPoint(dim, dim);
384 jacPointToPoint.UnitMatrix();
389 rep->
setTime(reference, trk->getTimeSeed());
390 rep->
setPosMom(reference, trk->getStateSeed());
393 reference.extrapolateToPlane(firstPlane);
395 double sumTrackLen = 0;
397 TMatrixDSym noise; TVectorD deltaState;
400 for (
int ipoint_meas = 0; ipoint_meas < npoints_meas; ipoint_meas++) {
402 TrackPoint* point_meas = trk->getPointWithMeasurement(ipoint_meas);
406 TVector3 trackDir = rep->
getDir(reference);
408 double trackMomMag = rep->
getMomMag(reference);
410 double particleCharge = rep->
getCharge(reference);
412 double particleMass = rep->
getMass(reference);
414 double trackLen = 0., scatTheta = 0., scatSMean = 0., scatDeltaS = 0.;
416 double theta1 = 0., theta2 = 0., s1 = 0., s2 = 0.;
418 TMatrixD jacMeas2Scat(dim, dim);
419 jacMeas2Scat.UnitMatrix();
423 if (ipoint_meas >= npoints_meas - 1) {
443 TVector3 segmentEntry = refCopy.getPos();
445 TVector3 segmentExit = refCopy.getPos();
447 getScattererFromMatList(trackLen,
457 s1 = 0.; s2 = scatSMean + scatDeltaS * scatDeltaS / (scatSMean - s1);
458 theta1 =
sqrt(scatTheta * scatTheta * scatDeltaS * scatDeltaS / (scatDeltaS * scatDeltaS + (scatSMean - s1) * (scatSMean - s1)));
459 theta2 =
sqrt(scatTheta * scatTheta * (scatSMean - s1) * (scatSMean - s1) / (scatDeltaS * scatDeltaS + (scatSMean - s1) * (scatSMean - s1)));
462 if (m_segmentController)
463 m_segmentController->controlTrackSegment(segmentEntry, segmentExit, scatTheta,
this);
466 if (m_enableScatterers && !m_enableIntermediateScatterer) {
469 }
else if (!m_enableScatterers) {
477 if (theta1 > scatEpsilon)
487 if (theta2 > scatEpsilon) {
497 scattp->setSortingParameter(point_meas->getSortingParameter() + s2);
502 for (
unsigned int itp = 0; itp < trk->getNumPoints(); itp++) {
503 if (trk->getPoint(itp) == point_meas) {
520 cout <<
" ERROR: The extrapolation to measurement point " << (ipoint_meas + 2) <<
" stepped back by " << nextStep <<
"cm !!! Track will be cut before this point." << endl;
531 cout <<
" ERROR: The extrapolation to measurement point " << (ipoint_meas + 2) <<
" stepped back by " << nextStep <<
"cm !!! Track will be cut before this point." << endl;
537 sumTrackLen += trackLen;
GblTrajectory definition.
unsigned int fit(double &Chi2, int &Ndf, double &lostWeight, std::string optionList="")
Perform fit of (valid) trajectory.
bool isValid() const
Retrieve validity of trajectory.
unsigned int getNumPoints() const
Retrieve number of point from trajectory.
virtual SharedPlanePtr constructPlane(const StateOnPlane &state) const =0
Construct (virtual) detector plane (use state's AbsTrackRep).
Abstract base class for a track representation.
virtual void setTime(StateOnPlane &state, double time) const =0
Set time at which the state was defined.
virtual unsigned int getDim() const =0
Get the dimension of the state vector used by the track representation.
double getMass(const StateOnPlane &state) const
Get tha particle mass in GeV/c^2.
virtual double getMomMag(const StateOnPlane &state) const =0
get the magnitude of the momentum in GeV.
virtual double getCharge(const StateOnPlane &state) const =0
Get the (fitted) charge of a state.
virtual void setPosMom(StateOnPlane &state, const TVector3 &pos, const TVector3 &mom) const =0
Set position and momentum of state.
virtual void getForwardJacobianAndNoise(TMatrixD &jacobian, TMatrixDSym &noise, TVectorD &deltaState) const =0
Get the jacobian and noise matrix of the last extrapolation.
virtual double extrapolateBy(StateOnPlane &state, double step, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state by step (cm) and returns the extrapolation length and, via reference,...
virtual double extrapolateToPlane(StateOnPlane &state, const genfit::SharedPlanePtr &plane, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state to plane, and returns the extrapolation length and, via reference,...
TVector3 getDir(const StateOnPlane &state) const
Get the direction vector of a state.
virtual std::vector< genfit::MatStep > getSteps() const =0
Get stepsizes and material properties of crossed materials of the last extrapolation.
TVector3 getFieldVal(const TVector3 &position)
This does NOT use the cache!
static FieldManager * getInstance()
Get singleton instance.
FitStatus for use with GblFitter.
Collects information needed and produced by a GblFitter/GBL and is specific to one AbsTrackRep of the...
void recalculateJacobian(GblFitterInfo *prevFitterInfo)
Re-extrapolates between prevFitterInfo and this point using forward state to update the Jacobian (if ...
void updateFitResults(gbl::GblTrajectory &traj)
Update fitter info from GBL fit results.
gbl::GblPoint constructGblPoint()
Collect all data and create a GblPoint.
void setJacobian(TMatrixD jacobian)
Set the Jacobian for further GblPoint construction.
TrackSegmentController for use with GblFitter.
A state with arbitrary dimension defined in a DetPlane.
Object containing AbsMeasurement and AbsFitterInfo objects.
AbsFitterInfo * getFitterInfo(const AbsTrackRep *rep=nullptr) const
Get fitterInfo for rep. Per default, use cardinal rep.
void setFitterInfo(genfit::AbsFitterInfo *fitterInfo)
Takes Ownership.
Collection of TrackPoint objects, AbsTrackRep objects and FitStatus objects.
bool sort()
Sort TrackPoint and according to their sorting parameters.
void insertPoint(TrackPoint *point, int id=-1)
Insert TrackPoint BEFORE TrackPoint with position id, if id >= 0.
const MeasuredStateOnPlane & getFittedState(int id=0, const AbsTrackRep *rep=nullptr, bool biased=true) const
Shortcut to get FittedStates.
double sqrt(double a)
sqrt for double
Namespace for the general broken lines package.
Defines for I/O streams used for error and debug printing.
std::shared_ptr< genfit::DetPlane > SharedPlanePtr
Shared Pointer to a DetPlane.
Simple struct containing MaterialProperties and stepsize in the material.