CMS 3D CMS Logo

MixEvtVtxGenerator.cc
Go to the documentation of this file.
1 #ifndef HI_MixEvtVtxGenerator_H
2 #define HI_MixEvtVtxGenerator_H
3 /*
4 */
7 
15 
19 
22 
23 #include "TMatrixD.h"
24 #include <iostream>
25 
26 using namespace edm;
27 using namespace std;
28 
29 
30 namespace HepMC {
31  class FourVector ;
32 }
33 
34 
36 {
37  public:
38 
39  // ctor & dtor
40  explicit MixEvtVtxGenerator( const edm::ParameterSet& );
41  virtual ~MixEvtVtxGenerator();
42 
43  virtual void produce( edm::Event&, const edm::EventSetup& ) override;
44 
45  virtual HepMC::FourVector* getVertex(edm::Event&);
46  virtual HepMC::FourVector* getRecVertex(edm::Event&);
47 
48  protected:
49 
50  HepMC::FourVector* fVertex ;
51  TMatrixD *boost_;
52 
53  private :
54 
58 
60  std::vector<double> vtxOffset;
61  bool useCF_;
62 
63 };
64 
66  : fVertex(0), boost_(0),
67  useRecVertex(pset.exists("useRecVertex")?pset.getParameter<bool>("useRecVertex"):false)
68 
69 {
70  produces<edm::HepMCProduct>();
71  vtxOffset.resize(3);
72  if(pset.exists("vtxOffset")) vtxOffset=pset.getParameter< std::vector<double> >("vtxOffset");
73 
74  if(useRecVertex) useCF_ = 0;
75  else{
76  useCF_ = pset.getUntrackedParameter<bool>("useCF",false);
77  cfLabel = consumes<CrossingFrame<HepMCProduct> >(pset.getParameter<edm::InputTag>("mixLabel"));
78  }
79  signalLabel = consumes<HepMCProduct>(pset.getParameter<edm::InputTag>("signalLabel"));
80 
81 }
82 
84 {
85  delete fVertex ;
86  if (boost_ != 0 ) delete boost_;
87  // no need since now it's done in HepMCProduct
88  // delete fEvt ;
89 }
90 
91 HepMC::FourVector* MixEvtVtxGenerator::getVertex( Event& evt){
92 
93  HepMC::GenVertex* genvtx = 0;
94  const HepMC::GenEvent* inev = 0;
95 
96  if(useCF_){
98  evt.getByToken(cfLabel,cf);
100  if(mix.size() < 2){
101  throw cms::Exception("MatchVtx")
102  <<"Mixing has "<<mix.size()<<" sub-events, should have been at least 2"
103  <<endl;
104  }
105  const HepMCProduct& bkg = mix.getObject(1);
106  if(!(bkg.isVtxGenApplied())){
107  throw cms::Exception("MatchVtx")<<"Input background does not have smeared vertex!"<<endl;
108  }else{
109  inev = bkg.GetEvent();
110  }
111  }else{
113  evt.getByToken(signalLabel,input);
114  inev = input->GetEvent();
115  }
116 
117  genvtx = inev->signal_process_vertex();
118  if(!genvtx){
119  HepMC::GenEvent::particle_const_iterator pt=inev->particles_begin();
120  HepMC::GenEvent::particle_const_iterator ptend=inev->particles_end();
121  while(!genvtx || ( genvtx->particles_in_size() == 1 && pt != ptend ) ){
122  if(pt == ptend) cout<<"End reached, No Gen Vertex!"<<endl;
123  genvtx = (*pt)->production_vertex();
124  ++pt;
125  }
126  }
127  double aX,aY,aZ,aT;
128 
129  aX = genvtx->position().x();
130  aY = genvtx->position().y();
131  aZ = genvtx->position().z();
132  aT = genvtx->position().t();
133 
134  if(!fVertex){
135  fVertex = new HepMC::FourVector();
136  }
137  LogInfo("MatchVtx")<<" setting vertex "<<" aX "<<aX<<" aY "<<aY<<" aZ "<<aZ<<" aT "<<aT<<endl;
138  fVertex->set(aX,aY,aZ,aT);
139 
140 
141  return fVertex;
142 
143 }
144 
145 
146 HepMC::FourVector* MixEvtVtxGenerator::getRecVertex( Event& evt){
147 
149  evt.getByToken(hiLabel,input);
150 
151  double aX,aY,aZ;
152 
153  aX = input->begin()->position().x() + vtxOffset[0];
154  aY = input->begin()->position().y() + vtxOffset[1];
155  aZ = input->begin()->position().z() + vtxOffset[2];
156 
157  if(!fVertex) fVertex = new HepMC::FourVector();
158  fVertex->set(10.0*aX,10.0*aY,10.0*aZ,0.0); // HepMC positions in mm (RECO in cm)
159 
160  return fVertex;
161 
162 }
163 
165 {
166  Handle<HepMCProduct> HepUnsmearedMCEvt ;
167 
168  evt.getByToken( signalLabel, HepUnsmearedMCEvt ) ;
169 
170  // generate new vertex & apply the shift
171  //
172  if(HepUnsmearedMCEvt->isVtxGenApplied()) throw cms::Exception("MatchVtx")
173  <<"Signal HepMCProduct is not compatible for embedding - it's vertex is already smeared."
174  <<std::endl;
175  // Copy the HepMC::GenEvent
176  HepMC::GenEvent* genevt = new HepMC::GenEvent(*HepUnsmearedMCEvt->GetEvent());
177  std::unique_ptr<edm::HepMCProduct> HepMCEvt(new edm::HepMCProduct(genevt));
178  // generate new vertex & apply the shift
179  //
180  HepMCEvt->applyVtxGen( useRecVertex ? getRecVertex(evt) : getVertex(evt) ) ;
181 
182  // HepMCEvt->boostToLab( GetInvLorentzBoost(), "vertex" );
183  // HepMCEvt->boostToLab( GetInvLorentzBoost(), "momentum" );
184 
185  evt.put(std::move(HepMCEvt)) ;
186 
187  return ;
188 
189 }
190 
192 
193 #endif
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
virtual HepMC::FourVector * getRecVertex(edm::Event &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
std::vector< double > vtxOffset
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
MixEvtVtxGenerator(const edm::ParameterSet &)
virtual void produce(edm::Event &, const edm::EventSetup &) override
virtual HepMC::FourVector * getVertex(edm::Event &)
static std::string const input
Definition: EdmProvDump.cc:44
edm::EDGetTokenT< CrossingFrame< HepMCProduct > > cfLabel
edm::EDGetTokenT< reco::VertexCollection > hiLabel
edm::EDGetTokenT< HepMCProduct > signalLabel
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
T const * product() const
Definition: Handle.h:81
HepMC::FourVector * fVertex
bool isVtxGenApplied() const
Definition: HepMCProduct.h:40
return(e1-e2)*(e1-e2)+dp *dp
HLT enums.
def move(src, dest)
Definition: eostools.py:510