// $Id$ // // File: DTrackCandidate_factory.cc // Created: Mon Jul 18 15:23:04 EDT 2005 // Creator: davidl (on Darwin wire129.jlab.org 7.8.0 powerpc) // /// Factory to match track candidates generated by the FDC- and CDC-specific /// track finding codes. The code uses several algorithms in the following /// order of precedence: /// Method 1: Swims cdc candidates pointing in the forward direction /// and tries to match them geometrically to fdc candidates. /// Method 2: Projects the fdc candidates into the cdc region and tries /// to match to the cdc hits in each cdc candidate that points /// in the forward direction. /// Method 3: Similar to method 2, except that it projects a cdc /// candidate into the FDC region. /// Method 4: Matches left-over fdc candidates to other left-over fdc /// candidates. /// Method 5: Matches cdc candidates that do not point toward the cdc /// endplate that were not matched by previous methods to fdc /// candidates. /// Method 6: Tries to match stray unused cdc hits with fdc candidates /// Method 7: Alternate method to match leftover fdc candidates that /// have already been matched to other track candidates /// If a match is found, the code attempts to improve the track parameters by /// redoing the Riemann Circle fit with the additional hits. If an fdc /// candidate is matched to a cdc candidate, or several previously unused hits /// are added to the track, the code attempts to place the "start" position /// parameters of the track at a radius just outside the start counter. #include using namespace std; #include "DTrackCandidate_factory.h" #include "DTrackCandidate_factory_CDC.h" #include "DANA/DApplication.h" #include #include #include #include #define CUT 10. #define RADIUS_CUT 50.0 #define BEAM_VAR 0.0208 // cm^2 #define EPS 0.001 //------------------ // cdc_fdc_match //------------------ inline bool cdc_fdc_match(double p_fdc,double p_cdc,double dist){ //double frac=fabs(1.-p_cdc/p_fdc); //double frac2=fabs(1.-p_fdc/p_cdc); double p=p_fdc; if (p_cdc

