5 #include "TMatrixDSym.h" 6 #include "TMatrixDSymEigen.h" 17 vector<float> segradlengths;
18 vector<float> cumseglens;
19 breakTrajInSegments(traj, breakpoints, segradlengths, cumseglens);
27 for (
unsigned int p = 0;
p<segradlengths.size();
p++) {
28 linearRegression(traj, breakpoints[
p], breakpoints[p+1], pcdir1);
30 if (segradlengths[p]<-100. || segradlengths[p-1]<-100.) {
31 dtheta.push_back(-999.);
33 const double cosval = pcdir0.X()*pcdir1.X()+pcdir0.Y()*pcdir1.Y()+pcdir0.Z()*pcdir1.Z();
36 double dt = 1000.*acos(cosval);
45 vector<float> cumLenFwd;
46 vector<float> cumLenBwd;
47 for (
unsigned int i = 0; i<cumseglens.size()-2; i++) {
48 cumLenFwd.push_back(cumseglens[i]);
49 cumLenBwd.push_back(cumseglens.back()-cumseglens[i+2]);
52 const ScanResult fwdResult = doLikelihoodScan(dtheta, segradlengths, cumLenFwd,
true , pid, detAngResol);
53 const ScanResult bwdResult = doLikelihoodScan(dtheta, segradlengths, cumLenBwd,
false, pid, detAngResol);
56 fwdResult.
p,fwdResult.
pUnc,fwdResult.
logL,
57 bwdResult.
p,bwdResult.
pUnc,bwdResult.
logL,
58 segradlengths,dtheta);
61 void TrajectoryMCSFitter::breakTrajInSegments(
const recob::TrackTrajectory& traj, vector<size_t>&
breakpoints, vector<float>& segradlengths, vector<float>& cumseglens)
const {
64 auto const* _SCE = (applySCEcorr_ ? lar::providerFrom<spacecharge::SpaceChargeService>() : NULL);
66 const double trajlen = traj.
Length();
67 const double thisSegLen = (trajlen>(segLen_*minNSegs_) ? segLen_ : trajlen/
double(minNSegs_) );
70 constexpr
double lar_radl_inv = 1./14.0;
71 cumseglens.push_back(0.);
75 breakpoints.push_back(nextValid);
78 const double Position[3] = {pos0.X(), pos0.Y(), pos0.Z()};
83 pos0_offset = _SCE->GetCalPosOffsets(p0, tpcid.
TPC);
85 pos0.SetX(pos0.X() - pos0_offset.X());
86 pos0.SetY(pos0.Y() + pos0_offset.Y());
87 pos0.SetZ(pos0.Z() + pos0_offset.Z());
96 const double Position[3] = {pos1.X(), pos1.Y(), pos1.Z()};
101 pos1_offset = _SCE->GetCalPosOffsets(p1, tpcid.
TPC);
103 pos1.SetX(pos1.X() - pos1_offset.X());
104 pos1.SetY(pos1.Y() + pos1_offset.Y());
105 pos1.SetZ(pos1.Z() + pos1_offset.Z());
108 auto step = (pos1-pos0).
R();
109 thislen += dir0.Dot(pos1-pos0);
123 if (thislen>=(thisSegLen-segLenTolerance_)) {
124 breakpoints.push_back(nextValid);
125 if (npoints>=minHitsPerSegment_) segradlengths.push_back(thislen*lar_radl_inv);
126 else segradlengths.push_back(-999.);
127 cumseglens.push_back(totlen);
136 segradlengths.push_back(thislen*lar_radl_inv);
137 cumseglens.push_back(cumseglens.back()+thislen);
142 const TrajectoryMCSFitter::ScanResult TrajectoryMCSFitter::doLikelihoodScan(std::vector<float>& dtheta, std::vector<float>& seg_nradlengths, std::vector<float>& cumLen,
143 bool fwdFit,
int pid,
float pmin,
float pmax,
float pstep,
float detAngResol)
const {
147 std::vector<float> vlogL;
148 for (
float p_test = pmin; p_test <= pmax; p_test+=pstep) {
149 float logL = mcsLikelihood(p_test, detAngResol, dtheta, seg_nradlengths, cumLen, fwdFit, pid);
150 if (logL < best_logL) {
153 best_idx = vlogL.size();
155 vlogL.push_back(logL);
161 for (
int j=best_idx-1;j>=0;j--) {
162 float dLL = vlogL[j]-vlogL[best_idx];
164 lunc = (best_idx-j)*pstep;
171 if (best_idx<
int(vlogL.size()-1)) {
172 for (
unsigned int j=best_idx+1;j<vlogL.size();j++) {
173 float dLL = vlogL[j]-vlogL[best_idx];
175 runc = (j-best_idx)*pstep;
183 const TrajectoryMCSFitter::ScanResult TrajectoryMCSFitter::doLikelihoodScan(std::vector<float>& dtheta, std::vector<float>& seg_nradlengths, std::vector<float>& cumLen,
184 bool fwdFit,
int pid,
float detAngResol)
const {
187 const ScanResult& coarseRes = doLikelihoodScan(dtheta, seg_nradlengths, cumLen, fwdFit, pid, pMin_, pMax_, pStepCoarse_, detAngResol);
189 float pmax =
std::min(coarseRes.
p+fineScanWindow_,pMax_);
190 float pmin =
std::max(coarseRes.
p-fineScanWindow_,pMin_);
197 const ScanResult& refineRes = doLikelihoodScan(dtheta, seg_nradlengths, cumLen, fwdFit, pid, pmin, pmax, pStep_, detAngResol);
205 auto const* _SCE = (applySCEcorr_ ? lar::providerFrom<spacecharge::SpaceChargeService>() : NULL);
209 size_t nextValid = firstPoint;
212 while (nextValid<lastPoint) {
215 const double Position[3] = {tempP.X(), tempP.Y(), tempP.Z()};
220 tempP_offset = _SCE->GetCalPosOffsets(ptemp, tpcid.
TPC);
222 tempP.SetX(tempP.X() - tempP_offset.X());
223 tempP.SetY(tempP.Y() + tempP_offset.Y());
224 tempP.SetZ(tempP.Z() + tempP_offset.Z());
226 middlePointCalc.
add(tempP);
232 const double norm = 1./double(npoints);
237 nextValid = firstPoint;
238 while (nextValid<lastPoint) {
241 const double Position[3] = {
p.X(),
p.Y(),
p.Z()};
246 p_offset = _SCE->GetCalPosOffsets(point, tpcid.
TPC);
248 p.SetX(
p.X() - p_offset.X());
249 p.SetY(
p.Y() + p_offset.Y());
250 p.SetZ(
p.Z() + p_offset.Z());
252 const double xxw0 =
p.X()-avgpos.X();
253 const double yyw0 =
p.Y()-avgpos.Y();
254 const double zzw0 =
p.Z()-avgpos.Z();
255 m(0, 0) += xxw0*xxw0*
norm;
256 m(0, 1) += xxw0*yyw0*
norm;
257 m(0, 2) += xxw0*zzw0*
norm;
258 m(1, 0) += yyw0*xxw0*
norm;
259 m(1, 1) += yyw0*yyw0*
norm;
260 m(1, 2) += yyw0*zzw0*
norm;
261 m(2, 0) += zzw0*xxw0*
norm;
262 m(2, 1) += zzw0*yyw0*
norm;
263 m(2, 2) += zzw0*zzw0*
norm;
267 const TMatrixDSymEigen
me(m);
268 const auto& eigenval = me.GetEigenValues();
269 const auto& eigenvec = me.GetEigenVectors();
272 double maxeval = eigenval(0);
273 for (
int i=1; i<3; ++i) {
274 if (eigenval(i)>maxeval) {
276 maxeval = eigenval(i);
280 pcdir =
Vector_t(eigenvec(0, maxevalidx),eigenvec(1, maxevalidx),eigenvec(2, maxevalidx));
285 double TrajectoryMCSFitter::mcsLikelihood(
double p,
double theta0x, std::vector<float>& dthetaij, std::vector<float>& seg_nradl, std::vector<float>& cumLen,
bool fwd,
int pid)
const {
287 const int beg = (fwd ? 0 : (dthetaij.size()-1));
288 const int end = (fwd ? dthetaij.size() : -1);
289 const int incr = (fwd ? +1 : -1);
293 const double m = mass(pid);
294 const double m2 = m*
m;
295 const double Etot = sqrt(p*p + m2);
298 double const fixedterm = 0.5 * std::log( 2.0 *
M_PI );
300 for (
int i = beg; i !=
end; i+=incr ) {
306 const double Eij = GetE(Etot,cumLen[i],m);
313 const double pij = sqrt(Eij2 - m2);
314 const double beta = sqrt( 1. - ((m2)/(pij*pij + m2)) );
315 constexpr
double HighlandSecondTerm = 0.038;
316 const double tH0 = ( HighlandFirstTerm(pij) / (pij*
beta) ) * ( 1.0 + HighlandSecondTerm * std::log( seg_nradl[i] ) ) * sqrt( seg_nradl[i] );
317 const double rms = sqrt( 2.0*( tH0 * tH0 + theta0x * theta0x ) );
319 std::cout <<
" Error : RMS cannot be zero ! " <<
std::endl;
322 const double arg = dthetaij[i]/
rms;
323 result += ( std::log( rms ) + 0.5 * arg * arg + fixedterm);
328 double TrajectoryMCSFitter::energyLossLandau(
const double mass2,
const double e2,
const double x)
const {
332 if (x<=0.)
return 0.;
333 constexpr
double Iinv2 = 1./(188.E-6*188.E-6);
334 constexpr
double matConst = 1.4*18./40.;
335 constexpr
double me = 0.511;
336 constexpr
double kappa = 0.307075;
337 constexpr
double j = 0.200;
339 const double beta2 = (e2-mass2)/e2;
340 const double gamma2 = 1./(1.0 - beta2);
341 const double epsilon = 0.5*kappa*x*matConst/beta2;
343 return 0.001*epsilon*( log(2.*me*beta2*gamma2*epsilon*Iinv2) + j - beta2 );
346 double TrajectoryMCSFitter::energyLossBetheBloch(
const double mass,
const double e2)
const {
348 constexpr
double Iinv = 1./188.E-6;
349 constexpr
double matConst = 1.4*18./40.;
350 constexpr
double me = 0.511;
351 constexpr
double kappa = 0.307075;
353 const double beta2 = (e2-mass*mass)/e2;
354 const double gamma2 = 1./(1.0 - beta2);
355 const double massRatio = me/mass;
356 const double argument = (2.*me*gamma2*beta2*Iinv) * std::sqrt(1+2*std::sqrt(gamma2)*massRatio + massRatio*massRatio);
358 double dedx = kappa*matConst/beta2;
360 if (mass==0.0)
return(0.0);
361 if (argument <= exp(beta2)) {
364 dedx *= (log(argument)-beta2)*1.
E-3;
365 if (dedx<0.) dedx = 0.;
370 double TrajectoryMCSFitter::GetE(
const double initial_E,
const double length_travelled,
const double m)
const {
374 constexpr
double kcal = 0.002105;
375 return (initial_E - kcal*length_travelled);
379 const double step_size = length_travelled / nElossSteps_;
381 double current_E = initial_E;
382 const double m2 = m*
m;
384 for (
auto i = 0; i < nElossSteps_; ++i) {
386 double dedx = energyLossBetheBloch(m,current_E);
387 current_E -= (dedx * step_size);
390 current_E -= energyLossLandau(m2,current_E*current_E,step_size);
392 if ( current_E <= m ) {
end
while True: pbar.update(maxval-len(onlies[E][S])) #print iS, "/", len(onlies[E][S]) found = False for...
Utilities to extend the interface of geometry vectors.
double beta(double KE, const simb::MCParticle *part)
double rms(sqlite3 *db, std::string const &table_name, std::string const &column_name)
T DirectionAtPoint(unsigned int p) const
Direction at point p. Use e.g. as:
recob::tracking::Vector_t Vector_t
Helper class to compute the middle point in a point set.
size_t LastValidPoint() const
Returns the index of the last valid point in the trajectory.
bool isValid
Whether this ID points to a valid element.
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Vector_t
Type for representation of momenta in 3D space.
void add(Point const &p)
Accumulates a point.
art framework interface to geometry description
geo::TPCID FindTPCAtPosition(double const worldLoc[3]) const
Returns the ID of the TPC at specified location.
double Length(size_t startAt=0) const
Returns the approximate length of the trajectory.
static constexpr double m2
A trajectory in space reconstructed from hits.
T LocationAtPoint(unsigned int p) const
Position at point p. Use e.g. as:
static int max(int a, int b)
The data type to uniquely identify a TPC.
auto norm(Vector const &v)
Return norm of the specified vector.
Class storing the result of the Maximum Likelihood fit of Multiple Coulomb Scattering angles between ...
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
geo::Point_t middlePoint() const
size_t FirstValidPoint() const
Returns the index of the first valid point in the trajectory.
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
static constexpr size_t InvalidIndex
Value returned on failed index queries.
size_t NextValidPoint(size_t index) const
Returns the index of the next valid point in the trajectory.
TPCID_t TPC
Index of the TPC within its cryostat.
Vector_t StartDirection() const
Returns the direction of the trajectory at the first point.
recob::tracking::Vector_t Vector_t
QTextStream & endl(QTextStream &s)