61 const std::vector<DTRecHit1DPair>& pairs){
64 vector<DTHitPairForFit*> hitsForFit =
initHits(sl, pairs);
66 vector<DTSegmentCand*> candidates =
buildSegments(sl, hitsForFit);
68 vector<DTSegmentCand*>::const_iterator cand=candidates.begin();
69 while (cand<candidates.end()) {
76 cout<<
"Reconstructed 2D segments "<<*segment<<endl;
91 vector<DTHitPairForFit*>
93 const std::vector<DTRecHit1DPair>& hits){
95 vector<DTHitPairForFit*>
result;
96 for (vector<DTRecHit1DPair>::const_iterator
hit=hits.begin();
103 vector<DTSegmentCand*>
105 const std::vector<DTHitPairForFit*>& hits){
107 typedef vector<DTHitPairForFit*> hitCont;
108 typedef hitCont::const_iterator hitIter;
109 vector<DTSegmentCand*>
result;
113 cout <<
"buildSegments: " << sl->
id() <<
" nHits " << hits.size() << endl;
114 for (hitIter
hit=hits.begin();
hit!=hits.end(); ++
hit)
cout << **
hit<<
" wire: " << (*hit)->
id() <<
" DigiTime: " << (*hit)->digiTime() << endl;
122 cout <<
"Warning: this SuperLayer " << sl->
id() <<
" has too many hits : "
124 cout <<
"Skipping segment reconstruction... " << endl;
131 if ((sl->
id()).superlayer()==2)
136 vector<AssPoint> usedHits;
140 for (hitCont::const_iterator firstHit=hits.begin(); firstHit!=hits.end();
142 for (hitCont::const_reverse_iterator lastHit=hits.rbegin();
143 (*lastHit)!=(*firstHit); ++lastHit) {
146 if (!
geometryFilter((*firstHit)->id(),(*lastHit)->id()))
continue;
150 for (hitCont::const_iterator tmpHit=firstHit+1; (*tmpHit)!=(*lastHit); tmpHit++)
152 && (
geometryFilter((*tmpHit)->id(),(*firstHit)->id()))) hitsForFit.push_back(*tmpHit);
154 for (
int firstLR=0; firstLR<2; ++firstLR) {
155 for (
int lastLR=0; lastLR<2; ++lastLR) {
165 float DAlpha=fabs(gvec.
theta()-gvecIP.
theta());
166 if (DAlpha>DAlphaMax)
continue;
169 cout <<
"Selected hit pair:" << endl;
170 cout <<
" First " << *(*firstHit) <<
" Layer Id: " << (*firstHit)->id().layerId() <<
" Side: " << firstLR <<
" DigiTime: " << (*firstHit)->digiTime() << endl;
171 cout <<
" Last " << *(*lastHit) <<
" Layer Id: " << (*lastHit)->id().layerId() <<
" Side: " << lastLR <<
" DigiTime: " << (*lastHit)->digiTime() << endl;
174 vector<AssPoint> assHits;
176 assHits.push_back(
AssPoint(*firstHit,codes[firstLR]));
177 assHits.push_back(
AssPoint(*lastHit,codes[lastLR]));
181 addHits(sl,assHits,hitsForFit,result, usedHits);
189 cout <<
"Result (before cleaning): " << result.size() << endl;
190 for (vector<DTSegmentCand*>::const_iterator seg=result.begin(); seg!=result.end(); ++seg)
191 cout << *(*seg) << endl;
195 cout <<
"Result (after cleaning): " << result.size() << endl;
196 for (vector<DTSegmentCand*>::const_iterator seg=result.begin();
197 seg!=result.end(); ++seg)
198 cout << *(*seg) << endl;
208 vector<DTSegmentCand*> &
result, vector<AssPoint>& usedHits) {
210 typedef vector<DTHitPairForFit*> hitCont;
211 double chi2l,chi2r,t0_corrl,t0_corrr;
212 bool foundSomething =
false;
215 cout <<
"DTMeantimerPatternReco::addHit " << endl <<
" Picked " << assHits.size() <<
" hits, " << hits.size() <<
" left." << endl;
217 if (assHits.size()+hits.size()<
maxfound)
return;
220 for (hitCont::const_iterator
hit=hits.begin();
hit!=hits.end(); ++
hit) {
222 cout <<
" Trying B: " << **
hit<<
" wire: " << (*hit)->
id() << endl;
225 bool left_ok=
fitWithT0(assHits, chi2l, t0_corrl,0);
229 bool right_ok=
fitWithT0(assHits, chi2r, t0_corrr,0);
233 int nHits=assHits.size()+1;
234 cout <<
" Left: t0_corr = " << t0_corrl <<
"ns chi2/nHits = " << chi2l <<
"/" << nHits <<
" ok: " << left_ok << endl;
235 cout <<
" Right: t0_corr = " << t0_corrr <<
"ns chi2/nHits = " << chi2r <<
"/" << nHits <<
" ok: " << right_ok << endl;
238 if (!left_ok && !right_ok)
continue;
240 foundSomething =
true;
244 for (hitCont::const_iterator tmpHit=
hit+1; tmpHit!=hits.end(); tmpHit++)
245 if (
geometryFilter((*tmpHit)->id(),(*hit)->id())) hitsForFit.push_back(*tmpHit);
247 reverse(hitsForFit.begin(),hitsForFit.end());
250 if (assHits.size()>3 && left_ok && right_ok) {
251 if (chi2l<chi2r-0.1) right_ok=
false;
else
252 if (chi2r<chi2l-0.1) left_ok=
false;
256 addHits(sl,assHits,hitsForFit,result,usedHits);
261 addHits(sl,assHits,hitsForFit,result,usedHits);
266 if (foundSomething)
return;
269 if (assHits.size()<
maxfound)
return;
277 pointsSet.insert(assHits.begin(),assHits.end());
282 for (vector<AssPoint>::const_iterator hiti = assHits.begin()+1; hiti != assHits.end()-1; hiti++)
283 usedHits.push_back(*hiti);
287 if (
debug)
cout <<
"Segment built: " << *seg<< endl;
289 result.push_back(seg);
290 if (
debug)
cout <<
" Result is now " << result.size() << endl;
292 if (
debug)
cout <<
" Exists - skipping" << endl;
304 const int layerLowerCut[4]={0,-1,-2,-2};
305 const int layerUpperCut[4]={0, 2, 2, 3};
316 if (!deltaLayer)
return false;
320 cout <<
"*** WARNING! DT Layer numbers differ by more than 3! for hits: " << endl;
321 cout <<
" " << first << endl;
322 cout <<
" " << second << endl;
327 int deltaWire=first.
wire()-second.
wire();
328 if (second.
layerId().
layer()%2==0) deltaWire=-deltaWire;
329 if ((deltaWire<=layerLowerCut[deltaLayer]) || (deltaWire>=layerUpperCut[deltaLayer]))
return false;
338 typedef vector < pair<double,double> > hitCoord;
339 double a,
b,coordError,
x,
y;
340 double sx=0,sy=0,sxy=0,sxx=0,ssx=0,ssy=0,
s=0,ss=0;
341 int leftHits=0,rightHits=0;
343 if (assHits.size()<3)
return false;
346 coordError=((*(assHits.begin())).first)->localPositionError().xx();
348 for (vector<AssPoint>::const_iterator assHit=assHits.begin(); assHit!=assHits.end(); ++assHit) {
349 if (coordError!=((*(assHits.begin())).first)->localPositionError().xx())
350 cout <<
" DTMeantimerPatternReco: Warning! Hit errors different within one segment!" << endl;
352 x=((*assHit).first)->localPosition((*assHit).second).z();
353 y=((*assHit).first)->localPosition((*assHit).second).x();
374 if (fitdebug &&
debug)
375 cout <<
" DTMeantimerPatternReco::fitWithT0 Left hits: " << leftHits <<
" Right hits: " << rightHits << endl;
377 if (leftHits && rightHits) {
379 double delta = ss*ss*sxx+
s*sx*sx+
s*ssx*ssx-
s*
s*sxx-2*ss*sx*ssx;
383 a=(ssy*
s*ssx+sxy*ss*ss+sy*sx*
s-sy*ss*ssx-ssy*sx*ss-sxy*
s*
s)/delta;
384 b=(ssx*sy*ssx+sxx*ssy*ss+sx*sxy*
s-sxx*sy*
s-ssx*sxy*ss-sx*ssy*ssx)/delta;
385 t0_corr=(ssx*
s*sxy+sxx*ss*sy+sx*sx*ssy-sxx*
s*ssy-sx*ss*sxy-ssx*sx*sy)/delta;
388 double d =
s*sxx - sx*sx;
389 b = (sxx*sy- sx*sxy)/ d;
390 a = (
s*sxy - sx*sy) / d;
396 double chi,chi2_not0;
400 for (vector<AssPoint>::const_iterator assHit=assHits.begin(); assHit!=assHits.end(); ++assHit) {
401 x=((*assHit).first)->localPosition((*assHit).second).z();
402 y=((*assHit).first)->localPosition((*assHit).second).x();
405 chi2_not0+=chi*chi/coordError;
410 chi2+=chi*chi/coordError;
414 if (assHits.size()<4) {
417 if (chi2<200.)
return true;
423 if (
debug && fitdebug)
425 cout <<
" t0_corr = " << t0_corr <<
"ns chi2/nHits = " << chi2 <<
"/" << assHits.size() << endl;
427 if (((chi2/(assHits.size()-2)<
theMaxChi2) && (t0_corr<theMaxT0) && (t0_corr>
theMinT0)) || (assHits.size()==4)) {
428 cout <<
" a = " << a <<
" b = " << b << endl;
429 for (vector<AssPoint>::const_iterator assHit=assHits.begin(); assHit!=assHits.end(); ++assHit) {
430 x=((*assHit).first)->localPosition((*assHit).second).z();
431 y=((*assHit).first)->localPosition((*assHit).second).x();
432 cout <<
" z= " << x <<
" x= " <<
y;
434 else cout <<
" x_corr= " << y-t0_corr*0.00543;
435 cout <<
" seg= " << a*x+b << endl;
440 if ((chi2/(assHits.size()-2)<
theMaxChi2) && (t0_corr<theMaxT0) && (t0_corr>
theMinT0))
return true;
447 double s=0,sx=0,sy=0,
x,
y;
451 if (hits.size()==0)
return;
453 if (hits.size()==1) {
454 b=(*(hits.begin())).second;
456 for (
unsigned int i = 0;
i != hits.size();
i++) {
467 double d = s*sxx - sx*sx;
468 b = (sxx*sy- sx*sxy)/ d;
469 a = (s*sxy - sx*sy) / d;
476 for (vector<DTSegmentCand*>::iterator cand=cands.begin();
477 cand!=cands.end(); ++cand) {
478 if (*(*cand)==*seg)
return false;
479 if (((*cand)->nHits()>=seg->
nHits()) && ((*cand)->chi2ndof()<seg->
chi2ndof()))
480 if ((*cand)->nSharedHitPairs(*seg)>int(seg->
nHits()-2))
return false;
T getParameter(std::string const &) const
std::pair< DTHitPairForFit *, DTEnums::DTCellSide > AssPoint
T getUntrackedParameter(std::string const &, T const &) const
void update(DTRecSegment4D *seg, const bool calcT0=false) const
recompute hits position and refit the segment4D
virtual double chi2ndof() const
the chi2/NDOF of the fit
bool checkDoubleCandidates(std::vector< DTSegmentCand * > &segs, DTSegmentCand *seg)
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
DTCellSide
Which side of the DT cell.
virtual void setChi2(double &chi2)
set chi2
virtual ~DTMeantimerPatternReco()
Destructor.
int layer() const
Return the layer number.
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
DTSegmentCleaner * theCleaner
Geom::Theta< T > theta() const
U second(std::pair< T, U > const &p)
DTSuperLayerId id() const
Return the DetId of this SL.
virtual void setES(const edm::EventSetup &setup)
void setES(const edm::EventSetup &setup)
set the setup
unsigned int theMaxAllowedHits
std::vector< DTHitPairForFit * > initHits(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
int superLayer() const
Return the superlayer number.
bool fit(DTSegmentCand *seg) const
bool fitWithT0(const std::vector< AssPoint > &assHits, double &chi2, double &t0_corr, const bool fitdebug)
int wire() const
Return the wire number.
std::set< AssPoint, AssPointLessZ > AssPointCont
void rawFit(double &a, double &b, const std::vector< std::pair< double, double > > &hits)
virtual unsigned int nHits() const
std::vector< DTSegmentCand * > buildSegments(const DTSuperLayer *sl, const std::vector< DTHitPairForFit * > &hits)
bool geometryFilter(const DTWireId first, const DTWireId second) const
virtual edm::OwnVector< DTSLRecSegment2D > reconstruct(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
this function is called in the producer
std::vector< DTSegmentCand * > clean(std::vector< DTSegmentCand * > inputCands) const
do the cleaning
void addHits(const DTSuperLayer *sl, std::vector< AssPoint > &assHits, const std::vector< DTHitPairForFit * > &hits, std::vector< DTSegmentCand * > &result, std::vector< AssPoint > &usedHits)
DTLayerId layerId() const
Return the corresponding LayerId.
DTMeantimerPatternReco(const edm::ParameterSet &pset)
Constructor.
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
edm::ESHandle< DTGeometry > theDTGeometry
DTSegmentUpdator * theUpdator