DFDCPseudo[0]->DFDCWire->layer int layer1 = 100; // defaults just in case there is a segment with no hits int layer2 = 100; if(segment1->hits.size()>0)layer1=segment1->hits[0]->wire->layer; if(segment2->hits.size()>0)layer2=segment2->hits[0]->wire->layer; return layer1 < layer2; } //------------------ // CDCHitSortByLayerincreasing //------------------ inline bool CDCHitSortByLayerincreasing(const DCDCTrackHit* const &hit1, const DCDCTrackHit* const &hit2) { // Used to sort CDC hits by layer (ring) with innermost layer hits first // if same ring, sort by wire number if(hit1->wire->ring == hit2->wire->ring) { if(hit1->wire->straw == hit2->wire->straw) return (hit1->dE > hit2->dE); return hit1->wire->straw < hit2->wire->straw; } return hit1->wire->ring < hit2->wire->ring; } //------------------ // FDCHitSortByLayerincreasing //------------------ inline bool FDCHitSortByLayerincreasing(const DFDCPseudo* const &hit1, const DFDCPseudo* const &hit2) { // Used to sort FDC hits by layer with most upstream layer hits first // if same layer, sort by wire number if(hit1->wire->layer == hit2->wire->layer) { if(hit1->wire->wire == hit2->wire->wire) return (hit1->dE > hit2->dE); return (hit1->wire->wire < hit2->wire->wire); } return hit1->wire->layer < hit2->wire->layer; } //------------------ // init //------------------ jerror_t DTrackCandidate_factory::init(void) { MAX_NUM_TRACK_CANDIDATES = 20; return NOERROR; } //------------------ // brun //------------------ jerror_t DTrackCandidate_factory::brun(JEventLoop* eventLoop,int runnumber){ DApplication* dapp=dynamic_cast(eventLoop->GetJApplication()); const DGeometry *dgeom = dapp->GetDGeometry(runnumber); bfield = dapp->GetBfield(runnumber); FactorForSenseOfRotation=(bfield->GetBz(0.,0.,65.)>0.)?-1.:1.; dIsNoFieldFlag = (dynamic_cast(bfield) != NULL); if(dIsNoFieldFlag) { //Setting this flag makes it so that JANA does not delete the objects in _data. This factory will manage this memory. //This is all of these pointers are just copied from the "StraightLine" factory, and should not be re-deleted. SetFactoryFlag(NOT_OBJECT_OWNER); } else ClearFactoryFlag(NOT_OBJECT_OWNER); //This factory will create it's own objects. // Get the position of the exit of the CDC endplate from DGeometry double endplate_z,endplate_dz,endplate_rmin; dgeom->GetCDCEndplate(endplate_z,endplate_dz,endplate_rmin,endplate_rmax); cdc_endplate.SetZ(endplate_z+endplate_dz); dParticleID = NULL; eventLoop->GetSingle(dParticleID); JCalibration *jcalib = dapp->GetJCalibration(runnumber); map targetparms; if (jcalib->Get("TARGET/target_parms",targetparms)==false){ TARGET_Z = targetparms["TARGET_Z_POSITION"]; } else{ dgeom->GetTargetZ(TARGET_Z); } // Initialize the stepper stepper=new DMagneticFieldStepper(bfield); stepper->SetStepSize(1.0); DEBUG_HISTS=false; gPARMS->SetDefaultParameter("TRKFIND:DEBUG_HISTS",DEBUG_HISTS); if(DEBUG_HISTS){ dapp->Lock(); /* match_center_dist2=(TH2F*)gROOT->FindObject("match_center_dist2"); if (!match_center_dist2){ match_center_dist2=new TH2F("match_center_dist2","larger radius vs matching distance squared between two circle centers",100,0,100.,100,0,100); match_center_dist2->SetYTitle("r_{c} [cm]"); match_center_dist2->SetXTitle("(#Deltad)^{2} [cm^{2}]"); } */ match_dist=(TH2F*)gROOT->FindObject("match_dist"); if (!match_dist){ match_dist=new TH2F("match_dist","Matching distance", 120,0.,60.,500,0,25.); match_dist->SetXTitle("r (cm)"); match_dist->SetYTitle("#Deltar (cm)"); } match_dist_vs_p=(TH2F*)gROOT->FindObject("match_dist_vs_p"); if (!match_dist_vs_p) { match_dist_vs_p=new TH2F("match_dist_vs_p","Matching distance vs p", 50,0.,7.,100,0,25.); match_dist_vs_p->SetYTitle("#Deltar (cm)"); match_dist_vs_p->SetXTitle("p (GeV/c)"); } dapp->Unlock(); } gPARMS->SetDefaultParameter("TRKFIND:MAX_NUM_TRACK_CANDIDATES", MAX_NUM_TRACK_CANDIDATES); MIN_NUM_HITS=6; gPARMS->SetDefaultParameter("TRKFIND:MIN_NUM_HITS", MIN_NUM_HITS); DEBUG_LEVEL=0; gPARMS->SetDefaultParameter("TRKFIND:DEBUG_LEVEL", DEBUG_LEVEL); return NOERROR; } //------------------ // erun //------------------ jerror_t DTrackCandidate_factory::erun(void) { if (stepper) delete stepper; return NOERROR; } //------------------ // erun //------------------ jerror_t DTrackCandidate_factory::fini(void) { if (stepper) delete stepper; return NOERROR; } //------------------ // evnt //------------------ jerror_t DTrackCandidate_factory::evnt(JEventLoop *loop, int eventnumber) { if(dIsNoFieldFlag) { //Clear previous objects: //JANA doesn't do it because NOT_OBJECT_OWNER was set //It DID delete them though, in the "StraightLine" factory _data.clear(); vector locStraightLineCandidates; loop->Get(locStraightLineCandidates, "StraightLine"); for(size_t loc_i = 0; loc_i < locStraightLineCandidates.size(); ++loc_i) _data.push_back(const_cast(locStraightLineCandidates[loc_i])); return NOERROR; } // Clear private vectors cdctrackcandidates.clear(); fdctrackcandidates.clear(); trackcandidates.clear(); mycdchits.clear(); // Get the track candidates from the CDC and FDC candidate finders loop->Get(cdctrackcandidates, "CDC"); loop->Get(fdctrackcandidates, "FDCCathodes"); // List of cdc hits loop->Get(mycdchits); // Vectors for keeping track of cdc matches and projections to the endplate vector cdc_forward_ids; vector cdc_backward_ids; vector cdc_endplate_projections; // Vector to keep track of cdc hits used in candidates vectorused_cdc_hits(mycdchits.size()); // vector to keep track of the matching status of each fdc candidate vectorforward_matches(fdctrackcandidates.size()); // Normal vector for CDC endplate DVector3 norm(0,0,1); //Loop through the list of CDC candidates, flagging those that point toward // the FDC. The others are put in the final candidate list. for(unsigned int i=0; imomentum(); DVector3 pos=srccan->position(); double theta=mom.Theta(); // Propagate track to CDC endplate if (theta0){ cdc_forward_ids.push_back(i); //ProjectHelixToZ(cdc_endplate.z(),srccan->charge(),mom,pos); stepper->SetCharge(srccan->charge()); stepper->SwimToPlane(pos,mom,cdc_endplate,norm,NULL); cdc_endplate_projections.push_back(pos); } else{ cdc_backward_ids.push_back(i); } } // Variables for candidate number accounting int num_forward_cdc_cands_remaining=cdc_forward_ids.size(); int num_fdc_cands_remaining=fdctrackcandidates.size(); // Loop through the list of FDC candidates looking for matches between the // CDC and the FDC in the transition region. if (num_forward_cdc_cands_remaining>0){ for(unsigned int i=0; isegments; fdccan->GetT(segments); sort(segments.begin(), segments.end(), SegmentSortByLayerincreasing); if (segments[0]->package!=0) continue; // Flag for matching status bool got_match=false; if (DEBUG_LEVEL>0) _DBG_ << "Attempting FDC/CDC matching method #1..." <0) _DBG_ << "... matched to FDC candidate #" << i <0) _DBG_ << "Attempting FDC/CDC matching method #2..." <0) _DBG_ << "... matched to FDC candidate #" << i <cdc_forward_matches(cdc_forward_ids.size()); for (unsigned int i=0;i0) _DBG_ << "Attempting FDC/CDC matching method #3..." <0) _DBG_ << "... matched to CDC candidate #" << i <0){ for (unsigned int j=0;jcdchits; cdccan->GetT(cdchits); sort(cdchits.begin(), cdchits.end(), CDCHitSortByLayerincreasing); can->setMomentum(cdccan->momentum()); can->setPosition(cdccan->position()); can->setCharge(cdccan->charge()); for (unsigned int n=0;nused_cdc_indexes[n]]=1; can->AddAssociatedObject(cdchits[n]); } can->chisq=cdccan->chisq; can->Ndof=cdccan->Ndof; trackcandidates.push_back(can); } } for (unsigned int j=0;j0 && MatchMethod8(cdccan,forward_matches,used_cdc_hits)==true){ num_fdc_cands_remaining--; } else{ DTrackCandidate *can = new DTrackCandidate; can->setMomentum(cdccan->momentum()); can->setPosition(cdccan->position()); can->setCharge(cdccan->charge()); // Get the cdc hits and add them to the candidate vectorcdchits; cdccan->GetT(cdchits); sort(cdchits.begin(), cdchits.end(), CDCHitSortByLayerincreasing); for (unsigned int n=0;nused_cdc_indexes[n]]=1; can->AddAssociatedObject(cdchits[n]); } can->chisq=cdccan->chisq; can->Ndof=cdccan->Ndof; trackcandidates.push_back(can); } } unsigned int num_unmatched_cdcs=0; for (unsigned int i=0;i1){ for (unsigned int j=0;j0){ for (unsigned int j=0;jsegments; fdccan->GetT(segments); if (segments.size()>1){ // Create new candidate for output DTrackCandidate *can = new DTrackCandidate; can->setMomentum(fdccan->momentum()); can->setPosition(fdccan->position()); can->setCharge(fdccan->charge()); can->chisq=fdccan->chisq; can->Ndof=fdccan->Ndof; for (unsigned int m=0;mhits.size();n++){ const DFDCPseudo *fdchit=segments[m]->hits[n]; can->AddAssociatedObject(fdchit); } } trackcandidates.push_back(can); num_fdc_cands_remaining--; forward_matches[j]=1; } } } } // Try to match remaining fdc candidates to track candidates that have // track segments that have already been matched together if (num_fdc_cands_remaining>0){ for (unsigned int i=0;isegments; srccan->GetT(segments); // redo circle fits for segments, forcing the circles to go through (0,0) if (MatchMethod9(j,srccan,segments[0],fdctrackcandidates, forward_matches)){ if (DEBUG_LEVEL>0) _DBG_ << "Matched FDC segments using method #9" << endl; num_fdc_cands_remaining-=2; } // Redo circle fit assuming track is almost straight else if (MatchMethod10(j,srccan,segments[0],fdctrackcandidates, forward_matches)){ if (DEBUG_LEVEL>0) _DBG_ << "Matched FDC segments using method #10" << endl; num_fdc_cands_remaining-=2; } else { // Not much we can do here -- add to the final list of candidates DTrackCandidate *can = new DTrackCandidate; can->Ndof=srccan->Ndof; can->chisq=srccan->chisq; can->setCharge(srccan->charge()); can->setMomentum(srccan->momentum()); can->setPosition(srccan->position()); for (unsigned int n=0;nhits.size();n++){ const DFDCPseudo *fdchit=segments[0]->hits[n]; can->AddAssociatedObject(fdchit); } trackcandidates.push_back(can); } } } } // Only output the candidates that have at least a minimum number of hits for (unsigned int i=0;imyfdchits; candidate->GetT(myfdchits); vectormycdchits; candidate->GetT(mycdchits); if (int(mycdchits.size()+myfdchits.size())>=MIN_NUM_HITS){ _data.push_back(candidate); } else delete candidate; } if((int(_data.size()) > MAX_NUM_TRACK_CANDIDATES) && (MAX_NUM_TRACK_CANDIDATES >= 0)) { if (DEBUG_LEVEL>0) _DBG_ << "Number of candidates = " << _data.size() << " > " << MAX_NUM_TRACK_CANDIDATES << " --- skipping track fitting! " << endl; for(size_t loc_i = 0; loc_i < _data.size(); ++loc_i) delete _data[loc_i]; _data.clear(); } // Set CDC ring & FDC plane hit patterns for(size_t loc_i = 0; loc_i < _data.size(); ++loc_i) { vector locCDCTrackHits; _data[loc_i]->Get(locCDCTrackHits); vector locFDCPseudos; _data[loc_i]->Get(locFDCPseudos); _data[loc_i]->dCDCRings = dParticleID->Get_CDCRingBitPattern(locCDCTrackHits); _data[loc_i]->dFDCPlanes = dParticleID->Get_FDCPlaneBitPattern(locFDCPseudos); } return NOERROR; } // Obtain position and momentum at the exit of a given package using the // helical track model. jerror_t DTrackCandidate_factory::GetPositionAndMomentum( const DFDCSegment *segment, DVector3 &pos, DVector3 &mom){ // Position of track segment at last hit plane of package double x=segment->xc+segment->rc*cos(segment->Phi1); double y=segment->yc+segment->rc*sin(segment->Phi1); double z=segment->hits[0]->wire->origin.z(); // Track parameters //double kappa=segment->q/(2.*segment->rc); double phi0=segment->phi0; double tanl=segment->tanl; double z0=segment->z_vertex; // Useful intermediate variables double cosp=cos(phi0); double sinp=sin(phi0); double sperp=(z-z0)/tanl; //double twoks=2.*kappa*sperp; double twoks=FactorForSenseOfRotation*segment->q*sperp/segment->rc; double sin2ks=sin(twoks); double cos2ks=cos(twoks); // Get Bfield double B=fabs(bfield->GetBz(x,y,z)); // Momentum double pt=0.003*B*segment->rc; double px=pt*(cosp*cos2ks-sinp*sin2ks); double py=pt*(sinp*cos2ks+cosp*sin2ks); double pz=pt*tanl; pos.SetXYZ(x,y,z); mom.SetXYZ(px,py,pz); return NOERROR; } // Get the position and momentum at a fixed radius from the beam line jerror_t DTrackCandidate_factory::GetPositionAndMomentum(DHelicalFit &fit, double Bz, DVector3 &pos, DVector3 &mom){ double r2=90.0; double xc=fit.x0; double yc=fit.y0; double rc=fit.r0; double rc2=rc*rc; double xc2=xc*xc; double yc2=yc*yc; double xc2_plus_yc2=xc2+yc2; double a=(r2-xc2_plus_yc2-rc2)/(2.*rc); double temp1=yc*sqrt(xc2_plus_yc2-a*a); double temp2=xc*a; double cosphi_plus=(temp2+temp1)/xc2_plus_yc2; double cosphi_minus=(temp2-temp1)/xc2_plus_yc2; // Direction tangent and transverse momentum double tanl=fit.tanl; double pt=0.003*Bz*rc; if(!isfinite(temp1) || !isfinite(temp2)){ // We did not find an intersection between the two circles, so return // an error. The values of mom and pos are not changed. // _DBG_ << endl; return VALUE_OUT_OF_RANGE; } double phi_plus=acos(cosphi_plus); double phi_minus=acos(cosphi_minus); double x_plus=xc+rc*cosphi_plus; double x_minus=xc+rc*cosphi_minus; double y_plus=yc+rc*sin(phi_plus); double y_minus=yc+rc*sin(phi_minus); // if the resulting radial position on the circle from the fit does not agree // with the radius to which we are matching, we have the wrong sign for phi+ // or phi- double r2_plus=x_plus*x_plus+y_plus*y_plus; double r2_minus=x_minus*x_minus+y_minus*y_minus; if (fabs(r2-r2_plus)>EPS){ phi_plus*=-1.; y_plus=yc+rc*sin(phi_plus); } if (fabs(r2-r2_minus)>EPS){ phi_minus*=-1.; y_minus=yc+rc*sin(phi_minus); } if (fit.h<0.){ phi_minus=M_PI-phi_minus; pos.SetXYZ(x_minus,y_minus,0.); // z will be filled later mom.SetXYZ(pt*sin(phi_minus),pt*cos(phi_minus),pt*tanl); } else{ phi_plus*=-1.; //phi_plus=M_PI-phi_plus; pos.SetXYZ(x_plus,y_plus,0.); // z will be filled later mom.SetXYZ(pt*sin(phi_plus),pt*cos(phi_plus),pt*tanl); } return NOERROR; } // Get the position and momentum at a fixed radius from the beam line jerror_t DTrackCandidate_factory::GetPositionAndMomentum(DHelicalFit &fit, double Bz, const DVector3 &origin, DVector3 &pos, DVector3 &mom){ double r2=90.0; double xc=fit.x0; double yc=fit.y0; double rc=fit.r0; double rc2=rc*rc; double xc2=xc*xc; double yc2=yc*yc; double xc2_plus_yc2=xc2+yc2; double a=(r2-xc2_plus_yc2-rc2)/(2.*rc); double temp1=yc*sqrt(xc2_plus_yc2-a*a); double temp2=xc*a; double cosphi_plus=(temp2+temp1)/xc2_plus_yc2; double cosphi_minus=(temp2-temp1)/xc2_plus_yc2; // Direction tangent and transverse momentum double tanl=fit.tanl; double pt=0.003*Bz*rc; if(!isfinite(temp1) || !isfinite(temp2)){ // We did not find an intersection between the two circles, so return // an error. The values of mom and pos are not changed. // _DBG_ << endl; return VALUE_OUT_OF_RANGE; } double phi_plus=acos(cosphi_plus); double phi_minus=acos(cosphi_minus); double x_plus=xc+rc*cosphi_plus; double x_minus=xc+rc*cosphi_minus; double y_plus=yc+rc*sin(phi_plus); double y_minus=yc+rc*sin(phi_minus); // if the resulting radial position on the circle from the fit does not agree // with the radius to which we are matching, we have the wrong sign for phi+ // or phi- double r2_plus=x_plus*x_plus+y_plus*y_plus; double r2_minus=x_minus*x_minus+y_minus*y_minus; if (fabs(r2-r2_plus)>EPS){ phi_plus*=-1.; y_plus=yc+rc*sin(phi_plus); } if (fabs(r2-r2_minus)>EPS){ phi_minus*=-1.; y_minus=yc+rc*sin(phi_minus); } // Choose phi- or phi+ depending on proximity to one of the cdc hits double xwire=origin.x(); double ywire=origin.y(); double dx=x_minus-xwire; double dy=y_minus-ywire; double d2_minus=dx*dx+dy*dy; dx=x_plus-xwire; dy=y_plus-ywire; double d2_plus=dx*dx+dy*dy; if (d2_plus>d2_minus){ fit.h=-1.; phi_minus=M_PI-phi_minus; pos.SetXYZ(x_minus,y_minus,0.); // z will be filled later mom.SetXYZ(pt*sin(phi_minus),pt*cos(phi_minus),pt*tanl); } else{ fit.h=1.; phi_plus*=-1.; pos.SetXYZ(x_plus,y_plus,0.); // z will be fillted later mom.SetXYZ(pt*sin(phi_plus),pt*cos(phi_plus),pt*tanl); } return NOERROR; } // Find the position along a helical path at the z-position z void DTrackCandidate_factory::ProjectHelixToZ(const double z,const double q, const DVector3 &mom, DVector3 &pos){ double pt=mom.Perp(); double phi=mom.Phi(); double sinphi=sin(phi); double cosphi=cos(phi); double tanl=tan(M_PI_2-mom.Theta()); double x0=pos.X(); double y0=pos.Y(); double z0=pos.Z(); double B=fabs(bfield->GetBz(x0,y0,z0)); double twokappa=FactorForSenseOfRotation*0.003*B*q/pt; double one_over_twokappa=1./twokappa; double sperp=(z-z0)/tanl; double twoks=twokappa*sperp; double sin2ks=sin(twoks); double one_minus_cos2ks=1.-cos(twoks); double x=x0+(cosphi*sin2ks-sinphi*one_minus_cos2ks)*one_over_twokappa; double y=y0+(sinphi*sin2ks+cosphi*one_minus_cos2ks)*one_over_twokappa; pos.SetXYZ(x,y,z); } // Routine to do a crude match between cdc wires and a helical approximation to // the trajectory double DTrackCandidate_factory::DocaToHelix(const DCDCTrackHit *hit, double q, const DVector3 &pos, const DVector3 &mom){ double pt=mom.Perp(); double phi=mom.Phi(); double sinphi=sin(phi); double cosphi=cos(phi); double tanl=tan(M_PI_2-mom.Theta()); double x0=pos.X(); double y0=pos.Y(); double z0=pos.Z(); double B=fabs(bfield->GetBz(x0,y0,z0)); double twokappa=FactorForSenseOfRotation*0.003*B*q/pt; double one_over_twokappa=1./twokappa; double sperp=0; DVector3 origin=hit->wire->origin; double z0w=origin.z(); DVector3 dir=(1./hit->wire->udir.z())*hit->wire->udir; DVector3 wirepos; double old_doca2=1e8; double doca2=old_doca2; double z=z0; double x=x0,y=y0; while (z>50. && z<170. && x<60. && y<60.){ old_doca2=doca2; sperp-=1.; double twoks=twokappa*sperp; double sin2ks=sin(twoks); double one_minus_cos2ks=1.-cos(twoks); x=x0+(cosphi*sin2ks-sinphi*one_minus_cos2ks)*one_over_twokappa; y=y0+(sinphi*sin2ks+cosphi*one_minus_cos2ks)*one_over_twokappa; z=z0+sperp*tanl; wirepos=origin+(z-z0w)*dir; double dxw=x-wirepos.x(); double dyw=y-wirepos.y(); doca2=dxw*dxw+dyw*dyw; if (doca2>old_doca2){ break; } } return old_doca2; } // Find the particle charge given the helical fit result and an fdc hit // on the track double DTrackCandidate_factory::GetSenseOfRotation(DHelicalFit &fit, const DFDCPseudo *fdchit, const DVector3 &pos ){ // Get circle parameters double rc=fit.r0; double xc=fit.x0; double yc=fit.y0; // Compute phi rotation from "vertex" to fdc hit double dphi=(fdchit->wire->origin.Z()-pos.z())/(rc*fit.tanl); double Phi1=atan2(pos.Y()-yc,pos.X()-xc); // Positive and negative changes in phi double phiplus=Phi1+dphi; double phiminus=Phi1-dphi; DVector2 plus(xc+rc*cos(phiplus),yc+rc*sin(phiplus)); DVector2 minus(xc+rc*cos(phiminus),yc+rc*sin(phiminus)); // Compute differences double d2plus=(plus-fdchit->xy).Mod2(); double d2minus=(minus-fdchit->xy).Mod2(); if (d2minussegments, vectorcdchits, double &Bz){ unsigned int num_hits=0; // Initialize Bz Bz=0.; // Add the cdc axial wires to the list of hits to use in the fit for (unsigned int k=0;kis_stereo==false){ double cov=0.213; //guess const DVector3 origin=cdchits[k]->wire->origin; double x=origin.x(),y=origin.y(),z=origin.z(); fit.AddHitXYZ(x,y,z,cov,cov,0.,true); Bz+=bfield->GetBz(x,y,z); num_hits++; } } // Add the FDC hits and estimate Bz for (unsigned int k=0;khits.size();n++){ const DFDCPseudo *fdchit=segments[k]->hits[n]; fit.AddHit(fdchit); //bfield->GetField(x,y,z,Bx,By,Bz); //Bz_avg-=Bz; Bz+=bfield->GetBz(fdchit->xy.X(),fdchit->xy.Y(),fdchit->wire->origin.z()); } num_hits+=segments[k]->hits.size(); } Bz=fabs(Bz)/double(num_hits); // Fit the points to a circle return (fit.FitCircleRiemann(segments[0]->rc)); } // Method to match FDC and CDC candidates that projects the CDC candidate to the // cdc endplate, where the FDC candidates are reported. bool DTrackCandidate_factory::MatchMethod1(const DTrackCandidate *fdccan, vector &cdc_forward_ids, vector&cdc_endplate_projections, vector&used_cdc_hits){ // Magnitude of the momentum of the potential cdc match double p_cdc=0.; // Momentum and position vectors for the FDC candidate DVector3 mom=fdccan->momentum(); DVector3 pos=fdccan->position(); double p_fdc=mom.Mag(); ProjectHelixToZ(cdc_endplate.z(),fdccan->charge(),mom,pos); // loop over the cdc candidates looking for the smallest distance // between the cdc and fdc projections to the end plate double diff_min=1000.; // candidate matching difference unsigned int jmin=0; double radius=0.; for (unsigned int j=0;jmomentum().Mag(); } } if (DEBUG_HISTS){ match_dist->Fill(radius,diff_min); match_dist_vs_p->Fill(p_fdc,diff_min); } if (cdc_fdc_match(p_fdc,p_cdc,diff_min)){ // Get the segment data vectorsegments; fdccan->GetT(segments); // The following code snippet is intended to prevent matching a cdc // track candidate to a low momentum particle curling from the cdc endplate if (segments.size()>1){ double theta=180./M_PI*mom.Theta(); if (p_fdc<0.3 && theta<5.) return false; } // JANA does not maintain the order that the segments were added // as associated objects. Therefore, we need to reorder them here // so segment[0] is the most upstream segment. sort(segments.begin(), segments.end(), SegmentSortByLayerincreasing); unsigned int cdc_index=cdc_forward_ids[jmin]; const DTrackCandidate *cdccan=cdctrackcandidates[cdc_index]; // Get the associated cdc hits vectorcdchits; cdccan->GetT(cdchits); // Sort CDC hits by layer sort(cdchits.begin(), cdchits.end(), CDCHitSortByLayerincreasing); // Abort if the track would have to pass from one quadrant into another // quadrant (this is more likely two roughly back-to-back tracks rather than // a single track). const DVector2 fdc_hit_pos=segments[0]->hits[0]->xy; const DVector3 cdc_wire_origin=cdchits[0]->wire->origin; if (cdc_wire_origin.x()*fdc_hit_pos.X()<0. || cdc_wire_origin.y()*fdc_hit_pos.Y()<0.){ if (DEBUG_LEVEL>1) _DBG_ << "Skipping match of potential back-to-back tracks." <momentum().Theta(); fit.tanl=tan(M_PI_2-theta); if (GetPositionAndMomentum(fit,Bz_avg,cdchits[0]->wire->origin, pos,mom)==NOERROR){ // FDC hit const DFDCPseudo *fdchit=segments[0]->hits[0]; // Find the z-position at the new position in x and y DVector2 xy0(pos.X(),pos.Y()); double tworc=2.*fit.r0; double ratio=(fdchit->xy-xy0).Mod()/tworc; double sperp=(ratio<1.)?tworc*asin(ratio):tworc*M_PI_2; // Create new track candidate object DTrackCandidate *can = new DTrackCandidate; can->used_cdc_indexes=cdccan->used_cdc_indexes; // Add cdc and fdc hits to the track as associated objects for (unsigned int m=0;mhits.size();n++){ can->AddAssociatedObject(segments[m]->hits[n]); } } for (unsigned int n=0;nused_cdc_indexes[n]]=1; can->AddAssociatedObject(cdchits[n]); } pos.SetZ(fdchit->wire->origin.z()-sperp*fit.tanl); can->chisq=fit.chisq; can->Ndof=fit.ndof; can->setCharge(FactorForSenseOfRotation*fit.h); can->setMomentum(mom); can->setPosition(pos); trackcandidates.push_back(can); // Remove the CDC candidate from the id list because we // found a match cdc_forward_ids.erase(cdc_forward_ids.begin()+jmin); if (DEBUG_LEVEL>0) _DBG_ << ".. matched to CDC candidate #" << cdc_index < &cdc_forward_ids, vector&used_cdc_hits ){ // loop over the forward-pointing cdc candidates for (unsigned int j=0;jcdchits; cdccan->GetT(cdchits); sort(cdchits.begin(),cdchits.end(),CDCHitSortByLayerincreasing); // Variables to count the number of matching hits and the match fraction unsigned int num_match=0; unsigned int num_cdc=0; // Get the track parameters from the fdc candidate DVector3 pos=fdccan->position(); DVector3 mom=fdccan->momentum(); double q=fdccan->charge(); // loop over the cdc hits and count hits that agree with a projection of // the helix into the cdc for (unsigned int m=0;m0.01) num_match++; if (DEBUG_LEVEL>1) _DBG_ << "CDC s: " << cdchits[m]->wire->straw << " r: " << cdchits[m]->wire->ring << " prob: " << prob << endl; num_cdc++; } if (num_match>=3 && double(num_match)/double(num_cdc)>0.33){ // Put the fdc candidate in the combined list DTrackCandidate *can = new DTrackCandidate; can->setCharge(q); // copy the list of cdc indices over from the cdc candidate can->used_cdc_indexes=cdccan->used_cdc_indexes; // Get the segment data vectorsegments; fdccan->GetT(segments); // JANA does not maintain the order that the segments were added // as associated objects. Therefore, we need to reorder them here // so segment[0] is the most upstream segment. sort(segments.begin(), segments.end(), SegmentSortByLayerincreasing); // Abort if the track would have to pass from one quadrant into the opposite // quadrant (this is more likely two roughly back-to-back tracks rather than // a single track). const DVector2 fdc_hit_pos=segments[0]->hits[0]->xy; const DVector3 cdc_wire_origin=cdchits[0]->wire->origin; if (cdc_wire_origin.x()*fdc_hit_pos.X()<0. && cdc_wire_origin.y()*fdc_hit_pos.Y()<0.){ if (DEBUG_LEVEL>1) _DBG_ << "Skipping match of potential back-to-back tracks." <hits.size();n++){ can->AddAssociatedObject(segments[m]->hits[n]); num_fdc_hits++; } } // Add the CDC hits to the track candidate for (unsigned int m=0;mused_cdc_indexes[m]]=1; can->AddAssociatedObject(cdchits[m]); } // Average Bz double Bz_avg=0.; // Redo circle fit with additional hits DHelicalFit fit; if (DoRefit(fit,segments,cdchits,Bz_avg)==NOERROR){ // Determine the polar angle double theta=fdccan->momentum().Theta(); fit.tanl=tan(M_PI_2-theta); if (GetPositionAndMomentum(fit,Bz_avg,cdchits[0]->wire->origin, pos,mom)==NOERROR){ // FDC hit const DFDCPseudo *fdchit=segments[0]->hits[0]; // Find the z-position at the new position in x and y DVector2 xy0(pos.X(),pos.Y()); double tworc=2.*fit.r0; double ratio=(fdchit->xy-xy0).Mod()/tworc; double sperp=(ratio<1.)?tworc*asin(ratio):tworc*M_PI_2; pos.SetZ(fdchit->wire->origin.z()-sperp*fit.tanl); can->setCharge(FactorForSenseOfRotation*fit.h); } else{ //Set the charge for the stepper stepper->SetCharge(q); // Swim to a radius just outside the start counter stepper->SwimToRadius(pos,mom,9.5,NULL); } can->chisq=fit.chisq; can->Ndof=fit.ndof; can->setMomentum(mom); can->setPosition(pos); } else{ //_DBG_ << endl; can->Ndof=fdccan->Ndof; can->chisq=fdccan->chisq; can->setMomentum(fdccan->momentum()); can->setPosition(fdccan->position()); } trackcandidates.push_back(can); // Remove the CDC candidate from the list cdc_forward_ids.erase(cdc_forward_ids.begin()+j); if (DEBUG_LEVEL>0) _DBG_ << ".. matched to CDC candidate #" < &forward_matches, vector&used_cdc_hits ){ // Get hits already linked to this candidate from associated objects vectorcdchits; cdccan->GetT(cdchits); sort(cdchits.begin(), cdchits.end(), CDCHitSortByLayerincreasing); // loop over fdc candidates for (unsigned int k=0;kmomentum().Mag(); double p_cdc=cdccan->momentum().Mag(); if (p_fdc/p_cdc>0.5){ // Get the segment data vectorsegments; fdccan->GetT(segments); sort(segments.begin(), segments.end(), SegmentSortByLayerincreasing); // Abort if the track would have to pass from one quadrant into the opposite // quadrant (this is more likely two roughly back-to-back tracks rather than // a single track). const DVector2 fdc_hit_pos=segments[0]->hits[0]->xy; const DVector3 cdc_wire_origin=cdchits[0]->wire->origin; if (cdc_wire_origin.x()*fdc_hit_pos.X()<0. && cdc_wire_origin.y()*fdc_hit_pos.Y()<0.){ if (DEBUG_LEVEL>1) _DBG_ << "Skipping match of potential back-to-back tracks." <momentum(); DVector3 pos=cdccan->position(); double q=cdccan->charge(); // Try to match unmatched fdc candidates int num_match=0; int num_hits=0; for (unsigned int m=0;mhits.size();n++){ unsigned int ind=segments[m]->hits.size()-1-n; const DFDCPseudo *hit=segments[m]->hits[ind]; ProjectHelixToZ(hit->wire->origin.z(),q,mom,pos); // difference DVector2 XY=hit->xy; double dx=XY.X()-pos.x(); double dy=XY.Y()-pos.y(); double dr2=dx*dx+dy*dy; double variance=1.0; double prob = isfinite(dr2) ? TMath::Prob(dr2/variance,1):0.0; if (prob>0.01) num_match++; num_hits++; } if (num_match>=3 || (num_match>0 && double(num_match)/double(num_hits)>0.33)){ DTrackCandidate *can = new DTrackCandidate; can->setCharge(cdccan->charge()); can->used_cdc_indexes=cdccan->used_cdc_indexes; // mark the fdc track candidate as matched forward_matches[k]=1; // Add cdc hits to the candidate for (unsigned int m=0;mused_cdc_indexes[m]]=1; can->AddAssociatedObject(cdchits[m]); } // Add fdc hits to the candidate for (unsigned int m=0;mhits.size();n++){ can->AddAssociatedObject(segments[m]->hits[n]); } } // average Bz double Bz_avg=0.; // Redo helical fit with all available hits DHelicalFit fit; if (DoRefit(fit,segments,cdchits,Bz_avg)==NOERROR){ // Determine the polar angle double theta=fdccan->momentum().Theta(); fit.tanl=tan(M_PI_2-theta); // Redo the line fit fit.FitLineRiemann(); // Set the charge fit.FindSenseOfRotation(); can->setCharge(FactorForSenseOfRotation*fit.h); if (GetPositionAndMomentum(fit,Bz_avg,cdchits[0]->wire->origin, pos,mom)==NOERROR){ // FDC hit const DFDCPseudo *fdchit=segments[0]->hits[0]; // Find the z-position at the new position in x and y DVector2 xy0(pos.X(),pos.Y()); double tworc=2.*fit.r0; double ratio=(fdchit->xy-xy0).Mod()/tworc; double sperp=(ratio<1.)?tworc*asin(ratio):tworc*M_PI_2; pos.SetZ(fdchit->wire->origin.z()-sperp*fit.tanl); } else{ mom=fdccan->momentum(); pos=fdccan->position(); } can->chisq=fit.chisq; can->Ndof=fit.ndof; can->setMomentum(mom); can->setPosition(pos); } else{ //_DBG_ << endl; can->Ndof=cdccan->Ndof; can->chisq=cdccan->chisq; can->setMomentum(cdccan->momentum()); can->setPosition(cdccan->position()); } trackcandidates.push_back(can); if (DEBUG_LEVEL>0) _DBG_ << ".. matched to CDC candidate #" << k < &forward_matches, int &num_fdc_cands_remaining){ // Get segments already linked to this candidate from associated objects vectorsrc_segments; srccan->GetT(src_segments); if (src_segments.size()==1) return false; sort(src_segments.begin(), src_segments.end(), SegmentSortByLayerincreasing); if (DEBUG_LEVEL>0) _DBG_ << "Attempting matching method #4..." <package; unsigned int pack1_first=src_segments[0]->package; for (unsigned int k=0;ksegments; fdccan->GetT(segments); sort(segments.begin(), segments.end(), SegmentSortByLayerincreasing); const DFDCPseudo *firsthit=segments[0]->hits[0]; unsigned int pack2_first=segments[0]->package; unsigned int pack2_last=segments[segments.size()-1]->package; if (pack2_first>pack1_last || pack2_lastmomentum(); DVector3 pos=srccan->position(); double q=srccan->charge(); DVector3 norm(0.,0.,1.); if (pack2_lastSetCharge(q); stepper->SwimToPlane(pos,mom,firsthit->wire->origin,norm,NULL); double dx=pos.x()-firsthit->xy.X(); double dy=pos.y()-firsthit->xy.Y(); double d2=dx*dx+dy*dy; double variance=1.; double prob=TMath::Prob(d2/variance,1); unsigned int num_match=(prob>0.01)?1:0; // Try to match unmatched fdc candidates for (unsigned int i=1;ihits.size();i++){ const DFDCPseudo *hit=segments[0]->hits[i]; ProjectHelixToZ(hit->wire->origin.z(),q,mom,pos); DVector2 XY=hit->xy; double dx=XY.X()-pos.x(); double dy=XY.Y()-pos.y(); double dr2=dx*dx+dy*dy; double variance=1.0; double prob = isfinite(dr2) ? TMath::Prob(dr2/variance,1):0.0; if (prob>0.01) num_match++; } if (num_match>=3){ forward_matches[k]=1; // Create new candidate for output DTrackCandidate *can = new DTrackCandidate; // average Bz double Bz_avg=0.; unsigned int num_hits=0; // Create a new DHelicalFit object for fitting combined data DHelicalFit fit; // Add hits to the fit object and also to the track candidate // itself as associated objects for (unsigned int m=0;mhits.size();n++){ const DFDCPseudo *fdchit=src_segments[m]->hits[n]; fit.AddHit(fdchit); can->AddAssociatedObject(fdchit); //bfield->GetField(x,y,z,Bx,By,Bz); //Bz_avg-=Bz; Bz_avg+=bfield->GetBz(fdchit->xy.X(),fdchit->xy.Y(), fdchit->wire->origin.z()); } num_hits+=src_segments[m]->hits.size(); } for (unsigned int m=0;mhits.size();n++){ const DFDCPseudo *fdchit=segments[m]->hits[n]; fit.AddHit(fdchit); can->AddAssociatedObject(fdchit); //bfield->GetField(x,y,z,Bx,By,Bz); //Bz_avg-=Bz; Bz_avg+=bfield->GetBz(fdchit->xy.X(),fdchit->xy.Y(), fdchit->wire->origin.z()); } num_hits+=segments[m]->hits.size(); } // Fit the points to a circle if (fit.FitCircleRiemann(segments[0]->rc)==NOERROR){ // Compute new transverse momentum Bz_avg=fabs(Bz_avg)/double(num_hits); // Guess for theta and z from input candidates double theta=fdccan->momentum().Theta(); fit.tanl=tan(M_PI_2-theta); fit.z_vertex=srccan->position().Z(); // Redo line fit fit.FitLineRiemann(); // Guess charge from fit fit.h=GetSenseOfRotation(fit,firsthit,srccan->position()); can->setCharge(FactorForSenseOfRotation*fit.h); // put z position just upstream of the first hit in z const DHFHit_t *myhit=fit.GetHits()[0]; pos.SetXYZ(myhit->x,myhit->y,myhit->z); GetPositionAndMomentum(myhit->z-1.,fit,Bz_avg,pos,mom); can->setPosition(pos); can->setMomentum(mom); can->chisq=fit.chisq; can->Ndof=fit.ndof; trackcandidates.push_back(can); num_fdc_cands_remaining--; if (DEBUG_LEVEL>0) _DBG_ << "Found a match using method #4" <&cdchits, vector &forward_matches ){ // Loop over the fdc track candidates for (unsigned int k=0;kposition(); DVector3 mom=fdccan->momentum(); double q=fdccan->charge(); for (unsigned int m=0;m0.01) num_match++; if (DEBUG_LEVEL>1) _DBG_ << "CDC s: " << cdchits[m]->wire->straw << " r: " << cdchits[m]->wire->ring << " prob: " << prob << endl; num_cdc++; } if (num_match>=3 && double(num_match)/double(num_cdc)>0.33){ // Get the segment data vectorsegments; fdccan->GetT(segments); // JANA does not maintain the order that the segments were added // as associated objects. Therefore, we need to reorder them here // so segment[0] is the most upstream segment. sort(segments.begin(), segments.end(), SegmentSortByLayerincreasing); // Abort if the track would have to pass from one quadrant into the opposite // quadrant (this is more likely two roughly back-to-back tracks rather than // a single track). const DVector2 fdc_hit_pos=segments[0]->hits[0]->xy; const DVector3 cdc_wire_origin=cdchits[0]->wire->origin; if (cdc_wire_origin.x()*fdc_hit_pos.X()<0. && cdc_wire_origin.y()*fdc_hit_pos.Y()<0.){ if (DEBUG_LEVEL>1) _DBG_ << "Skipping match of potential back-to-back tracks." <0) _DBG_ << "... matched to FDC candidate #" << k <hits.size();n++){ can->AddAssociatedObject(segments[m]->hits[n]); } } // Average Bz double Bz_avg=0.; // Instantiate the helical fitter to do the refit DHelicalFit fit; if (DoRefit(fit,segments,cdchits,Bz_avg)==NOERROR){ // Determine the polar angle double theta=fdccan->momentum().Theta(); fit.tanl=tan(M_PI_2-theta); if (GetPositionAndMomentum(fit,Bz_avg,cdchits[0]->wire->origin, pos,mom)==NOERROR){ // FDC hit const DFDCPseudo *fdchit=segments[0]->hits[0]; // Find the z-position at the new position in x and y DVector2 xy0(pos.X(),pos.Y()); double tworc=2.*fit.r0; double ratio=(fdchit->xy-xy0).Mod()/tworc; double sperp=(ratio<1.)?tworc*asin(ratio):tworc*M_PI_2; pos.SetZ(fdchit->wire->origin.z()-sperp*fit.tanl); can->setCharge(FactorForSenseOfRotation*fit.h); } can->chisq=fit.chisq; can->Ndof=fit.ndof; can->setMomentum(mom); can->setPosition(pos); } return true; } // check for match } // check that the candidate has not already been matched } // loop over fdc candidates return false; } // Method to match FDC candidates with stray CDC hits unassociated with track // candidates void DTrackCandidate_factory::MatchMethod6(DTrackCandidate *can, vector&segments, vector&used_cdc_hits, unsigned int &num_unmatched_cdcs ){ // Input candidate parameters DVector3 pos=can->position(); DVector3 mom=can->momentum(); double q=can->charge(); // Variables to keep track of cdc hits unsigned int num_match_cdc=0; unsigned int num_axial=0; int id_for_smallest_dr=-1; int old_ring=-1; vectormatched_ids; double dr2_min=1e6; for (unsigned int k=0;kwire->ring!=old_ring){ if (id_for_smallest_dr>=0){ if (!used_cdc_hits[id_for_smallest_dr]){ num_match_cdc++; num_unmatched_cdcs--; used_cdc_hits[id_for_smallest_dr]=1; can->used_cdc_indexes.push_back(id_for_smallest_dr); can->AddAssociatedObject(mycdchits[id_for_smallest_dr]); if (mycdchits[id_for_smallest_dr]->is_stereo==false){ num_axial++; matched_ids.push_back(id_for_smallest_dr); } } } // Reset matching variables dr2_min=1e6; id_for_smallest_dr=-1; } if (!used_cdc_hits[k]){ double variance=1.0; // Use a helical approximation to the track to match both axial and // stereo wires double dr2=DocaToHelix(mycdchits[k],q,pos,mom); double prob=isfinite(dr2) ? TMath::Prob(dr2/variance,1):0.0; if (prob>0.5){ // Look for closest match if (dr2wire->ring; } // Grab the last potential match if (id_for_smallest_dr>=0){ if (!used_cdc_hits[id_for_smallest_dr]){ num_match_cdc++; num_unmatched_cdcs--; used_cdc_hits[id_for_smallest_dr]=1; can->used_cdc_indexes.push_back(id_for_smallest_dr); can->AddAssociatedObject(mycdchits[id_for_smallest_dr]); if (mycdchits[id_for_smallest_dr]->is_stereo==false){ num_axial++; matched_ids.push_back(id_for_smallest_dr); } } } // We found some matched hits. Update the start position of the candidate. if (num_match_cdc>0){ const DFDCPseudo *firsthit=segments[0]->hits[0]; int axial_id=-1; if (num_axial>0) axial_id=matched_ids[0]; // Compute the average Bz double Bz_avg=0.,denom=0.; for (unsigned int m=0;mhits.size();n++){ const DFDCPseudo *fdchit=segments[m]->hits[n]; double x=fdchit->xy.X(); double y=fdchit->xy.Y(); double z=fdchit->wire->origin.z(); Bz_avg+=bfield->GetBz(x,y,z); denom+=1.; } } Bz_avg=fabs(Bz_avg)/denom; // Store the current track parameters in the DHelicalFit class DHelicalFit fit; fit.r0=mom.Perp()/(0.003*Bz_avg); fit.phi=mom.Phi(); fit.h=q*FactorForSenseOfRotation; fit.x0=pos.x()-fit.h*fit.r0*sin(fit.phi); fit.y0=pos.y()+fit.h*fit.r0*cos(fit.phi); fit.tanl=tan(M_PI_2-mom.Theta()); fit.z_vertex=0; // this will be changed later // Update the position and momentum of the candidate based on the most // recent fit results UpdatePositionAndMomentum(can,firsthit,fit,Bz_avg,axial_id); if (DEBUG_LEVEL>0) _DBG_ << "... matched stray CDC hits ..." << endl; cout<< "... matched stray CDC hits ..." < &forward_matches, int &num_fdc_cands_remaining){ if (DEBUG_LEVEL>0) _DBG_ << "Attempting matching method #7..." <fdchits; srccan->GetT(fdchits); if (fdchits.size()==0) return false; sort(fdchits.begin(),fdchits.end(),FDCHitSortByLayerincreasing); unsigned int pack1_last=(fdchits[fdchits.size()-1]->wire->layer-1)/6; for (unsigned int k=0;ksegments; fdccan->GetT(segments); sort(segments.begin(), segments.end(), SegmentSortByLayerincreasing); const DFDCPseudo *firsthit=segments[0]->hits[0]; unsigned int pack2_first=segments[0]->package; if (pack2_first-pack1_last==1){ // match in adjacent packages // Momentum and position vectors for the input candidate DVector3 mom=srccan->momentum(); DVector3 pos=srccan->position(); double q=srccan->charge(); DVector3 norm(0.,0.,1.); // Swim to the first hit in the candidate to which we wish to match // using the stepper stepper->SetCharge(q); stepper->SwimToPlane(pos,mom,firsthit->wire->origin,norm,NULL); double dx=pos.x()-firsthit->xy.X(); double dy=pos.y()-firsthit->xy.Y(); double d2=dx*dx+dy*dy; double variance=1.; double prob=TMath::Prob(d2/variance,1); unsigned int num_match=(prob>0.01)?1:0; // Try to match unmatched fdc candidates for (unsigned int i=1;ihits.size();i++){ const DFDCPseudo *hit=segments[0]->hits[i]; ProjectHelixToZ(hit->wire->origin.z(),q,mom,pos); DVector2 XY=hit->xy; double dx=XY.X()-pos.x(); double dy=XY.Y()-pos.y(); double dr2=dx*dx+dy*dy; double variance=1.0; double prob = isfinite(dr2) ? TMath::Prob(dr2/variance,1):0.0; if (prob>0.01) num_match++; } if (num_match>=3){ forward_matches[k]=1; num_fdc_cands_remaining--; // average Bz double Bz=0.; unsigned int num_hits=0; // Create a new DHelicalFit object for fitting combined data DHelicalFit fit; // Get the cdc hits associated with this track candidate and add the // axial straws to the list of hits to use in the circle fit vectorcdchits; srccan->GetT(cdchits); for (unsigned int i=0;iis_stereo==false){ double cov=0.213; //guess const DVector3 origin=cdchits[i]->wire->origin; double x=origin.x(),y=origin.y(),z=origin.z(); fit.AddHitXYZ(x,y,z,cov,cov,0.,true); Bz+=bfield->GetBz(x,y,z); num_hits++; } } // Add the FDC hits from the existing track to this list for (unsigned int i=0;iwire->origin.z(); Bz+=bfield->GetBz(fdchits[i]->xy.X(),fdchits[i]->xy.Y(),zhit); num_hits++; if (zhitwire->origin.z()){ firsthit=fdchits[i]; } } // Add the new set of FDC hits; also add them as associated objects // to the track for (unsigned int m=0;mhits.size();n++){ const DFDCPseudo *fdchit=segments[m]->hits[n]; fit.AddHit(fdchit); double zhit=fdchit->wire->origin.z(); Bz+=bfield->GetBz(fdchit->xy.X(),fdchit->xy.Y(),zhit); srccan->AddAssociatedObject(fdchit); } num_hits+=segments[m]->hits.size(); } Bz=fabs(Bz)/double(num_hits); // Redo the circle fit with the extra hits double rc=srccan->momentum().Perp()/(0.003*Bz); if (fit.FitCircleRiemann(rc)==NOERROR){ // Determine the polar angle double theta=srccan->momentum().Theta(); fit.tanl=tan(M_PI_2-theta); if (cdchits.size()>0){ if (GetPositionAndMomentum(fit,Bz,cdchits[0]->wire->origin, pos,mom)==NOERROR){ // Find the z-position at the new position in x and y DVector2 xy0(pos.X(),pos.Y()); double tworc=2.*fit.r0; double ratio=(firsthit->xy-xy0).Mod()/tworc; double sperp=(ratio<1.)?tworc*asin(ratio):tworc*M_PI_2; pos.SetZ(firsthit->wire->origin.z()-sperp*fit.tanl); } } else{ // put z position just upstream of the first hit in z const DHFHit_t *myhit=fit.GetHits()[0]; pos.SetXYZ(myhit->x,myhit->y,myhit->z); GetPositionAndMomentum(myhit->z-1.,fit,Bz,pos,mom); } srccan->chisq=fit.chisq; srccan->Ndof=fit.ndof; srccan->setMomentum(mom); srccan->setPosition(pos); } if (DEBUG_LEVEL>0) _DBG_ << "Found a match using method #7" < &forward_matches, vector&used_cdc_hits ){ // Get the cdc hits associated with this track candidate vectorcdchits; cdccan->GetT(cdchits); sort(cdchits.begin(), cdchits.end(), CDCHitSortByLayerincreasing); unsigned int last_index=cdchits.size()-1; if (cdchits[last_index]->is_stereo==true && cdchits[last_index]->wire->ring<13){ DHelicalFit fit; // Assume the outermost stereo hit occurs near the end of the straw DVector3 origin=cdchits[last_index]->wire->origin; DVector3 dir=cdchits[last_index]->wire->udir; DVector3 wirepos=origin+(75./dir.z())*dir; double cov=0.2; double x=wirepos.x(),y=wirepos.y(),z=wirepos.z(); fit.AddHitXYZ(x,y,z,cov,cov,0,true); double Bz=bfield->GetBz(x,y,z); int num_hits=1; // Add the cdc axial wires to the list of hits to use in the fit for (unsigned int k=0;kis_stereo==false){ cov=0.2; //guess origin=cdchits[k]->wire->origin; double x=origin.x(),y=origin.y(),z=origin.z(); fit.AddHitXYZ(x,y,z,cov,cov,0.,true); Bz+=bfield->GetBz(x,y,z); num_hits++; } } Bz=fabs(Bz)/num_hits; // Fit the points to a circle assuming that the circle goes through the // origin. fit.FitCircle(); // Determine the tangent of the dip angle double tworc=2.*fit.r0; double sperp=tworc*asin(wirepos.Perp()/tworc); fit.tanl=(wirepos.z()-TARGET_Z)/sperp; // Find sense of rotation (proportional to the charge) fit.GuessChargeFromCircleFit(); double q=fit.h*FactorForSenseOfRotation; // Find the momentum of the particle and the position just outside the // start counter DVector3 pos,mom; if (GetPositionAndMomentum(fit,Bz,cdchits[0]->wire->origin,pos,mom) ==NOERROR){ // Find the z-position at the new position in x and y DVector2 xy0(pos.X(),pos.Y()); double tworc=2.*fit.r0; double ratio=(wirepos-pos).Perp()/tworc; double sperp=(ratio<1.)?tworc*asin(ratio):tworc*M_PI_2; pos.SetZ(wirepos.z()-sperp*fit.tanl); // Set the charge for the stepper stepper->SetCharge(q); // Vector normal to the fdc planes DVector3 norm(0.,0.,1.); // Loop over the fdc candidates looking for a match for (unsigned int k=0;ksegments; fdccan->GetT(segments); sort(segments.begin(), segments.end(), SegmentSortByLayerincreasing); const DFDCPseudo *firsthit=segments[0]->hits[0]; // Make copies of the momentum and position vectors for input to // the swimmer DVector3 my_mom(mom); DVector3 my_pos(pos); // Swim to the first hit in the candidate to which we wish to match // using the stepper stepper->SwimToPlane(my_pos,my_mom,firsthit->wire->origin,norm,NULL); // Match based on proximity of projection to hit double dx=my_pos.x()-firsthit->xy.X(); double dy=my_pos.y()-firsthit->xy.Y(); double d2=dx*dx+dy*dy; double variance=1.; double prob=TMath::Prob(d2/variance,1); unsigned int num_match=(prob>0.01)?1:0; // Try to match further hits in most upstream segment for (unsigned int i=1;ihits.size();i++){ const DFDCPseudo *hit=segments[0]->hits[i]; ProjectHelixToZ(hit->wire->origin.z(),q,my_mom,my_pos); DVector2 XY=hit->xy; double dx=XY.X()-my_pos.x(); double dy=XY.Y()-my_pos.y(); double dr2=dx*dx+dy*dy; double variance=1.0; double prob = isfinite(dr2) ? TMath::Prob(dr2/variance,1):0.0; if (prob>0.01) num_match++; } if (num_match>=3){ forward_matches[k]=1; // Keep track of magnetic field in FDC double Bz_fdc=0; unsigned int num_hits_fdc=0; // Drop the first cdc hit in the list, since it isn't an // axial straw so the actual x,y position of the wire // depends on z, which we don't really know yet. fit.PruneHit(0); // Add the fdc hits to the fit class for (unsigned int m=0;mhits.size();n++){ const DFDCPseudo *fdchit=segments[m]->hits[n]; fit.AddHit(fdchit); Bz_fdc+=bfield->GetBz(fdchit->xy.X(),fdchit->xy.Y(), fdchit->wire->origin.z()); } num_hits_fdc+=segments[m]->hits.size(); } // Fake point at origin fit.AddHitXYZ(0.,0.,TARGET_Z,BEAM_VAR,BEAM_VAR,0.,true); // Fit the points to a circle if (fit.FitCircleRiemann(segments[0]->rc)==NOERROR){ // Use the fdc track candidate to get tanl double theta=fdccan->momentum().Theta(); fit.tanl=tan(M_PI_2-theta); Bz=0.5*(Bz+fabs(Bz_fdc)/num_hits_fdc); if (GetPositionAndMomentum(fit,Bz,cdchits[0]->wire->origin, my_pos,my_mom)==NOERROR){ // FDC hit const DFDCPseudo *fdchit=segments[0]->hits[0]; // Find the z-position at the new position in x and y DVector2 xy0(my_pos.X(),my_pos.Y()); double tworc=2.*fit.r0; double ratio=(fdchit->xy-xy0).Mod()/tworc; double sperp=(ratio<1.)?tworc*asin(ratio):tworc*M_PI_2; my_pos.SetZ(fdchit->wire->origin.z()-sperp*fit.tanl); DTrackCandidate *can = new DTrackCandidate; can->setCharge(q); can->setMomentum(my_mom); can->setPosition(my_pos); can->chisq=fit.chisq; can->Ndof=fit.ndof; // Add hits as associated objects for (unsigned int m=0;mhits.size();n++){ const DFDCPseudo *fdchit=segments[m]->hits[n]; can->AddAssociatedObject(fdchit); } } for (unsigned int n=0;nused_cdc_indexes[n]]=1; can->AddAssociatedObject(cdchits[n]); } trackcandidates.push_back(can); if (DEBUG_LEVEL>0) _DBG_ << "Matched using Method #8" <&cands, vector &forward_matches){ double q=srccan->charge(); int pack1=segment->package; // Get hits from segment and redo fit forcing the circle to go through (0,0) DHelicalFit fit1; for (unsigned int n=0;nhits.size();n++){ const DFDCPseudo *hit=segment->hits[n]; fit1.AddHit(hit); } fit1.FitCircle(); // Sense of rotation fit1.h=q*FactorForSenseOfRotation; // Use the fdc track candidate to get tanl double theta=srccan->momentum().Theta(); fit1.tanl=tan(M_PI_2-theta); // Find the magnetic field at the first hit in the first segment double x=segment->hits[0]->xy.X(); double y=segment->hits[0]->xy.Y(); double z=segment->hits[0]->wire->origin.z(); double Bz=fabs(bfield->GetBz(x,y,z)); // Get position and momentum for the first segment DVector3 mypos,mymom; GetPositionAndMomentum(z,fit1,Bz,mypos,mymom); // Loop over the rest of the fdc track candidates, skipping those that have // already been used for (unsigned int k=src_index+1;ksegments2; can2->GetT(segments2); int pack2=segments2[0]->package; if (abs(pack1-pack2)>0){ // Get hits from the second segment and redo fit forcing circle // to go through (0,0) DHelicalFit fit2; for (unsigned int n=0;nhits.size();n++){ const DFDCPseudo *hit=segments2[0]->hits[n]; fit2.AddHit(hit); } fit2.FitCircle(); // Match using centers of circles double dx=fit1.x0-fit2.x0; double dy=fit1.y0-fit2.y0; double circle_center_diff2=dx*dx+dy*dy; double got_match=false; if (circle_center_diff2<4.0) got_match=true; // try another matching method here if got_match==false else{ // Sense of rotation double q2=can2->charge(); fit2.h=q2*FactorForSenseOfRotation; // Use the fdc track candidate to get tanl theta=can2->momentum().Theta(); fit2.tanl=tan(M_PI_2-theta); // Try to match segments by swimming through the field got_match=MatchMethod11(q,mypos,mymom,fit2,segment,segments2[0]); } if (got_match){ forward_matches[k]=1; forward_matches[src_index]=1; // Create a new DTrackCandidate for output DTrackCandidate *can = new DTrackCandidate; // variables for finding double Bz=0; int num_hits=0; // Add hits from first segment as associated objects for (unsigned int n=0;nhits.size();n++){ const DFDCPseudo *fdchit=segment->hits[n]; can->AddAssociatedObject(fdchit); Bz+=bfield->GetBz(fdchit->xy.X(),fdchit->xy.Y(), fdchit->wire->origin.z()); num_hits++; } // Add the hits from the second segment to the first fit object // and refit the circle. Also add them as associated objects // to the new candidate. for (unsigned int n=0;nhits.size();n++){ const DFDCPseudo *hit=segments2[0]->hits[n]; fit1.AddHit(hit); can->AddAssociatedObject(hit); Bz+=bfield->GetBz(hit->xy.X(),hit->xy.Y(), hit->wire->origin.z()); num_hits++; } Bz=fabs(Bz)/double(num_hits); // Initialize variables needed for output DVector3 mom=srccan->momentum(); DVector3 pos=srccan->position(); can->Ndof=srccan->Ndof; can->chisq=srccan->chisq; if (fit1.FitCircleRiemann(fit1.r0)==NOERROR){ can->Ndof=fit1.ndof; can->chisq=fit1.chisq; // Redo line fit fit1.FitLineRiemann(); // Guess charge from fit fit1.h=GetSenseOfRotation(fit1,segments2[0]->hits[0], srccan->position()); q=FactorForSenseOfRotation*fit1.h; // put z position just upstream of the first hit in z const DHFHit_t *myhit=fit1.GetHits()[0]; pos.SetXYZ(myhit->x,myhit->y,myhit->z); GetPositionAndMomentum(myhit->z-1.,fit1,Bz,pos,mom); } can->setCharge(q); can->setPosition(pos); can->setMomentum(mom); trackcandidates.push_back(can); return true; } // circle center match } // different packages? } // already matched? } // loop over tracks return false; } // Method to try to match unmatched FDC segments by assuming that the tracks // are sufficiently stiff we are better served by somethign closer to a // "straight-line" approximation bool DTrackCandidate_factory::MatchMethod10(unsigned int src_index, const DTrackCandidate *srccan, const DFDCSegment *segment, vector&cands, vector &forward_matches){ double q=srccan->charge(); int pack1=segment->package; // Get hits from segment and redo fit forcing the circle to go through (0,0) DHelicalFit fit1; for (unsigned int n=0;nhits.size();n++){ const DFDCPseudo *hit=segment->hits[n]; fit1.AddHit(hit); } fit1.FitCircleStraightTrack(); // Sense of rotation fit1.h=q*FactorForSenseOfRotation; // Use the fdc track candidate to get tanl double theta=srccan->momentum().Theta(); fit1.tanl=tan(M_PI_2-theta); // Find the magnetic field at the first hit in the first segment double x=segment->hits[0]->xy.X(); double y=segment->hits[0]->xy.Y(); double z=segment->hits[0]->wire->origin.z(); double Bz=fabs(bfield->GetBz(x,y,z)); // Get position and momentum for the first segment DVector3 mypos,mymom; GetPositionAndMomentum(z,fit1,Bz,mypos,mymom); // Loop over the rest of the fdc track candidates, skipping those that have // already been used for (unsigned int k=src_index+1;ksegments2; can2->GetT(segments2); int pack2=segments2[0]->package; if (abs(pack1-pack2)>0){ // Get hits from the second segment and redo fit forcing circle // to go through (0,0) DHelicalFit fit2; for (unsigned int n=0;nhits.size();n++){ const DFDCPseudo *hit=segments2[0]->hits[n]; fit2.AddHit(hit); } fit2.FitCircleStraightTrack(); // Sense of rotation double q2=can2->charge(); fit2.h=q2*FactorForSenseOfRotation; // Use the fdc track candidate to get tanl theta=can2->momentum().Theta(); fit2.tanl=tan(M_PI_2-theta); // Try to match segments by swimming through the field if (MatchMethod11(q,mypos,mymom,fit2,segment,segments2[0])){ forward_matches[k]=1; forward_matches[src_index]=1; // Create a new DTrackCandidate for output DTrackCandidate *can = new DTrackCandidate; // variables for finding double Bz=0; int num_hits=0; // Add hits from first segment as associated objects for (unsigned int n=0;nhits.size();n++){ const DFDCPseudo *fdchit=segment->hits[n]; can->AddAssociatedObject(fdchit); Bz+=bfield->GetBz(fdchit->xy.X(),fdchit->xy.Y(), fdchit->wire->origin.z()); num_hits++; } // Add the hits from the second segment to the first fit object // and refit the circle. Also add them as associated objects // to the new candidate. for (unsigned int n=0;nhits.size();n++){ const DFDCPseudo *hit=segments2[0]->hits[n]; fit1.AddHit(hit); can->AddAssociatedObject(hit); Bz+=bfield->GetBz(hit->xy.X(),hit->xy.Y(), hit->wire->origin.z()); num_hits++; } Bz=fabs(Bz)/double(num_hits); // Initialize variables needed for output DVector3 mom=srccan->momentum(); DVector3 pos=srccan->position(); double q=srccan->charge(); can->Ndof=srccan->Ndof; can->chisq=srccan->chisq; if (fit1.FitCircleRiemann(fit1.r0)==NOERROR){ can->Ndof=fit1.ndof; can->chisq=fit1.chisq; // Redo line fit fit1.FitLineRiemann(); // Guess charge from fit fit1.h=GetSenseOfRotation(fit1,segments2[0]->hits[0], srccan->position()); q=FactorForSenseOfRotation*fit1.h; // put z position just upstream of the first hit in z const DHFHit_t *myhit=fit1.GetHits()[0]; pos.SetXYZ(myhit->x,myhit->y,myhit->z); GetPositionAndMomentum(myhit->z-1.,fit1,Bz,pos,mom); } can->setCharge(q); can->setPosition(pos); can->setMomentum(mom); trackcandidates.push_back(can); return true; } // minimum number of matching hits } // different packages? } // already matched? } // loop over tracks return false; } // Swims from one FDC segment to another segment looking for a match. If this // fails, swims from the second second to the first in the opposite direction. bool DTrackCandidate_factory::MatchMethod11(double q,DVector3 &mypos, DVector3 &mymom, DHelicalFit &fit2, const DFDCSegment *segment1, const DFDCSegment *segment2 ){ // Package numbers int pack1=segment1->package; int pack2=segment2->package; // Set direction of propagation for traversing from segment1 to segment2 if ((pack20.) || (pack1>pack2 && mymom.z()<0.)) mymom=(-1.)*mymom; // Find the magnetic field at the first hit of the second segment double x=segment2->hits[0]->xy.X(); double y=segment2->hits[0]->xy.Y(); double z=segment2->hits[0]->wire->origin.z(); double Bz=fabs(bfield->GetBz(x,y,z)); // Get position and momentum for the second segment DVector3 mypos2,mymom2; GetPositionAndMomentum(z,fit2,Bz,mypos2,mymom2); if (pack2>pack1) mymom2=(-1.)*mymom2; // Swim from the first segment to the second segment using the stepper const DFDCPseudo *secondhit=segment2->hits[0]; DVector3 norm(0.,0.,1.); stepper->SetCharge(q); stepper->SwimToPlane(mypos,mymom,secondhit->wire->origin,norm,NULL); // Look for match at first hit double dx=mypos.x()-secondhit->xy.X(); double dy=mypos.y()-secondhit->xy.Y(); double d2=dx*dx+dy*dy; double variance=1.; double prob=TMath::Prob(d2/variance,1); unsigned int num_match=(prob>0.01)?1:0; // Try to match more hits for (unsigned int i=1;ihits.size();i++){ const DFDCPseudo *hit=segment2->hits[i]; ProjectHelixToZ(hit->wire->origin.z(),q,mymom,mypos); DVector2 XY=hit->xy; double dx=XY.X()-mypos.x(); double dy=XY.Y()-mypos.y(); double dr2=dx*dx+dy*dy; double variance=1.0; double prob = isfinite(dr2) ? TMath::Prob(dr2/variance,1):0.0; if (prob>0.01) num_match++; } if (num_match>=3){ return true; } // Try swimming in opposite direction secondhit=segment1->hits[0]; double q2=fit2.h*FactorForSenseOfRotation; stepper->SetCharge(q2); stepper->SwimToPlane(mypos2,mymom2,secondhit->wire->origin,norm,NULL); // Look for match at first hit dx=mypos2.x()-secondhit->xy.X(); dy=mypos2.y()-secondhit->xy.Y(); d2=dx*dx+dy*dy; prob=TMath::Prob(d2/variance,1); num_match=(prob>0.01)?1:0; // Try to match more hits for (unsigned int i=1;ihits.size();i++){ const DFDCPseudo *hit=segment1->hits[i]; ProjectHelixToZ(hit->wire->origin.z(),q2,mymom2,mypos2); DVector2 XY=hit->xy; double dx=XY.X()-mypos2.x(); double dy=XY.Y()-mypos2.y(); double dr2=dx*dx+dy*dy; double variance=1.0; double prob = isfinite(dr2) ? TMath::Prob(dr2/variance,1):0.0; if (prob>0.01) num_match++; } if (num_match>=3){ return true; } return false; } // Update the momentum and position entries for the candidate based on the // results of the most recent helical fit void DTrackCandidate_factory::UpdatePositionAndMomentum(DTrackCandidate *can, const DFDCPseudo *fdchit, DHelicalFit &fit, double Bz_avg, int axial_id){ DVector3 pos=can->position(),mom=can->momentum(); double q=can->charge(); if (GetPositionAndMomentum(fit,Bz_avg, pos,mom)==NOERROR){ // Find the z-position at the new position in x and y DVector2 xy0(pos.X(),pos.Y()); double tworc=2.*fit.r0; double ratio=(fdchit->xy-xy0).Mod()/tworc; double sperp=(ratio<1.)?tworc*asin(ratio):tworc*M_PI_2; pos.SetZ(fdchit->wire->origin.z()-sperp*fit.tanl); // update the track parameters can->setMomentum(mom); can->setPosition(pos); } else if (axial_id>=0){ // if the previous attempt to get the position at r^2=90 did not work, // then the circle after the refit does not intercept this radius. // Use these previous estimates before the refit for the track parameters // and the inner-most axial wire to estimate the position. fit.r0=mom.Perp()/(0.003*Bz_avg); fit.phi=mom.Phi(); fit.h=q*FactorForSenseOfRotation; fit.x0=pos.x()-fit.h*fit.r0*sin(fit.phi); fit.y0=pos.y()+fit.h*fit.r0*cos(fit.phi); fit.tanl=tan(M_PI_2-mom.Theta()); fit.z_vertex=0; // this will be changed later if (GetPositionAndMomentum(fit,Bz_avg,mycdchits[axial_id]->wire->origin,pos, mom)==NOERROR){ // Find the z-position at the new position in x and y DVector2 xy0(pos.X(),pos.Y()); double tworc=2.*fit.r0; double ratio=(fdchit->xy-xy0).Mod()/tworc; double sperp=(ratio<1.)?tworc*asin(ratio):tworc*M_PI_2; pos.SetZ(fdchit->wire->origin.z()-sperp*fit.tanl); // update the track parameters can->setMomentum(mom); can->setPosition(pos); } } } // Routine to return momentum and position given the helical parameters and the // z-component of the magnetic field jerror_t DTrackCandidate_factory::GetPositionAndMomentum(double z,DHelicalFit &fit, double Bz,DVector3 &pos, DVector3 &mom){ double xc=fit.x0; double yc=fit.y0; double rc=fit.r0; // Position double phi1=atan2(pos.y()-yc,pos.x()-xc); double q_over_rc_tanl=FactorForSenseOfRotation*fit.h/(rc*fit.tanl); double dphi_s=(pos.z()-z)*q_over_rc_tanl; double dphi1=phi1-dphi_s;// was - double x=xc+rc*cos(dphi1); double y=yc+rc*sin(dphi1); pos.SetXYZ(x,y,z); dphi1*=-1.; if (fit.h<0) dphi1+=M_PI; // Momentum double pt=0.003*fabs(Bz)*rc; double px=pt*sin(dphi1); double py=pt*cos(dphi1); double pz=pt*fit.tanl; mom.SetXYZ(px,py,pz); return NOERROR; }