Skip to content
Snippets Groups Projects
Commit 9305eea2 authored by Stephen A. Wood's avatar Stephen A. Wood
Browse files

Copy DC tracks into the podd THaTrack list. Add H.tr.* to root file.

parent 33a2b4f4
No related branches found
No related tags found
No related merge requests found
......@@ -4,6 +4,7 @@ block H.dc.*
block H.hod.*
block H.cal.*
block H.aero.*
block H.tr.*
block g.evtyp
# TDC hits per paddle
......
......@@ -459,7 +459,7 @@ Int_t THcDC::ApplyCorrections( void )
}
//_____________________________________________________________________________
Int_t THcDC::CoarseTrack( TClonesArray& /* tracks */ )
Int_t THcDC::CoarseTrack( TClonesArray& tracks )
{
// Calculation of coordinates of particle track cross point with scint
// plane in the detector coordinate system. For this, parameters of track
......@@ -475,7 +475,26 @@ Int_t THcDC::CoarseTrack( TClonesArray& /* tracks */ )
}
// Now link the stubs between chambers
LinkStubs();
if(fNDCTracks > 0) TrackFit();
if(fNDCTracks > 0) {
TrackFit();
// Copy tracks into podd tracks list
for(Int_t itrack=0;itrack<fNDCTracks;itrack++) {
THaTrack* theTrack = NULL;
theTrack = AddTrack(tracks, 0.0, 0.0, 0.0, 0.0); // Leaving off trackID
// Should we add stubs with AddCluster?
THcDCTrack *tr = static_cast<THcDCTrack*>( fDCTracks->At(itrack));
theTrack->SetD(tr->GetX(), tr->GetY(), tr->GetXP(), tr->GetYP());
theTrack->SetFlag((UInt_t) 0);
Int_t nhits=tr->GetNHits();
// Need to look at how engine does chi2 and track selection. Reduced?
theTrack->SetChi2(tr->GetChisq(),nhits-4); // Nconstraints - Nparameters
// CalcFocalPlaneCoords. Aren't our tracks already in focal plane coords
// We should have some kind of track ID so that the THaTrack can be
// associate back with the DC track
}
}
// Check for internal TrackFit errors
// Histogram the focal plane tracks
// Histograms made in h_fill_dc_fp_hist
......
......@@ -40,6 +40,7 @@ public:
Double_t GetY() const {return fY_fp;}
Double_t GetXP() const {return fXp_fp;}
Double_t GetYP() const {return fYp_fp;}
Double_t GetChisq() const {return fChi2_fp;}
void SetNFree(Int_t nfree) {fNfree = nfree;}
void SetCoord(Int_t ip, Double_t coord) {fCoords[ip] = coord;}
void SetResidual(Int_t ip, Double_t coord) {fResiduals[ip] = coord;}
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment