// $Id: DTrackWireBased_factory.cc 5612 2009-10-15 20:51:25Z staylor $ // // File: DTrackWireBased_factory.cc // Created: Wed Sep 3 09:33:40 EDT 2008 // Creator: davidl (on Darwin harriet.jlab.org 8.11.1 i386) // #include #include #include using namespace std; #include "DTrackWireBased_factory.h" #include #include #include #include using namespace jana; //------------------ // CDCSortByRincreasing //------------------ bool CDCSortByRincreasing(const DCDCTrackHit* const &hit1, const DCDCTrackHit* const &hit2) { // use the ring number to sort by R(decreasing) and then straw(increasing) if(hit1->wire->ring == hit2->wire->ring){ return hit1->wire->straw < hit2->wire->straw; } return hit1->wire->ring < hit2->wire->ring; } //------------------ // FDCSortByZincreasing //------------------ bool FDCSortByZincreasing(const DFDCPseudo* const &hit1, const DFDCPseudo* const &hit2) { // use the layer number to sort by Z(decreasing) and then wire(increasing) if(hit1->wire->layer == hit2->wire->layer){ return hit1->wire->wire < hit2->wire->wire; } return hit1->wire->layer < hit2->wire->layer; } //------------------ // count_common_members //------------------ template static unsigned int count_common_members(vector &a, vector &b) { unsigned int n=0; for(unsigned int i=0; iSetDefaultParameter("TRKFIT:DEBUG_LEVEL",DEBUG_LEVEL); gPARMS->SetDefaultParameter("TRKFIT:MOMENTUM_CUT_FOR_DEDX",MOMENTUM_CUT_FOR_DEDX); return NOERROR; } //------------------ // brun //------------------ jerror_t DTrackWireBased_factory::brun(jana::JEventLoop *loop, int runnumber) { // Get pointer to DTrackFitter object that actually fits a track vector fitters; loop->Get(fitters); if(fitters.size()<1){ _DBG_<<"Unable to get a DTrackFitter object! NO Charged track fitting will be done!"<(fitters[0]); // Warn user if something happened that caused us NOT to get a fitter object pointer if(!fitter){ _DBG_<<"Unable to get a DTrackFitter object! NO Charged track fitting will be done!"<SetDefaultParameter("TRKFIT:MASS_HYPOTHESES", MASS_HYPOTHESES); // Parse MASS_HYPOTHESES string to make list of masses to try if(MASS_HYPOTHESES.length()>0){ string &str = MASS_HYPOTHESES; unsigned int cutAt; while( (cutAt = str.find(",")) != (unsigned int)str.npos ){ if(cutAt > 0)mass_hypotheses.push_back(atof(str.substr(0,cutAt).c_str())); str = str.substr(cutAt+1); } if(str.length() > 0)mass_hypotheses.push_back(atof(str.c_str())); }else{ mass_hypotheses.push_back(0.0); // If empty string is specified, assume they want massless particle } return NOERROR; } //------------------ // evnt //------------------ jerror_t DTrackWireBased_factory::evnt(JEventLoop *loop, int eventnumber) { if(!fitter)return NOERROR; // Get candidates and hits vector candidates; loop->Get(candidates); // Deallocate some reference trajectories occasionally unsigned int rts_to_keep = 5; if(candidates.size()>rts_to_keep)rts_to_keep=candidates.size(); for(unsigned int i=rts_to_keep; iGetDMagneticFieldMap())); DReferenceTrajectory *rt = rtv[_data.size()]; // Loop over potential particle masses until one is found that gives a chisq/Ndof<3.0 // If none does, then use the one with the smallest chisq DTrackWireBased *best_track = NULL; double best_fom = 0.0; for(unsigned int j=0; j1){_DBG__;_DBG_<<"---- Starting wire based fit for candidate "<SetFitType(DTrackFitter::kWireBased); DTrackFitter::fit_status_t status = fitter->FindHitsAndFitTrack(*candidate, rt, loop, mass_hypotheses[j]); DTrackWireBased *dtrack = NULL; switch(status){ case DTrackFitter::kFitNotDone: _DBG_<<"Fitter returned kFitNotDone. This should never happen!!"<Ndof < 1){ if(DEBUG_LEVEL>1)_DBG_<<"-- new track with mass "<1)_DBG_<<"-- first successful fit this candidate with mass: "< best_fom){ if(DEBUG_LEVEL>1)_DBG_<<"-- new best track with mass "<1)_DBG_<<"-- keeping best track with mass "<mass()<<" (old chisq/Ndof="<<(best_track->chisq/best_track->Ndof)<<" , new chisq/Ndof="<<(dtrack->chisq/dtrack->Ndof)<<") (old fom="<2)_DBG_<<"adding wire-based track for candidate "<GetDMagneticFieldMap())); DReferenceTrajectory *rt = rtv[_data.size()]; DTrackWireBased *track = new DTrackWireBased; // Copy over DKinematicData part DKinematicData *track_kd = track; *track_kd = fitter->GetFitParameters(); rt->SetMass(track_kd->mass()); rt->Swim(track->position(), track->momentum(), track->charge()); track->rt = rt; track->chisq = fitter->GetChisq(); track->Ndof = fitter->GetNdof(); track->candidateid = candidate->id; // Add hits used as associated objects vector cdchits = fitter->GetCDCFitHits(); vector fdchits = fitter->GetFDCFitHits(); sort(cdchits.begin(), cdchits.end(), CDCSortByRincreasing); sort(fdchits.begin(), fdchits.end(), FDCSortByZincreasing); for(unsigned int i=0; iAddAssociatedObject(cdchits[i]); for(unsigned int i=0; iAddAssociatedObject(fdchits[i]); // Add DTrackCandidate as associated object track->AddAssociatedObject(candidate); return track; } //------------------ // GetFOM //------------------ // Uses dEdx from the track to provide a measure of the figure of merit for a track for a given mass // hypothesis double DTrackWireBased_factory::GetFOM(DTrackWireBased *dtrack) { double dedx,mean_path_length,p_avg; unsigned int num_hits=0; if (fitter->GetdEdx(dtrack->rt,dedx,mean_path_length,p_avg,num_hits) ==NOERROR){ dtrack->setdEdx(dedx); double dedx_sigma=fitter->GetdEdxSigma(num_hits,mean_path_length); double dedx_most_probable=fitter->GetdEdx(p_avg,dtrack->rt->GetMass(),mean_path_length); //figure of merit return ( dedx_sigma/fabs(dedx/dedx_most_probable-1.) ); } // If we got here, GetdEdx failed for this track dtrack->setdEdx(0.); return 0.; } //------------------ // GetRangeOutFOM //------------------ double DTrackWireBased_factory::GetRangeOutFOM(DTrackWireBased *dtrack) { /// Calculate a figure of merit for the track ranging out within the /// detector. If the particle does range out (lose all of its energy) /// then the value returned would be essentially zero. If it does not, /// then the value will be greater than zero. /// /// The FOM is ratio of the total pathlength of the track to the /// pathlength to the outermost wire associated with the track. /// Therefore, FOM=1 corresponds to a track that goes just /// as far after it hit the last wire as before, before finally /// hitting the BCAL, FCAL, etc... // We want the pathlength to the last wire that is on this track // (as determined by the DTrackHitSelector??? class). Since the // tracks can curl back in toward the beamline, we have to check every // wire to see what the pathlength to it is // Need the reference trajectory to find pathlengths _DBG__; DReferenceTrajectory *rt = const_cast(dtrack->rt); // Look first at FDC hits _DBG__; vector fdchits; dtrack->Get(fdchits); const DCoordinateSystem *outermost_wire = NULL; double s_to_outermost_wire=-1.0; _DBG__; for(unsigned int i=0; iFindClosestSwimStep(fdchits[i]->wire); if(step->s > s_to_outermost_wire){ s_to_outermost_wire = step->s; outermost_wire = fdchits[i]->wire; } } // Check CDC if no FDC wire was found _DBG__; if(!outermost_wire){ vector cdchits; dtrack->Get(cdchits); for(unsigned int i=0; iFindClosestSwimStep(cdchits[i]->wire); if(step->s > s_to_outermost_wire){ s_to_outermost_wire = step->s; outermost_wire = cdchits[i]->wire; } } } // Make sure *a* wire was found. (This is just a dummy check) _DBG__; if(!outermost_wire){ _DBG_<<"ERROR: No outermost wire found for track!"<rt->swim_steps[dtrack->rt->Nswim_steps-1].s; // Calculate figure of merit double fom = (total_s - s_to_outermost_wire)/s_to_outermost_wire; _DBG__; return fom; } //------------------ // FilterDuplicates //------------------ void DTrackWireBased_factory::FilterDuplicates(void) { /// Look through all current DTrackWireBased objects and remove any /// that have all of their hits in common with another track if(_data.size()==0)return; if(DEBUG_LEVEL>2)_DBG_<<"Looking for clones of wire-based tracks ..."< indexes_to_delete; for(unsigned int i=0; i<_data.size()-1; i++){ DTrackWireBased *dtrack1 = _data[i]; vector cdchits1; vector fdchits1; dtrack1->Get(cdchits1); dtrack1->Get(fdchits1); for(unsigned int j=i+1; j<_data.size(); j++){ DTrackWireBased *dtrack2 = _data[j]; vector cdchits2; vector fdchits2; dtrack2->Get(cdchits2); dtrack2->Get(fdchits2); // Count number of cdc and fdc hits in common unsigned int Ncdc = count_common_members(cdchits1, cdchits2); unsigned int Nfdc = count_common_members(fdchits1, fdchits2); if(Ncdc!=cdchits1.size() && Ncdc!=cdchits2.size())continue; if(Nfdc!=fdchits1.size() && Nfdc!=fdchits2.size())continue; unsigned int total = Ncdc + Nfdc; unsigned int total1 = cdchits1.size()+fdchits1.size(); unsigned int total2 = cdchits2.size()+fdchits2.size(); if(total!=total1 && total!=total2)continue; if(total12)_DBG_<<"Found "< new_data; for(unsigned int i=0; i<_data.size(); i++){ if(indexes_to_delete.find(i)==indexes_to_delete.end()){ new_data.push_back(_data[i]); }else{ delete _data[i]; if(DEBUG_LEVEL>1)_DBG_<<"Deleting clone wire-based track "<