Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DTrackWireBased_factory_StraightLine.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: DTrackWireBased_factory_StraightLine.cc
4 // Created: Wed Mar 13 10:00:25 EDT 2019
5 // Creator: staylor (on Linux ifarm1402.jlab.org 3.10.0-327.el7.x86_64 x86_64)
6 //
7 
8 
9 #include <iostream>
10 #include <iomanip>
11 using namespace std;
12 
14 #include <CDC/DCDCTrackHit.h>
15 #include <FDC/DFDCPseudo.h>
16 using namespace jana;
17 
18 //------------------
19 // init
20 //------------------
22 {
23  CDC_MATCH_CUT=2.;
24  gPARMS->SetDefaultParameter("TRKFIT:CDC_MATCH_CUT",CDC_MATCH_CUT);
25  FDC_MATCH_CUT=1.25;
26  gPARMS->SetDefaultParameter("TRKFIT:FDC_MATCH_CUT",FDC_MATCH_CUT);
27 
28  return NOERROR;
29 }
30 
31 //------------------
32 // brun
33 //------------------
34 jerror_t DTrackWireBased_factory_StraightLine::brun(jana::JEventLoop *eventLoop, int32_t runnumber)
35 {
36  // Get the particle ID algorithms
37  eventLoop->GetSingle(dPIDAlgorithm);
38 
39  // Get pointer to TrackFinder object
40  vector<const DTrackFinder *> finders;
41  eventLoop->Get(finders);
42 
43  if(finders.size()<1){
44  _DBG_<<"Unable to get a DTrackFinder object!"<<endl;
45  return RESOURCE_UNAVAILABLE;
46  }
47 
48  // Drop the const qualifier from the DTrackFinder pointer
49  finder = const_cast<DTrackFinder*>(finders[0]);
50 
51  // Get pointer to DTrackFitter object that actually fits a track
52  vector<const DTrackFitter *> fitters;
53  eventLoop->Get(fitters,"StraightTrack");
54  if(fitters.size()<1){
55  _DBG_<<"Unable to get a DTrackFitter object!"<<endl;
56  return RESOURCE_UNAVAILABLE;
57  }
58 
59  // Drop the const qualifier from the DTrackFitter pointer
60  fitter = const_cast<DTrackFitter*>(fitters[0]);
61 
62  return NOERROR;
63 }
64 
65 //------------------
66 // evnt
67 //------------------
68 jerror_t DTrackWireBased_factory_StraightLine::evnt(JEventLoop *loop, uint64_t eventnumber)
69 {
70  // Get candidates
71  vector<const DTrackCandidate*> candidates;
72  loop->Get(candidates);
73 
74  // Get hits
75  vector<const DCDCTrackHit *>cdchits;
76  loop->Get(cdchits);
77  vector<const DFDCPseudo *>fdchits;
78  loop->Get(fdchits);
79 
80  for (unsigned int i=0;i<candidates.size();i++){
81  // Reset the fitter
82  fitter->Reset();
83 
84  const DTrackCandidate *cand=candidates[i];
85  DVector3 pos=cand->position();
86  DVector3 dir=cand->momentum();
87  // Select hits that belong to the track
88  for (unsigned int j=0;j<cdchits.size();j++){
89  double d=finder->FindDoca(pos,dir,cdchits[j]->wire->origin,
90  cdchits[j]->wire->udir);
91  if (d<CDC_MATCH_CUT) fitter->AddHit(cdchits[j]);
92  }
93  for (unsigned int i=0;i<fdchits.size();i++){
94  double pz=dir.z();
95  double tx=dir.x()/pz;
96  double ty=dir.y()/pz;
97  double dz=fdchits[i]->wire->origin.z()-pos.z();
98  DVector2 predpos(pos.x()+tx*dz,pos.y()+ty*dz);
99  DVector2 diff=predpos-fdchits[i]->xy;
100  if (diff.Mod()<FDC_MATCH_CUT) fitter->AddHit(fdchits[i]);
101  }
102 
103  // Fit the track using the list of hits we gathered above
104  if (fitter->FitTrack(pos,dir,1.,0.,0.)==DTrackFitter::kFitSuccess){
105  // Make a new wire-based track
107  *static_cast<DTrackingData*>(track) = fitter->GetFitParameters();
108 
109  track->chisq = fitter->GetChisq();
110  track->Ndof = fitter->GetNdof();
111  track->setPID(PiPlus);
112  track->FOM = TMath::Prob(track->chisq, track->Ndof);
113  track->pulls =std::move(fitter->GetPulls());
114  track->extrapolations=std::move(fitter->GetExtrapolations());
115  track->IsSmoothed = fitter->GetIsSmoothed();
116 
117  // candidate id
118  track->candidateid=i+1;
119 
120  // Add hits used as associated objects
121  vector<const DCDCTrackHit*> cdchits_on_track = fitter->GetCDCFitHits();
122  vector<const DFDCPseudo*> fdchits_on_track = fitter->GetFDCFitHits();
123 
124  // ...also find minimum drift time...
125  double tmin=1e6;
126  DetectorSystem_t detector=SYS_NULL;
127  for (unsigned int k=0;k<cdchits_on_track.size();k++){
128  if (cdchits_on_track[k]->tdrift<tmin){
129  tmin=cdchits_on_track[k]->tdrift;
130  detector=SYS_CDC;
131  }
132  track->AddAssociatedObject(cdchits_on_track[k]);
133  }
134  for (unsigned int k=0;k<fdchits_on_track.size();k++){
135  if (fdchits_on_track[k]->time<tmin){
136  tmin=fdchits_on_track[k]->time;
137  detector=SYS_FDC;
138  }
139  track->AddAssociatedObject(fdchits_on_track[k]);
140  }
141  track->setT0(tmin,10.0,detector);
142 
143  track->dCDCRings = dPIDAlgorithm->Get_CDCRingBitPattern(cdchits_on_track);
144  track->dFDCPlanes = dPIDAlgorithm->Get_FDCPlaneBitPattern(fdchits_on_track);
145 
146  _data.push_back(track);
147  }
148  }
149 
150  return NOERROR;
151 }
152 
153 //------------------
154 // erun
155 //------------------
157 {
158  return NOERROR;
159 }
160 
161 //------------------
162 // fini
163 //------------------
165 {
166  return NOERROR;
167 }
168 
map< DetectorSystem_t, vector< DTrackFitter::Extrapolation_t > > extrapolations
int GetNdof(void) const
Definition: DTrackFitter.h:155
jerror_t brun(jana::JEventLoop *eventLoop, int32_t runnumber)
Called everytime a new run number is detected.
Definition: track.h:16
vector< pull_t > & GetPulls(void)
Definition: DTrackFitter.h:162
float chisq
Chi-squared for the track (not chisq/dof!)
The DTrackFitter class is a base class for different charged track fitting algorithms. It does not actually fit the track itself, but provides the interface and some common support features most algorthims will need to implement.
Definition: DTrackFitter.h:61
TVector2 DVector2
Definition: DVector2.h:9
TVector3 DVector3
Definition: DVector3.h:14
Definition: GlueX.h:16
const DVector3 & position(void) const
bool GetIsSmoothed(void) const
Definition: DTrackFitter.h:160
void setT0(double at0, double at0_err, DetectorSystem_t at0_detector)
Definition: DTrackingData.h:64
DetectorSystem_t
Definition: GlueX.h:15
unsigned int dFDCPlanes
Definition: GlueX.h:17
jerror_t init(void)
Called once at program start.
void Reset(void)
Definition: DTrackFitter.cc:94
jerror_t evnt(jana::JEventLoop *eventLoop, uint64_t eventnumber)
Called every event.
vector< DTrackFitter::pull_t > pulls
Holds pulls used in chisq calc. (not including off-diagonals)
const map< DetectorSystem_t, vector< Extrapolation_t > > & GetExtrapolations(void) const
Definition: DTrackFitter.h:163
TLatex tx
double GetChisq(void) const
Definition: DTrackFitter.h:154
Definition: GlueX.h:18
jerror_t erun(void)
Called everytime run number changes, provided brun has been called.
#define _DBG_
Definition: HDEVIO.h:12
jerror_t fini(void)
Called after last event of last event source has been processed.
const vector< const DCDCTrackHit * > & GetCDCFitHits(void) const
Definition: DTrackFitter.h:139
&lt;A href=&quot;index.html#legend&quot;&gt; &lt;IMG src=&quot;CORE.png&quot; width=&quot;100&quot;&gt; &lt;/A&gt;
int Ndof
Number of degrees of freedom in the fit.
void setPID(Particle_t locPID)
const vector< const DFDCPseudo * > & GetFDCFitHits(void) const
Definition: DTrackFitter.h:140
void AddHit(const DCDCTrackHit *cdchit)
unsigned int dCDCRings
const DVector3 & momentum(void) const
fit_status_t FitTrack(const DVector3 &pos, const DVector3 &mom, double q, double mass, double t0=QuietNaN, DetectorSystem_t t0_det=SYS_NULL)
const DTrackFitter * fitter
oid_t candidateid
which DTrackCandidate this came from
TDirectory * dir
Definition: bcal_hist_eff.C:31
const DTrackingData & GetFitParameters(void) const
Definition: DTrackFitter.h:153