CMS 3D CMS Logo

BeamHaloProducer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <ctime>
3 
8 
12 
15 
16 using namespace edm;
17 using namespace std;
18 
19 #include "HepMC/IO_HEPEVT.h"
20 #include "HepMC/HEPEVT_Wrapper.h"
21 // #include "HepMC/ConvertHEPEVT.h"
22 // #include "HepMC/CBhepevt.h"
23 #include "HepMC/WeightContainer.h"
24 
25 #define KI_BHG_INIT ki_bhg_init_
26 extern "C" {
27 void KI_BHG_INIT(long& seed);
28 }
29 
30 #define BHSETPARAM bhsetparam_
31 extern "C" {
32 void BHSETPARAM(int* iparam, float* fparam, const char* cparam, int length);
33 }
34 
35 #define KI_BHG_FILL ki_bhg_fill_
36 extern "C" {
37 void KI_BHG_FILL(int& iret, float& weight);
38 }
39 
40 #define KI_BHG_STAT ki_bhg_stat_
41 extern "C" {
42 void KI_BHG_STAT(int& iret);
43 }
44 
45 // HepMC::ConvertHEPEVT conv;
46 //include "HepMC/HEPEVT_Wrapper.h"
47 static HepMC::HEPEVT_Wrapper wrapper;
48 static HepMC::IO_HEPEVT conv;
49 
51  int iret = 0;
52  call_ki_bhg_stat(iret);
53 }
54 
55 BeamHaloProducer::BeamHaloProducer(const ParameterSet& pset) : evt(nullptr), isInitialized_(false) {
56  int iparam[8];
57  float fparam[4];
58  std::string cparam;
59  // -- from bhgctrl.inc
60  iparam[0] = pset.getUntrackedParameter<int>("GENMOD");
61  iparam[1] = pset.getUntrackedParameter<int>("LHC_B1");
62  iparam[2] = pset.getUntrackedParameter<int>("LHC_B2");
63  iparam[3] = pset.getUntrackedParameter<int>("IW_MUO");
64  iparam[4] = pset.getUntrackedParameter<int>("IW_HAD");
65  iparam[5] = 9999999;
66  iparam[6] = pset.getUntrackedParameter<int>("OFFSET", 0);
67  iparam[7] = pset.getUntrackedParameter<int>("shift_bx");
68 
69  fparam[0] = (float)pset.getUntrackedParameter<double>("EG_MIN");
70  fparam[1] = (float)pset.getUntrackedParameter<double>("EG_MAX");
71 
72  fparam[2] = (float)pset.getUntrackedParameter<double>("BXNS");
73  fparam[3] = (float)pset.getUntrackedParameter<double>("W0", 1.0);
74 
75  cparam = pset.getUntrackedParameter<std::string>("G3FNAME", "input.txt");
76  call_bh_set_parameters(iparam, fparam, cparam);
77 
78  produces<HepMCProduct>("unsmeared");
79  produces<GenEventInfoProduct>();
80  produces<GenRunInfoProduct, Transition::EndRun>();
81 
82  usesResource("BeamHaloProducer");
83 
84  cout << "BeamHaloProducer: starting event generation ... " << endl;
85 }
86 
88 
89 void BeamHaloProducer::setRandomEngine(CLHEP::HepRandomEngine* v) { _BeamHalo_randomEngine = v; }
90 
92  if (!isInitialized_) {
93  isInitialized_ = true;
94  RandomEngineSentry<BeamHaloProducer> randomEngineSentry(this, lumi.index());
95 
96  // -- initialisation
97  long seed = 1; // This seed is not actually used
99  }
100 }
101 
103  RandomEngineSentry<BeamHaloProducer> randomEngineSentry(this, e.streamID());
104 
105  // cout << "in produce " << endl;
106 
107  // unique_ptr<HepMCProduct> bare_product(new HepMCProduct());
108 
109  // cout << "apres autoptr " << endl;
110 
111  int iret = 0;
112  float weight = 0;
113  call_ki_bhg_fill(iret, weight);
114 
115  // Throw an exception if call_ki_bhg_fill(...) fails. Use the EventCorruption
116  // exception since it maps onto SkipEvent which is what we want to do here.
117 
118  if (iret < 0)
120  << "BeamHaloProducer: function call_ki_bhg_fill returned " << iret << endl;
121 
122  // cout << "apres fortran " << endl;
123 
124  // HepMC::GenEvent* evt = conv.getGenEventfromHEPEVT();
125  // HepMC::GenEvent* evt = conv.read_next_event(); seems to be broken (?)
126  evt = new HepMC::GenEvent();
127 
128  for (int theindex = 1; theindex <= wrapper.number_entries(); theindex++) {
129  HepMC::GenVertex* Vtx = new HepMC::GenVertex(
130  HepMC::FourVector(wrapper.x(theindex), wrapper.y(theindex), wrapper.z(theindex), wrapper.t(theindex)));
131  HepMC::FourVector p(wrapper.px(theindex), wrapper.py(theindex), wrapper.pz(theindex), wrapper.e(theindex));
132  HepMC::GenParticle* Part = new HepMC::GenParticle(p, wrapper.id(theindex), wrapper.status(theindex));
133  Vtx->add_particle_out(Part);
134  evt->add_vertex(Vtx);
135  }
136 
137  evt->set_event_number(e.id().event());
138 
139  HepMC::WeightContainer& weights = evt->weights();
140  weights.push_back(weight);
141  // evt->print();
142  std::unique_ptr<HepMCProduct> CMProduct(new HepMCProduct());
143  if (evt)
144  CMProduct->addHepMCData(evt);
145  e.put(std::move(CMProduct), "unsmeared");
146 
147  unique_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(evt));
148  e.put(std::move(genEventInfo));
149 }
150 
152  // just create an empty product
153  // to keep the EventContent definitions happy
154  // later on we might put the info into the run info that this is a PGun
155  unique_ptr<GenRunInfoProduct> genRunInfo(new GenRunInfoProduct());
156  run.put(std::move(genRunInfo));
157 }
158 
159 bool BeamHaloProducer::call_bh_set_parameters(int* ival, float* fval, const std::string cval_string) {
160  BHSETPARAM(ival, fval, cval_string.c_str(), cval_string.length());
161  return true;
162 }
163 
165  KI_BHG_INIT(seed);
166  return true;
167 }
168 
170  KI_BHG_FILL(iret, weight);
171  return true;
172 }
173 
175  KI_BHG_STAT(iret);
176  return true;
177 }
GenEventInfoProduct
Definition: GenEventInfoProduct.h:17
edm::BeamHaloProducer::clear
void clear()
Definition: BeamHaloProducer.cc:87
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
funct::false
false
Definition: Factorize.h:29
edm::BeamHaloProducer::setRandomEngine
void setRandomEngine(CLHEP::HepRandomEngine *v)
Definition: BeamHaloProducer.cc:89
edm::BeamHaloProducer::~BeamHaloProducer
~BeamHaloProducer() override
Destructor.
Definition: BeamHaloProducer.cc:50
conv
static HepMC::IO_HEPEVT conv
Definition: BeamHaloProducer.cc:48
edm::LuminosityBlock
Definition: LuminosityBlock.h:50
edm::Run
Definition: Run.h:45
edm
HLT enums.
Definition: AlignableModifier.h:19
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
edm::BeamHaloProducer::beginLuminosityBlock
void beginLuminosityBlock(LuminosityBlock const &, EventSetup const &) override
Definition: BeamHaloProducer.cc:91
ZgammaFilter_cfi.HepMCProduct
HepMCProduct
Definition: ZgammaFilter_cfi.py:9
KI_BHG_INIT
#define KI_BHG_INIT
Definition: BeamHaloProducer.cc:25
gather_cfg.cout
cout
Definition: gather_cfg.py:144
edm::errors::EventCorruption
Definition: EDMException.h:43
wrapper
static HepMC::HEPEVT_Wrapper wrapper
Definition: BeamHaloProducer.cc:47
findQualityFiles.v
v
Definition: findQualityFiles.py:179
GenRunInfoProduct.h
PYR.h
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
edm::BeamHaloProducer::produce
void produce(Event &e, const EventSetup &es) override
Definition: BeamHaloProducer.cc:102
fileCollector.seed
seed
Definition: fileCollector.py:127
HLT_FULL_cff.weights
weights
Definition: HLT_FULL_cff.py:99178
Run.h
KI_BHG_FILL
#define KI_BHG_FILL
Definition: BeamHaloProducer.cc:35
GenRunInfoProduct
Definition: GenRunInfoProduct.h:8
edm::BeamHaloProducer::evt
HepMC::GenEvent * evt
Definition: BeamHaloProducer.h:45
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::ParameterSet
Definition: ParameterSet.h:47
GenEventInfoProduct.h
Event.h
edm::BeamHaloProducer::call_bh_set_parameters
bool call_bh_set_parameters(int *ival, float *fval, const std::string cval_string)
Definition: BeamHaloProducer.cc:159
BHSETPARAM
#define BHSETPARAM
Definition: BeamHaloProducer.cc:30
edm::BeamHaloProducer::endRunProduce
void endRunProduce(Run &r, const EventSetup &es) override
Definition: BeamHaloProducer.cc:151
edm::EventSetup
Definition: EventSetup.h:57
edm::BeamHaloProducer::BeamHaloProducer
BeamHaloProducer(const ParameterSet &)
Constructor.
Definition: BeamHaloProducer.cc:55
_BeamHalo_randomEngine
CLHEP::HepRandomEngine * _BeamHalo_randomEngine
Definition: PYR.cc:3
edm::BeamHaloProducer::call_ki_bhg_fill
bool call_ki_bhg_fill(int &iret, float &weight)
Definition: BeamHaloProducer.cc:169
edm::BeamHaloProducer::call_ki_bhg_init
bool call_ki_bhg_init(long &seed)
Definition: BeamHaloProducer.cc:164
GenParticle.GenParticle
GenParticle
Definition: GenParticle.py:18
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
writedatasetfile.run
run
Definition: writedatasetfile.py:27
RandomEngineSentry.h
genParticles2HepMC_cfi.genEventInfo
genEventInfo
Definition: genParticles2HepMC_cfi.py:6
edm::RandomEngineSentry
Definition: RandomEngineSentry.h:28
KI_BHG_STAT
#define KI_BHG_STAT
Definition: BeamHaloProducer.cc:40
Exception
Definition: hltDiff.cc:246
edm::BeamHaloProducer::call_ki_bhg_stat
bool call_ki_bhg_stat(int &iret)
Definition: BeamHaloProducer.cc:174
Exception.h
HepMCProduct.h
edm::Event
Definition: Event.h:73
lumi
Definition: LumiSectionData.h:20
edm::BeamHaloProducer::isInitialized_
bool isInitialized_
Definition: BeamHaloProducer.h:56
weight
Definition: weight.py:1
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
BeamHaloProducer.h
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37