CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/FastSimulation/EventProducer/src/FamosProducer.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/Event.h"
00002 #include "FWCore/Framework/interface/MakerMacros.h"
00003 #include "FWCore/Framework/interface/EventSetup.h"
00004 
00005 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00006 #include "SimDataFormats/CaloHit/interface/PCaloHitContainer.h"
00007 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00008 
00009 #include "FastSimDataFormats/NuclearInteractions/interface/FSimVertexTypeFwd.h"
00010 
00011 #include "DataFormats/Common/interface/Handle.h"
00012 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00013 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00014 #include "DataFormats/Provenance/interface/EventID.h"
00015 
00016 #include "FastSimulation/EventProducer/interface/FamosProducer.h"
00017 #include "FastSimulation/EventProducer/interface/FamosManager.h"
00018 #include "FastSimulation/Event/interface/FSimEvent.h"
00019 #include "FastSimulation/Event/interface/KineParticleFilter.h"
00020 #include "FastSimulation/Event/interface/PrimaryVertexGenerator.h"
00021 #include "FastSimulation/Calorimetry/interface/CalorimetryManager.h"
00022 #include "FastSimulation/TrajectoryManager/interface/TrajectoryManager.h"
00023 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
00024 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00025 
00026 #include "HepMC/GenEvent.h"
00027 #include "HepMC/GenVertex.h"
00028 
00029 #include <iostream>
00030 #include <memory>
00031 #include <vector>
00032 
00033 FamosProducer::FamosProducer(edm::ParameterSet const & p)      
00034 {    
00035 
00036     produces<edm::SimTrackContainer>();
00037     produces<edm::SimVertexContainer>();
00038     produces<FSimVertexTypeCollection>("VertexTypes");
00039     produces<edm::PSimHitContainer>("TrackerHits");
00040     produces<edm::PCaloHitContainer>("EcalHitsEB");
00041     produces<edm::PCaloHitContainer>("EcalHitsEE");
00042     produces<edm::PCaloHitContainer>("EcalHitsES");
00043     produces<edm::PCaloHitContainer>("HcalHits");
00044     // Temporary facility to allow for the crossing frame to work...
00045     simulateMuons = p.getParameter<bool>("SimulateMuons");
00046     if ( simulateMuons ) produces<edm::SimTrackContainer>("MuonSimTracks");
00047 
00048     // The generator input label
00049     theSourceLabel = p.getParameter<edm::InputTag>("SourceLabel");
00050     theGenParticleLabel = p.getParameter<edm::InputTag>("GenParticleLabel");
00051     theBeamSpotLabel = p.getParameter<edm::InputTag>("BeamSpotLabel");
00052 
00053     famosManager_ = new FamosManager(p);
00054   
00055 }
00056 
00057 FamosProducer::~FamosProducer() 
00058 { if ( famosManager_ ) delete famosManager_; }
00059 
00060 void
00061 FamosProducer::beginRun(edm::Run & run, const edm::EventSetup & es) {
00062   famosManager_->setupGeometryAndField(run,es);
00063 }
00064  
00065 void FamosProducer::endJob()
00066 { 
00067 }
00068  
00069 void FamosProducer::produce(edm::Event & iEvent, const edm::EventSetup & es)
00070 {
00071    using namespace edm;
00072 
00073   //  // The beam spot position
00074    edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00075    iEvent.getByLabel(theBeamSpotLabel,recoBeamSpotHandle); 
00076    math::XYZPoint BSPosition_ = recoBeamSpotHandle->position();
00077       
00078    const HepMC::GenEvent* myGenEvent = 0;
00079    FSimEvent* fevt = famosManager_->simEvent();
00080    //   fevt->setBeamSpot(BSPosition_);     
00081    
00082    // Get the generated event(s) from the edm::Event
00083    // 1. Check if a HepMCProduct exists
00084    //    a. Take the VtxSmeared if it exists
00085    //    b. Take the source  otherwise
00086    // 2. Otherwise go for the CandidateCollection
00087    
00088    Handle<HepMCProduct> theHepMCProduct;
00089    
00090    PrimaryVertexGenerator* theVertexGenerator = fevt->thePrimaryVertexGenerator();
00091    
00092    
00093    const reco::GenParticleCollection* myGenParticlesXF = 0;
00094    const reco::GenParticleCollection* myGenParticles = 0;
00095    const HepMC::GenEvent* thePUEvents = 0;
00096 
00097    Handle<CrossingFrame<HepMCProduct> > theHepMCProductCrossingFrame;
00098    bool isPileUpXF = iEvent.getByLabel("mixGenPU","generator",theHepMCProductCrossingFrame);
00099 
00100    if (isPileUpXF){// take the GenParticle from crossingframe event collection, if it exists 
00101      Handle<reco::GenParticleCollection> genEvtXF;
00102      bool genPartXF = iEvent.getByLabel("genParticlesFromMixingModule",genEvtXF);
00103      if(genPartXF) myGenParticlesXF = &(*genEvtXF);
00104    }
00105    else{// otherwise, use the old famos PU     
00106      // Get the generated signal event
00107      bool source = iEvent.getByLabel(theSourceLabel,theHepMCProduct);
00108      if ( source ) { 
00109        myGenEvent = theHepMCProduct->GetEvent();
00110        // First rotate in case of beam crossing angle (except if done already)
00111        if ( theVertexGenerator ) { 
00112          TMatrixD* boost = theVertexGenerator->boost();
00113          if ( boost ) theHepMCProduct->boostToLab(boost,"momentum");
00114        }          
00115        myGenEvent = theHepMCProduct->GetEvent();
00116      } 
00117 
00118      fevt->setBeamSpot(BSPosition_);
00119        
00120      //     In case there is no HepMCProduct, seek a genParticle Candidate Collection
00121      bool genPart = false;
00122      if ( !myGenEvent ) { 
00123        // Look for the particle CandidateCollection
00124        Handle<reco::GenParticleCollection> genEvt;
00125        genPart = iEvent.getByLabel(theGenParticleLabel,genEvt);
00126        if ( genPart ) myGenParticles = &(*genEvt);
00127      }
00128        
00129      if ( !myGenEvent && !genPart )
00130        std::cout << "There is no generator input for this event, under " 
00131                  << "any form (HepMCProduct, genParticles)" << std::endl
00132                  << "Please check SourceLabel or GenParticleLabel" << std::endl;
00133        
00134      // Get the pile-up events from the pile-up producer
00135      // There might be no pile-up events, by the way, in that case, just continue
00136      Handle<HepMCProduct> thePileUpEvents;
00137      bool isPileUp = iEvent.getByLabel("famosPileUp","PileUpEvents",thePileUpEvents);
00138      thePUEvents = isPileUp ? thePileUpEvents->GetEvent() : 0;
00139       
00140    }//end else
00141 
00142    // pass the event to the Famos Manager for propagation and simulation
00143    if (myGenParticlesXF) {
00144      famosManager_->reconstruct(myGenParticlesXF);
00145    } else {
00146      famosManager_->reconstruct(myGenEvent,myGenParticles,thePUEvents);
00147    }
00148 
00149    // Set the vertex back to the HepMCProduct (except if it was smeared already)
00150    if ( myGenEvent ) { 
00151      if ( theVertexGenerator ) { 
00152        HepMC::FourVector theVertex(
00153                                    (theVertexGenerator->X()-theVertexGenerator->beamSpot().X()+BSPosition_.X())*10.,
00154                                    (theVertexGenerator->Y()-theVertexGenerator->beamSpot().Y()+BSPosition_.Y())*10.,
00155                                    (theVertexGenerator->Z()-theVertexGenerator->beamSpot().Z()+BSPosition_.Z())*10.,
00156                                    0.);
00157        if ( fabs(theVertexGenerator->Z()) > 1E-10 ) theHepMCProduct->applyVtxGen( &theVertex );
00158      }
00159    }
00160    
00161    CalorimetryManager * calo = famosManager_->calorimetryManager();
00162    TrajectoryManager * tracker = famosManager_->trackerManager();
00163    
00164    // Save everything in the edm::Event
00165    std::auto_ptr<edm::SimTrackContainer> p1(new edm::SimTrackContainer);
00166    std::auto_ptr<edm::SimTrackContainer> m1(new edm::SimTrackContainer);
00167    std::auto_ptr<edm::SimVertexContainer> p2(new edm::SimVertexContainer);
00168    std::auto_ptr<FSimVertexTypeCollection> v1(new FSimVertexTypeCollection);
00169    std::auto_ptr<edm::PSimHitContainer> p3(new edm::PSimHitContainer);
00170    std::auto_ptr<edm::PCaloHitContainer> p4(new edm::PCaloHitContainer);
00171    std::auto_ptr<edm::PCaloHitContainer> p5(new edm::PCaloHitContainer);
00172    std::auto_ptr<edm::PCaloHitContainer> p6(new edm::PCaloHitContainer); 
00173    std::auto_ptr<edm::PCaloHitContainer> p7(new edm::PCaloHitContainer);
00174    
00175    fevt->load(*p1,*m1);
00176    fevt->load(*p2);
00177    fevt->load(*v1);
00178    //   fevt->print();
00179    tracker->loadSimHits(*p3);
00180 
00181    //   fevt->print();
00182 
00183    if ( calo ) {  
00184      calo->loadFromEcalBarrel(*p4);
00185      calo->loadFromEcalEndcap(*p5);
00186      calo->loadFromPreshower(*p6);
00187      calo->loadFromHcal(*p7);
00188      // update the muon SimTracks
00189      calo->loadMuonSimTracks(*m1);
00190    }
00191 
00192    // Write muon first, to allow tracking particles to work... (pending MixingModule fix)
00193    if ( simulateMuons ) iEvent.put(m1,"MuonSimTracks");
00194    iEvent.put(p1);
00195    iEvent.put(p2);
00196    iEvent.put(p3,"TrackerHits");
00197    iEvent.put(v1,"VertexTypes");
00198    iEvent.put(p4,"EcalHitsEB");
00199    iEvent.put(p5,"EcalHitsEE");
00200    iEvent.put(p6,"EcalHitsES");
00201    iEvent.put(p7,"HcalHits");
00202 
00203 }
00204 
00205 DEFINE_FWK_MODULE(FamosProducer);
00206