CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HectorProducer.cc
Go to the documentation of this file.
1 // Framework headers
4 // SimpleConfigurable replacement
6 
7 // HepMC headers
8 #include "HepMC/GenEvent.h"
9 
10 // Hector headers
13 
14 // SimDataFormats headers
17 
18 #include <iostream>
19 #include <memory>
20 
21 using std::cout;
22 using std::endl;
23 
25 
26 
27  // TransportHector
28 
29  m_InTag = parameters.getParameter<std::string>("HepMCProductLabel") ;
30  m_verbosity = parameters.getParameter<bool>("Verbosity");
31  m_FP420Transport = parameters.getParameter<bool>("FP420Transport");
32  m_ZDCTransport = parameters.getParameter<bool>("ZDCTransport");
33 
34  produces<edm::HepMCProduct>();
35  produces<edm::LHCTransportLinkContainer>();
36 
37  hector = new Hector(parameters,
38  m_verbosity,
39  m_FP420Transport,
40  m_ZDCTransport);
41 
42 }
43 
45 
46  if(m_verbosity) {
47  LogDebug("HectorSetup") << "Delete HectorProducer"
48  << "Number of events analysed: " << eventsAnalysed;
49  }
50 
51 }
52 
54 
55  using namespace edm;
56  using namespace std;
57 
59 
60  Handle<HepMCProduct> HepMCEvt;
61  iEvent.getByLabel( m_InTag, HepMCEvt ) ;
62 
63  if ( !HepMCEvt.isValid() )
64  {
65  throw cms::Exception("InvalidReference")
66  << "Invalid reference to HepMCProduct\n";
67  }
68 
69  if ( HepMCEvt.provenance()->moduleLabel() == "LHCTransport" )
70  {
71  throw cms::Exception("LogicError")
72  << "HectorTrasported HepMCProduce already exists\n";
73  }
74 
75  evt_ = new HepMC::GenEvent( *HepMCEvt->GetEvent() );
77  if(m_FP420Transport) {
78  hector->clear();
79  hector->add( evt_ ,es);
81  }
82  if(m_ZDCTransport) {
83  hector->clear();
84  hector->add( evt_ ,es);
85  hector->filterZDC();
86 
87  hector->clear();
88  hector->add( evt_ ,es);
89  hector->filterD1();
90  }
92  if (m_verbosity) {
93  evt_->print();
94  }
95 
96  auto_ptr<HepMCProduct> NewProduct(new HepMCProduct()) ;
97  NewProduct->addHepMCData( evt_ ) ;
98 
99  iEvent.put( NewProduct ) ;
100 
101  auto_ptr<LHCTransportLinkContainer> NewCorrespondenceMap(new edm::LHCTransportLinkContainer() );
103  (*NewCorrespondenceMap).swap(thisLink);
104 
105  if ( m_verbosity ) {
106  for ( unsigned int i = 0; i < (*NewCorrespondenceMap).size(); i++)
107  LogDebug("HectorEventProcessing") << "Hector correspondence table: " << (*NewCorrespondenceMap)[i];
108  }
109 
110  iEvent.put( NewCorrespondenceMap );
111 
112 }
113 
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
dictionary parameters
Definition: Parameters.py:2
std::string m_InTag
void produce(edm::Event &iEvent, const edm::EventSetup &es)
this method will do the user analysis
void filterFP420()
Definition: Hector.cc:249
HectorProducer(edm::ParameterSet const &p)
default constructor
Definition: Hector.h:40
void clear()
Definition: Hector.cc:160
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
int eventsAnalysed
just to count events that have been analysed
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
virtual ~HectorProducer()
default destructor
void add(const HepMC::GenEvent *ev, const edm::EventSetup &es)
Definition: Hector.cc:183
std::vector< LHCTransportLink > & getCorrespondenceMap()
Definition: Hector.h:72
HepMC::GenEvent * addPartToHepMC(HepMC::GenEvent *event)
Definition: Hector.cc:495
std::vector< LHCTransportLink > LHCTransportLinkContainer
tuple cout
Definition: gather_cfg.py:121
HepMC::GenEvent * evt_
void filterD1()
Definition: Hector.cc:404
void filterZDC()
Definition: Hector.cc:333
void clearApertureFlags()
Definition: Hector.cc:154