CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DuplicateTrackMerger.cc
Go to the documentation of this file.
1 
9 
13 
14 
18 
20 
21 #include "TrackMerger.h"
22 
23 #include <vector>
24 #include <algorithm>
25 #include <string>
26 #include <iostream>
27 #include <atomic>
28 
30 
31 using namespace reco;
32 namespace {
33  class DuplicateTrackMerger final : public edm::stream::EDProducer<> {
34  public:
36  explicit DuplicateTrackMerger(const edm::ParameterSet& iPara);
38  virtual ~DuplicateTrackMerger();
39 
40  using CandidateToDuplicate = std::vector<std::pair<int, int>>;
41 
42  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
43 
44 
45  private:
47  void produce( edm::Event &, const edm::EventSetup &) override;
48 
49  private:
51  const GBRForest* forest_;
52 
53 
55  std::string dbFileName_;
56  bool useForestFromDB_;
57  std::string forestLabel_;
58 
59 
63  double minDeltaR3d2_;
65  double minBDTG_;
67  double minpT2_;
69  double minP_;
71  float maxDCA2_;
73  float maxDPhi_;
75  float maxDLambda_;
77  float maxDdxy_;
79  float maxDdsz_;
81  float maxDQoP_;
82 
84 
86  TrackMerger merger_;
87  };
88  }
89 
99 #include "TFile.h"
100 
101 namespace {
102 
103 
104 void
106 {
108  desc.add<edm::InputTag>("source",edm::InputTag());
109  desc.add<double>("minDeltaR3d",-4.0);
110  desc.add<double>("minBDTG",-0.1);
111  desc.add<double>("minpT",0.2);
112  desc.add<double>("minP",0.4);
113  desc.add<double>("maxDCA",30.0);
114  desc.add<double>("maxDPhi",0.30);
115  desc.add<double>("maxDLambda",0.30);
116  desc.add<double>("maxDdsz",10.0);
117  desc.add<double>("maxDdxy",10.0);
118  desc.add<double>("maxDQoP",0.25);
119  desc.add<std::string>("forestLabel","MVADuplicate");
120  desc.add<std::string>("GBRForestFileName","");
121  desc.add<bool>("useInnermostState",true);
122  desc.add<std::string>("ttrhBuilderName","WithAngleAndTemplate");
123  descriptions.add("DuplicateTrackMerger", desc);
124 }
125 
126 
127 DuplicateTrackMerger::DuplicateTrackMerger(const edm::ParameterSet& iPara) : forest_(nullptr), merger_(iPara)
128 {
129 
130  trackSource_ = consumes<reco::TrackCollection>(iPara.getParameter<edm::InputTag>("source"));
131  minDeltaR3d2_ = iPara.getParameter<double>("minDeltaR3d"); minDeltaR3d2_*=std::abs(minDeltaR3d2_);
132  minBDTG_ = iPara.getParameter<double>("minBDTG");
133  minpT2_ = iPara.getParameter<double>("minpT"); minpT2_ *= minpT2_;
134  minP_ = iPara.getParameter<double>("minP");
135  maxDCA2_ = iPara.getParameter<double>("maxDCA"); maxDCA2_*=maxDCA2_;
136  maxDPhi_ = iPara.getParameter<double>("maxDPhi");
137  maxDLambda_ = iPara.getParameter<double>("maxDLambda");
138  maxDdsz_ = iPara.getParameter<double>("maxDdsz");
139  maxDdxy_ = iPara.getParameter<double>("maxDdxy");
140  maxDQoP_ = iPara.getParameter<double>("maxDQoP");
141 
142  produces<std::vector<TrackCandidate> >("candidates");
143  produces<CandidateToDuplicate>("candidateMap");
144 
145  forestLabel_ = iPara.getParameter<std::string>("forestLabel");
146 
147  dbFileName_ = iPara.getParameter<std::string>("GBRForestFileName");
148  useForestFromDB_ = dbFileName_.empty();
149 
150  /*
151  tmvaReader_ = new TMVA::Reader("!Color:Silent");
152  tmvaReader_->AddVariable("ddsz",&tmva_ddsz_);
153  tmvaReader_->AddVariable("ddxy",&tmva_ddxy_);
154  tmvaReader_->AddVariable("dphi",&tmva_dphi_);
155  tmvaReader_->AddVariable("dlambda",&tmva_dlambda_);
156  tmvaReader_->AddVariable("dqoverp",&tmva_dqoverp_);
157  tmvaReader_->AddVariable("delta3d_r",&tmva_d3dr_);
158  tmvaReader_->AddVariable("delta3d_z",&tmva_d3dz_);
159  tmvaReader_->AddVariable("outer_nMissingInner",&tmva_outer_nMissingInner_);
160  tmvaReader_->AddVariable("inner_nMissingOuter",&tmva_inner_nMissingOuter_);
161  tmvaReader_->BookMVA("BDTG",mvaFilePath);
162  */
163 
164 
165 }
166 
167 DuplicateTrackMerger::~DuplicateTrackMerger()
168 {
169 
170  if(!useForestFromDB_) delete forest_;
171 
172 }
173 
174 
175 #ifdef VI_STAT
176  struct Stat {
177  Stat() : maxCos(1.1), nCand(0),nLoop0(0) {}
178  ~Stat() {
179  std::cout << "Stats " << nCand << ' ' << nLoop0 << ' ' << maxCos << std::endl;
180  }
181  std::atomic<float> maxCos;
182  std::atomic<int> nCand, nLoop0;
183  };
184  Stat stat;
185 #endif
186 
187 template<typename T>
188 void update_maximum(std::atomic<T>& maximum_value, T const& value) noexcept
189 {
190  T prev_value = maximum_value;
191  while(prev_value < value &&
192  !maximum_value.compare_exchange_weak(prev_value, value))
193  ;
194 }
195 
196  template<typename T>
197 void update_minimum(std::atomic<T>& minimum_value, T const& value) noexcept
198 {
199  T prev_value = minimum_value;
200  while(prev_value > value &&
201  !minimum_value.compare_exchange_weak(prev_value, value))
202  ;
203 }
204 
205 
206 void DuplicateTrackMerger::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
207 {
208 
209  merger_.init(iSetup);
210 
211  if(!forest_){
212  if(useForestFromDB_){
213  edm::ESHandle<GBRForest> forestHandle;
214  iSetup.get<GBRWrapperRcd>().get(forestLabel_,forestHandle);
215  forest_ = forestHandle.product();
216  }else{
217  TFile gbrfile(dbFileName_.c_str());
218  forest_ = dynamic_cast<const GBRForest*>(gbrfile.Get(forestLabel_.c_str()));
219  }
220  }
221 
222  //edm::Handle<edm::View<reco::Track> >handle;
224  iEvent.getByToken(trackSource_,handle);
225  auto const & tracks = *handle;
226 
227  iSetup.get<IdealMagneticFieldRecord>().get(magfield_);
229  TSCPBuilderNoMaterial tscpBuilder;
230  auto out_duplicateCandidates = std::make_unique<std::vector<TrackCandidate>>();
231 
232  auto out_candidateMap = std::make_unique<CandidateToDuplicate>();
233  LogDebug("DuplicateTrackMerger") << "Number of tracks to be checked for merging: " << tracks.size();
234 
235  for(int i = 0; i < (int)tracks.size(); i++){
236  const reco::Track *rt1 = &tracks[i];
237  if(rt1->innerMomentum().perp2() < minpT2_)continue;
238  // if(rt1->innerMomentum().R() < minP_)continue;
239  for(int j = i+1; j < (int)tracks.size();j++){
240  const reco::Track *rt2 = &tracks[j];
241  if(rt1->charge() != rt2->charge())continue;
242  auto cosT = (*rt1).momentum().unit().Dot((*rt2).momentum().unit());
243  if (cosT<0.) continue;
244  if(rt2->innerMomentum().perp2() < minpT2_)continue;
245  // if(rt2->innerMomentum().R() < minP_)continue;
246  const reco::Track* t1,*t2;
247  if(rt1->outerPosition().perp2() < rt2->outerPosition().perp2()){
248  t1 = rt1;
249  t2 = rt2;
250  }else{
251  t1 = rt2;
252  t2 = rt1;
253  }
254  auto deltaR3d2 = (t1->outerPosition() - t2->innerPosition()).mag2();
255 
256  if(t1->outerPosition().perp2() > t2->innerPosition().perp2()) deltaR3d2 *= -1.0;
257  if(deltaR3d2 < minDeltaR3d2_)continue;
258 
259  FreeTrajectoryState fts1 = trajectoryStateTransform::outerFreeState(*t1, &*magfield_,false);
260  FreeTrajectoryState fts2 = trajectoryStateTransform::innerFreeState(*t2, &*magfield_,false);
261  GlobalPoint avgPoint((t1->outerPosition().x()+t2->innerPosition().x())*0.5,(t1->outerPosition().y()+t2->innerPosition().y())*0.5,(t1->outerPosition().z()+t2->innerPosition().z())*0.5);
262  TrajectoryStateClosestToPoint TSCP1 = tscpBuilder(fts1, avgPoint);
263  TrajectoryStateClosestToPoint TSCP2 = tscpBuilder(fts2, avgPoint);
264  if(!TSCP1.isValid())continue;
265  if(!TSCP2.isValid())continue;
266 
267  const FreeTrajectoryState ftsn1 = TSCP1.theState();
268  const FreeTrajectoryState ftsn2 = TSCP2.theState();
269 
270  if ( (ftsn2.position()-ftsn1.position()).mag2() > maxDCA2_ ) continue;
271 
272  auto qoverp1 = ftsn1.signedInverseMomentum();
273  auto qoverp2 = ftsn2.signedInverseMomentum();
274  float tmva_dqoverp_ = qoverp1-qoverp2;
275  if ( std::abs(tmva_dqoverp_) > maxDQoP_ ) continue;
276 
277 
278  //auto pp = [&](TrajectoryStateClosestToPoint const & ts) { std::cout << ' ' << ts.perigeeParameters().vector()[0] << '/' << ts.perigeeError().transverseCurvatureError();};
279  //if(qoverp1*qoverp2 <0) { std::cout << "charge different " << qoverp1 <<',' << qoverp2; pp(TSCP1); pp(TSCP2); std::cout << std::endl;}
280 
281  auto lambda1 = M_PI/2 - ftsn1.momentum().theta();
282  auto lambda2 = M_PI/2 - ftsn2.momentum().theta();
283  float tmva_dlambda_ = lambda1-lambda2;
284  if ( std::abs(tmva_dlambda_) > maxDLambda_ ) continue;
285 
286  auto phi1 = ftsn1.momentum().phi();
287  auto phi2 = ftsn2.momentum().phi();
288  float tmva_dphi_ = phi1-phi2;
289  if(std::abs(tmva_dphi_) > float(M_PI)) tmva_dphi_ = 2.f*float(M_PI) - std::abs(tmva_dphi_);
290  if (std::abs(tmva_dphi_) > maxDPhi_ ) continue;
291 
292  auto dxy1 = (-ftsn1.position().x() * ftsn1.momentum().y() + ftsn1.position().y() * ftsn1.momentum().x())/TSCP1.pt();
293  auto dxy2 = (-ftsn2.position().x() * ftsn2.momentum().y() + ftsn2.position().y() * ftsn2.momentum().x())/TSCP2.pt();
294  float tmva_ddxy_ = dxy1-dxy2;
295  if ( std::abs(tmva_ddxy_) > maxDdxy_ ) continue;
296 
297  auto dsz1 = ftsn1.position().z() * TSCP1.pt() / TSCP1.momentum().mag()
298  - (ftsn1.position().x() * ftsn1.momentum().y() + ftsn1.position().y() * ftsn1.momentum().x())/TSCP1.pt() * ftsn1.momentum().z()/ftsn1.momentum().mag();
299  auto dsz2 = ftsn2.position().z() * TSCP2.pt() / TSCP2.momentum().mag()
300  - (ftsn2.position().x() * ftsn2.momentum().y() + ftsn2.position().y() * ftsn2.momentum().x())/TSCP2.pt() * ftsn2.momentum().z()/ftsn2.momentum().mag();
301  float tmva_ddsz_ = dsz1-dsz2;
302  if ( std::abs(tmva_ddsz_) > maxDdsz_ ) continue;
303 
304  float tmva_d3dr_ = avgPoint.perp();
305  float tmva_d3dz_ = avgPoint.z();
306  float tmva_outer_nMissingInner_ = t2->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
307  float tmva_inner_nMissingOuter_ = t1->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_OUTER_HITS);
308 
309  float gbrVals_[9];
310  gbrVals_[0] = tmva_ddsz_;
311  gbrVals_[1] = tmva_ddxy_;
312  gbrVals_[2] = tmva_dphi_;
313  gbrVals_[3] = tmva_dlambda_;
314  gbrVals_[4] = tmva_dqoverp_;
315  gbrVals_[5] = tmva_d3dr_;
316  gbrVals_[6] = tmva_d3dz_;
317  gbrVals_[7] = tmva_outer_nMissingInner_;
318  gbrVals_[8] = tmva_inner_nMissingOuter_;
319 
320 
321  auto mvaBDTG = forest_->GetClassifier(gbrVals_);
322  if(mvaBDTG < minBDTG_)continue;
323 
324  // std::cout << "to merge " << mvaBDTG << ' ' << std::copysign(std::sqrt(std::abs(deltaR3d2)),deltaR3d2) << ' ' << tmva_dphi_ << ' ' << TSCP1.pt() <<'/'<<TSCP2.pt() << std::endl;
325 
326  out_duplicateCandidates->push_back(merger_.merge(*t1,*t2));
327  out_candidateMap->emplace_back(i,j);
328 
329 #ifdef VI_STAT
330  ++stat.nCand;
331  // auto cosT = float((*t1).momentum().unit().Dot((*t2).momentum().unit()));
332  if (cosT>0) update_minimum(stat.maxCos,float(cosT));
333  else ++stat.nLoop0;
334 #endif
335 
336  }
337  }
338  iEvent.put(std::move(out_duplicateCandidates),"candidates");
339  iEvent.put(std::move(out_candidateMap),"candidateMap");
340 
341 
342 }
343 
344 
345 
346 }
347 
350 
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
const FreeTrajectoryState & theState() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
static void fillDescriptions(ConfigurationDescriptions &descriptions)
T y() const
Definition: PV3DBase.h:63
#define noexcept
#define nullptr
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:65
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:55
int iEvent
Definition: GenABIO.cc:230
T mag() const
Definition: PV3DBase.h:67
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
T z() const
Definition: PV3DBase.h:64
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
def move
Definition: eostools.py:510
tuple handle
Definition: patZpeak.py:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
ParameterDescriptionBase * add(U const &iLabel, T const &value)
GlobalVector momentum() const
#define M_PI
GlobalPoint position() const
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:445
tuple tracks
Definition: testEve_cfg.py:39
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
int numberOfLostHits(HitCategory category) const
Definition: HitPattern.h:902
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:60
tuple cout
Definition: gather_cfg.py:145
int charge() const
track electric charge
Definition: TrackBase.h:562
FreeTrajectoryState innerFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
long double T
T x() const
Definition: PV3DBase.h:62
FreeTrajectoryState outerFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
double signedInverseMomentum() const