Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:56:27

0001 /** \file LaserAlignmentProducer.cc
0002  *  Producer to be used for the Simulation of the Laser Alignment System
0003  *  an empty MCHepEvent will be generated (needed by OscarProducer). The actual
0004  * simulation of the laser beams is done in the SimWatcher attached to
0005  * OscarProducer
0006  *
0007  *  $Date: 2011/09/16 06:23:27 $
0008  *  $Revision: 1.6 $
0009  *  \author Maarten Thomas
0010  */
0011 // system include files
0012 #include "FWCore/Framework/interface/Event.h"
0013 
0014 // user include files
0015 #include "Alignment/LaserAlignmentSimulation/plugins/LaserAlignmentProducer.h"
0016 #include "FWCore/Framework/interface/MakerMacros.h"
0017 
0018 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0019 
0020 //
0021 // constructors and destructor
0022 //
0023 LaserAlignmentProducer::LaserAlignmentProducer(const edm::ParameterSet &) : EDProducer(), theEvent(nullptr) {
0024   // register your products
0025   produces<edm::HepMCProduct>("unsmeared");
0026 
0027   // now do what ever other initialization is needed
0028 }
0029 
0030 LaserAlignmentProducer::~LaserAlignmentProducer() {
0031   // no need to cleanup theEvent since it's done in HepMCProduct
0032 }
0033 
0034 // ------------ method called to produce the event  ------------
0035 void LaserAlignmentProducer::produce(edm::Event &iEvent, const edm::EventSetup &) {
0036   // create the event
0037   theEvent = new HepMC::GenEvent();
0038 
0039   // create a primary vertex
0040   HepMC::GenVertex *theVtx = new HepMC::GenVertex(HepMC::FourVector(0., 0., 0.));
0041 
0042   // add a particle to the vertex; this is needed to avoid crashes in
0043   // OscarProducer. Use a electron neutrino, with zero energy and mass
0044   HepMC::GenParticle *theParticle = new HepMC::GenParticle(HepMC::FourVector(0., 0., 0., 0.), 12, 1);
0045 
0046   theVtx->add_particle_out(theParticle);
0047 
0048   // add the vertex to the event
0049   theEvent->add_vertex(theVtx);
0050 
0051   // set the event number
0052   theEvent->set_event_number(iEvent.id().event());
0053   // set the signal process id
0054   theEvent->set_signal_process_id(20);
0055 
0056   // create an empty output collection
0057   auto theOutput = std::make_unique<edm::HepMCProduct>();
0058   theOutput->addHepMCData(theEvent);
0059 
0060   // put the output to the event
0061   iEvent.put(std::move(theOutput));
0062 }
0063 
0064 // define this as a plug-in
0065 
0066 DEFINE_FWK_MODULE(LaserAlignmentProducer);