Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DTrackCandidate_factory_THROWN.cc
Go to the documentation of this file.
1 // $Id$
2 //
3 // File: DTrackCandidate_factory_THROWN.cc
4 // Created: Tue Dec 12 12:42:57 EST 2006
5 // Creator: davidl (on Darwin swire-b241.jlab.org 8.8.0 powerpc)
6 //
7 
8 #include <cmath>
9 using namespace std;
10 
11 #include <JANA/JEventLoop.h>
12 #include <DANA/DApplication.h>
13 
14 #include <CDC/DCDCTrackHit.h>
15 #include <FDC/DFDCPseudo.h>
16 
18 #include <TRACKING/DMCThrown.h>
20 #include <TRACKING/DTrackFitter.h>
22 #include <DRandom.h>
23 #include <DMatrix.h>
24 
25 //------------------
26 // Constructor
27 //------------------
29 {
30  fitter = NULL;
31  hitselector = NULL;
32  bfield = NULL;
33 }
34 
35 //------------------
36 // brun
37 //------------------
38 jerror_t DTrackCandidate_factory_THROWN::brun(jana::JEventLoop *loop, int32_t runnumber)
39 {
40  DApplication* dapp=dynamic_cast<DApplication*>(eventLoop->GetJApplication());
41  bfield = dapp->GetBfield(runnumber);
42 
43  // Get pointer to DTrackFitter object that actually fits a track
44  vector<const DTrackFitter *> fitters;
45  loop->Get(fitters);
46  if(fitters.size()<1){
47  _DBG_<<"Unable to get a DTrackFitter object! NO Charged track fitting will be done!"<<endl;
48  return RESOURCE_UNAVAILABLE;
49  }
50 
51  // Drop the const qualifier from the DTrackFitter pointer (I'm surely going to hell for this!)
52  fitter = const_cast<DTrackFitter*>(fitters[0]);
53 
54  // Warn user if something happened that caused us NOT to get a fitter object pointer
55  if(!fitter){
56  _DBG_<<"ERROR: Unable to get a DTrackFitter object! Chisq for DTrackCandidate:THROWN will NOT be calculated!"<<endl;
57  return RESOURCE_UNAVAILABLE;
58  }
59 
60  // Get pointer to DTrackHitSelector object
61  vector<const DTrackHitSelector *> hitselectors;
62  loop->Get(hitselectors);
63  if(hitselectors.size()<1){
64  _DBG_<<"ERROR: Unable to get a DTrackHitSelector object! NO DTrackCandidate:THROWN objects will be created!"<<endl;
65  return RESOURCE_UNAVAILABLE;
66  }
67  hitselector = hitselectors[0];
68 
69  // Get the particle ID algorithms
70  loop->GetSingle(dParticleID);
71 
72  return NOERROR;
73 }
74 
75 //------------------
76 // evnt
77 //------------------
78 jerror_t DTrackCandidate_factory_THROWN::evnt(JEventLoop *loop, uint64_t eventnumber)
79 {
80  vector<const DMCThrown*> mcthrowns;
81  vector<const DCDCTrackHit*> cdctrackhits;
82  vector<const DFDCPseudo*> fdcpseudos;
83  loop->Get(mcthrowns);
84  loop->Get(cdctrackhits);
85  loop->Get(fdcpseudos);
86 
87  for(unsigned int i=0; i< mcthrowns.size(); i++){
88  const DMCThrown *thrown = mcthrowns[i];
89  const DKinematicData *kd_thrown = thrown;
90 
91  if(fabs(thrown->charge())<1)continue;
92 
93  // First, copy over the DKinematicData part
94  DTrackCandidate *candidate = new DTrackCandidate;
95  DKinematicData *kd_candidate = candidate;
96  *kd_candidate = *kd_thrown;
97 
98  // Add DMCThrown as associated object
99  candidate->AddAssociatedObject(thrown);
100 
101  // We need to swim a reference trajectory here. To avoid the overhead
102  // of allocating/deallocating them every event, we keep a pool and
103  // re-use them. If the pool is not big enough, then add one to the
104  // pool.
105  unsigned int locNumInitialReferenceTrajectories = rt_pool.size();
106  if(rt_pool.size()<=_data.size()){
107  // This is a little ugly, but only gets called a few times throughout the life of the process
108  // Note: these never get deleted, even at the end of process.
109  rt_pool.push_back(new DReferenceTrajectory(bfield));
110  }
111  DReferenceTrajectory *rt = rt_pool[_data.size()];
112  if(locNumInitialReferenceTrajectories == rt_pool.size()) //didn't create a new one
113  rt->Reset();
114  rt->q = candidate->charge();
115  candidate->rt = rt;
116  DVector3 pos = candidate->position();
117  DVector3 mom = candidate->momentum();
118  rt->Swim(pos, mom, candidate->charge());
119 
120  // Find hits that should be on this track and add them as associated objects
121  vector<const DCDCTrackHit*> cdchits;
122  vector<const DFDCPseudo*> fdchits;
123  if(hitselector)hitselector->GetCDCHits(DTrackHitSelector::kHelical, rt, cdctrackhits, cdchits);
124  if(hitselector)hitselector->GetFDCHits(DTrackHitSelector::kHelical, rt, fdcpseudos, fdchits);
125  for(unsigned int i=0; i<cdchits.size(); i++)candidate->AddAssociatedObject(cdchits[i]);
126  for(unsigned int i=0; i<fdchits.size(); i++)candidate->AddAssociatedObject(fdchits[i]);
127 
128  // We want to get chisq and Ndof values for this track using the hits from above.
129  // We do this using the DTrackFitter object. This more or less guarantees that the
130  // chisq calculation is done in the same way as it is for track fitting. Note
131  // that no fitting is actually done here so this should be reasonably fast
132  if(fitter){
133  fitter->Reset();
134  fitter->AddHits(cdchits);
135  fitter->AddHits(fdchits);
136  double chisq;
137  int Ndof;
138  fitter->ChiSq(DTrackFitter::kTimeBased, rt, &chisq, &Ndof);
139  candidate->chisq = chisq;
140  candidate->Ndof = Ndof;
141  }else{
142  candidate->chisq = 0.0;
143  candidate->Ndof = 0;
144  }
145 
146  _data.push_back(candidate);
147  }
148 
149  // Set CDC ring & FDC plane hit patterns
150  for(size_t loc_i = 0; loc_i < _data.size(); ++loc_i)
151  {
152  vector<const DCDCTrackHit*> locCDCTrackHits;
153  _data[loc_i]->Get(locCDCTrackHits);
154 
155  vector<const DFDCPseudo*> locFDCPseudos;
156  _data[loc_i]->Get(locFDCPseudos);
157 
158  _data[loc_i]->dCDCRings = dParticleID->Get_CDCRingBitPattern(locCDCTrackHits);
159  _data[loc_i]->dFDCPlanes = dParticleID->Get_FDCPlaneBitPattern(locFDCPseudos);
160  }
161 
162  return NOERROR;
163 }
164 
DApplication * dapp
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
TVector3 DVector3
Definition: DVector3.h:14
void AddHits(vector< const DCDCTrackHit * > cdchits)
DMagneticFieldMap * GetBfield(unsigned int run_number=1)
const DVector3 & position(void) const
void Swim(const DVector3 &pos, const DVector3 &mom, double q=-1000.0, const TMatrixFSym *cov=NULL, double smax=2000.0, const DCoordinateSystem *wire=NULL)
virtual double ChiSq(fit_type_t fit_type, DReferenceTrajectory *rt, double *chisq_ptr=NULL, int *dof_ptr=NULL, vector< pull_t > *pulls_ptr=NULL)=0
void Reset(void)
Definition: DTrackFitter.cc:94
#define _DBG_
Definition: HDEVIO.h:12
const DReferenceTrajectory * rt
pointer to reference trjectory representing this track (if any)
double charge(void) const
&lt;A href=&quot;index.html#legend&quot;&gt; &lt;IMG src=&quot;CORE.png&quot; width=&quot;100&quot;&gt; &lt;/A&gt;
float chisq
Chi-squared for the track (not chisq/dof!)
const DVector3 & momentum(void) const
jerror_t evnt(jana::JEventLoop *eventLoop, uint64_t eventnumber)
Called every event.
const DTrackFitter * fitter
jerror_t brun(jana::JEventLoop *eventLoop, int32_t runnumber)
Called everytime a new run number is detected.
int Ndof
Number of degrees of freedom in the fit.