CMS 3D CMS Logo

TwoBodyDecayConstraintProducer.cc
Go to the documentation of this file.
1 
10 #include <memory>
11 
18 
20 
23 
25 
30 
31 // // debug
32 // #include <map>
33 // #include "TH1F.h"
34 // #include "TFile.h"
35 // #include "TLorentzVector.h"
36 // #include "Alignment/TwoBodyDecay/interface/TwoBodyDecayModel.h"
37 
39 public:
41  ~TwoBodyDecayConstraintProducer() override = default;
42 
43 private:
44  void produce(edm::StreamID streamid, edm::Event&, const edm::EventSetup&) const override;
45 
46  std::pair<bool, TrajectoryStateOnSurface> innermostState(const reco::TransientTrack& ttrack) const;
47  bool match(const TrajectoryStateOnSurface& newTsos, const TrajectoryStateOnSurface& oldTsos) const;
48 
51 
53 
54  const double primaryMass_;
55  const double primaryWidth_;
56  const double secondaryMass_;
57 
58  const double sigmaPositionCutValue_;
59  const double chi2CutValue_;
60  const double errorRescaleValue_;
61 
64 
65  // // debug
66  // std::map<std::string, TH1F*> histos_;
67 };
68 
70  : srcTag_(iConfig.getParameter<edm::InputTag>("src")),
71  bsSrcTag_(iConfig.getParameter<edm::InputTag>("beamSpot")),
72  tbdFitter_(iConfig),
73  primaryMass_(iConfig.getParameter<double>("primaryMass")),
74  primaryWidth_(iConfig.getParameter<double>("primaryWidth")),
75  secondaryMass_(iConfig.getParameter<double>("secondaryMass")),
76  sigmaPositionCutValue_(iConfig.getParameter<double>("sigmaPositionCut")),
77  chi2CutValue_(iConfig.getParameter<double>("chi2Cut")),
78  errorRescaleValue_(iConfig.getParameter<double>("rescaleError")) {
79  trackCollToken_ = consumes<reco::TrackCollection>(edm::InputTag(srcTag_));
80  bsToken_ = consumes<reco::BeamSpot>(edm::InputTag(bsSrcTag_));
81 
82  produces<std::vector<TrackParamConstraint> >();
83  produces<TrackParamConstraintAssociationCollection>();
84 
85  // //debug
86  // histos_["deltaEta1"] = new TH1F( "deltaEta1", "deltaEta1", 200, -1., 1. );
87  // histos_["deltaP1"] = new TH1F( "deltaP1", "deltaP1", 200, -50., 50. );
88 
89  // histos_["deltaEta2"] = new TH1F( "deltaEta2", "deltaEta2", 200, -1., 1. );
90  // histos_["deltaP2"] = new TH1F( "deltaP2", "deltaP2", 200, -50., 50. );
91 }
92 
95  const edm::EventSetup& iSetup) const {
96  using namespace edm;
97 
99  iEvent.getByToken(trackCollToken_, trackColl);
100 
102  iEvent.getByToken(bsToken_, beamSpot);
103 
104  ESHandle<MagneticField> magField;
105  iSetup.get<IdealMagneticFieldRecord>().get(magField);
106 
108  iEvent.getRefBeforePut<std::vector<TrackParamConstraint> >();
109  std::unique_ptr<std::vector<TrackParamConstraint> > pairs(new std::vector<TrackParamConstraint>);
110  std::unique_ptr<TrackParamConstraintAssociationCollection> output(
111  new TrackParamConstraintAssociationCollection(trackColl, rPairs));
112 
113  if (trackColl->size() == 2) {
116 
118  std::vector<reco::TransientTrack> ttracks(2);
119  ttracks[0] = reco::TransientTrack(reco::TrackRef(trackColl, 0), magField.product());
120  ttracks[0].setES(iSetup);
121  ttracks[1] = reco::TransientTrack(reco::TrackRef(trackColl, 1), magField.product());
122  ttracks[1].setES(iSetup);
123 
125  TwoBodyDecay tbd = tbdFitter_.estimate(ttracks, vm);
126 
127  if (!tbd.isValid() or (tbd.chi2() > chi2CutValue_))
128  return;
129 
131  std::pair<bool, TrajectoryStateOnSurface> oldInnermostState1 = innermostState(ttracks[0]);
132  std::pair<bool, TrajectoryStateOnSurface> oldInnermostState2 = innermostState(ttracks[1]);
133  if (!oldInnermostState1.second.isValid() || !oldInnermostState2.second.isValid())
134  return;
135 
137  TwoBodyDecayTrajectoryState::TsosContainer trackTsos(oldInnermostState1.second, oldInnermostState2.second);
138  TwoBodyDecayTrajectoryState tbdTrajState(trackTsos, tbd, secondaryMass_, magField.product(), true);
139  if (!tbdTrajState.isValid())
140  return;
141 
143  bool match1 = match(tbdTrajState.trajectoryStates(true).first, oldInnermostState1.second);
144  bool match2 = match(tbdTrajState.trajectoryStates(true).second, oldInnermostState2.second);
145  if (!match1 || !match2)
146  return;
147 
148  // re-scale error of constraintTsos
149  tbdTrajState.rescaleError(errorRescaleValue_);
150 
151  // produce constraint for first track
152  pairs->push_back(tbdTrajState.trajectoryStates(true).first);
153  output->insert(reco::TrackRef(trackColl, 0), edm::Ref<std::vector<TrackParamConstraint> >(rPairs, 0));
154 
155  // produce constraint for second track
156  pairs->push_back(tbdTrajState.trajectoryStates(true).second);
157  output->insert(reco::TrackRef(trackColl, 1), edm::Ref<std::vector<TrackParamConstraint> >(rPairs, 1));
158 
159  // // debug
160  // if ( tbd.isValid() ) {
161  // TwoBodyDecayModel model;
162  // const std::pair< AlgebraicVector, AlgebraicVector > fitMomenta = model.cartesianSecondaryMomenta( tbd );
163 
164  // TLorentzVector recoMomentum1( ttracks[0].track().px(), ttracks[0].track().py(), ttracks[0].track().pz(),
165  // sqrt((ttracks[0].track().p()*ttracks[0].track().p())+0.105658*0.105658) );
166  // TLorentzVector fitMomentum1( fitMomenta.first[0], fitMomenta.first[1], fitMomenta.first[2],
167  // sqrt( fitMomenta.first.normsq()+0.105658*0.105658) );
168  // histos_["deltaP1"]->Fill( recoMomentum1.P() - fitMomentum1.P() );
169  // histos_["deltaEta1"]->Fill( recoMomentum1.Eta() - fitMomentum1.Eta() );
170 
171  // TLorentzVector recoMomentum2( ttracks[1].track().px(), ttracks[1].track().py(), ttracks[1].track().pz(),
172  // sqrt((ttracks[1].track().p()*ttracks[1].track().p())+0.105658*0.105658) );
173  // TLorentzVector fitMomentum2( fitMomenta.second[0], fitMomenta.second[1], fitMomenta.second[2],
174  // sqrt( fitMomenta.second.normsq()+0.105658*0.105658) );
175  // histos_["deltaP2"]->Fill( recoMomentum2.P() - fitMomentum2.P() );
176  // histos_["deltaEta2"]->Fill( recoMomentum2.Eta() - fitMomentum2.Eta() );
177  // }
178  }
179 
180  iEvent.put(std::move(pairs));
181  iEvent.put(std::move(output));
182 }
183 
184 std::pair<bool, TrajectoryStateOnSurface> TwoBodyDecayConstraintProducer::innermostState(
185  const reco::TransientTrack& ttrack) const {
186  double outerR = ttrack.outermostMeasurementState().globalPosition().perp();
187  double innerR = ttrack.innermostMeasurementState().globalPosition().perp();
188  bool insideOut = (outerR > innerR);
189  return std::make_pair(insideOut, insideOut ? ttrack.innermostMeasurementState() : ttrack.outermostMeasurementState());
190 }
191 
193  const TrajectoryStateOnSurface& oldTsos) const {
194  LocalPoint lp1 = newTsos.localPosition();
195  LocalPoint lp2 = oldTsos.localPosition();
196 
197  double deltaX = lp1.x() - lp2.x();
198  double deltaY = lp1.y() - lp2.y();
199 
200  return ((fabs(deltaX) < sigmaPositionCutValue_) && (fabs(deltaY) < sigmaPositionCutValue_));
201 }
202 
203 //define this as a plug-in
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
T perp() const
Definition: PV3DBase.h:69
~TwoBodyDecayConstraintProducer() override=default
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
void produce(edm::StreamID streamid, edm::Event &, const edm::EventSetup &) const override
T y() const
Definition: PV3DBase.h:60
GlobalPoint globalPosition() const
TrajectoryStateOnSurface innermostMeasurementState() const
virtual const TwoBodyDecay estimate(const std::vector< reco::TransientTrack > &tracks, const TwoBodyDecayVirtualMeasurement &vm) const
std::pair< bool, TrajectoryStateOnSurface > innermostState(const reco::TransientTrack &ttrack) const
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::AssociationMap< edm::OneToOne< reco::TrackCollection, std::vector< TrackParamConstraint > > > TrackParamConstraintAssociationCollection
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
TwoBodyDecayConstraintProducer(const edm::ParameterSet &)
edm::EDGetTokenT< reco::BeamSpot > bsToken_
TrajectoryStateOnSurface outermostMeasurementState() const
RefProd< PROD > getRefBeforePut()
Definition: Event.h:156
void setES(const edm::EventSetup &es)
bool isValid(void) const
Definition: TwoBodyDecay.h:46
T const * product() const
Definition: Handle.h:69
std::pair< TrajectoryStateOnSurface, TrajectoryStateOnSurface > TsosContainer
double chi2(void) const
Definition: TwoBodyDecay.h:44
edm::EDGetTokenT< reco::TrackCollection > trackCollToken_
HLT enums.
T get() const
Definition: EventSetup.h:73
T x() const
Definition: PV3DBase.h:59
T const * product() const
Definition: ESHandle.h:86
bool match(const TrajectoryStateOnSurface &newTsos, const TrajectoryStateOnSurface &oldTsos) const
def move(src, dest)
Definition: eostools.py:511