CMS 3D CMS Logo

PyquenHadronizer.cc
Go to the documentation of this file.
1 /*
2  *
3  * Generates PYQUEN HepMC events
4  *
5 */
6 
7 #include <iostream>
8 #include <ctime>
9 
11 
16 
18 
23 
25 
26 #include "HepMC/IO_HEPEVT.h"
27 #include "HepMC/PythiaWrapper.h"
28 
29 using namespace gen;
30 using namespace edm;
31 using namespace std;
32 
33 HepMC::IO_HEPEVT hepevtio;
34 
37 
40  pset_(pset),
41  abeamtarget_(pset.getParameter<double>("aBeamTarget")),
42  angularspecselector_(pset.getParameter<int>("angularSpectrumSelector")),
43  bmin_(pset.getParameter<double>("bMin")),
44  bmax_(pset.getParameter<double>("bMax")),
45  bfixed_(pset.getParameter<double>("bFixed")),
46  cflag_(pset.getParameter<int>("cFlag")),
47  comenergy(pset.getParameter<double>("comEnergy")),
48  doquench_(pset.getParameter<bool>("doQuench")),
49  doradiativeenloss_(pset.getParameter<bool>("doRadiativeEnLoss")),
50  docollisionalenloss_(pset.getParameter<bool>("doCollisionalEnLoss")),
51  doIsospin_(pset.getParameter<bool>("doIsospin")),
52  protonSide_(pset.getUntrackedParameter<int>("protonSide", 0)),
53  embedding_(pset.getParameter<bool>("embeddingMode")),
54  evtPlane_(0),
55  nquarkflavor_(pset.getParameter<int>("qgpNumQuarkFlavor")),
56  qgpt0_(pset.getParameter<double>("qgpInitialTemperature")),
57  qgptau0_(pset.getParameter<double>("qgpProperTimeFormation")),
58  maxEventsToPrint_(pset.getUntrackedParameter<int>("maxEventsToPrint", 1)),
59  pythiaHepMCVerbosity_(pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity", false)),
60  pythiaPylistVerbosity_(pset.getUntrackedParameter<int>("pythiaPylistVerbosity", 0)),
61  pythia6Service_(new Pythia6Service(pset)),
62  filterType_(pset.getUntrackedParameter<string>("filterType", "None")) {
63  // Verbosity Level
64  // Valid PYLIST arguments are: 1, 2, 3, 5, 7, 11, 12, 13
65  LogDebug("PYLISTverbosity") << "Pythia PYLIST verbosity level = " << pythiaPylistVerbosity_ << endl;
66  // HepMC event verbosity Level
67  pythiaHepMCVerbosity_ = pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity", false);
68  LogDebug("HepMCverbosity") << "Pythia HepMC verbosity = " << pythiaHepMCVerbosity_ << endl;
69 
70  //Max number of events printed on verbosity level
71  maxEventsToPrint_ = pset.getUntrackedParameter<int>("maxEventsToPrint", 0);
72  LogDebug("Events2Print") << "Number of events to be printed = " << maxEventsToPrint_ << endl;
73 
74  if (embedding_) {
75  cflag_ = 0;
76  src_ = pset.getParameter<InputTag>("backgroundLabel");
77  }
79 
80  int cm = 1, va, vb, vc;
81  PYQVER(cm, va, vb, vc);
82 }
83 
84 //_____________________________________________________________________
86  // distructor
87  call_pystat(1);
88 
89  delete pythia6Service_;
90 }
91 
92 //_____________________________________________________________________
94 
95 //_____________________________________________________________________
97  HepMC::HeavyIon* hi = new HepMC::HeavyIon(1, // Ncoll_hard
98  -1, // Npart_proj
99  -1, // Npart_targ
100  1, // Ncoll
101  -1, // spectator_neutrons
102  -1, // spectator_protons
103  -1, // N_Nwounded_collisions
104  -1, // Nwounded_N_collisions
105  -1, // Nwounded_Nwounded_collisions
106  plfpar.bgen, // impact_parameter in [fm]
107  evtPlane_, // event_plane_angle
108  0, // eccentricity
109  -1 // sigma_inel_NN
110  );
111 
112  evt->set_heavy_ion(*hi);
113 
114  delete hi;
115 }
116 
117 //_____________________________________________________________________
120 
121  // Not possible to retrieve impact paramter and event plane info
122  // at this part, need to overwrite filter() in
123  // PyquenGeneratorFilter
124 
125  const edm::Event& e = getEDMEvent();
126 
127  if (embedding_) {
129  e.getByLabel(src_, input);
130  const HepMC::GenEvent* inev = input->GetEvent();
131  const HepMC::HeavyIon* hi = inev->heavy_ion();
132  if (hi) {
133  bfixed_ = hi->impact_parameter();
134  evtPlane_ = hi->event_plane_angle();
135  } else {
136  LogWarning("EventEmbedding") << "Background event does not have heavy ion record!";
137  }
138  }
139 
140  // Generate PYQUEN event
141  // generate single partonic PYTHIA jet event
142 
143  // Take into account whether it's a nn or pp or pn interaction
144  if (doIsospin_) {
145  string projN = "p";
146  string targN = "p";
147  if (protonSide_ != 1)
148  projN = nucleon();
149  if (protonSide_ != 2)
150  targN = nucleon();
151  call_pyinit("CMS", projN.data(), targN.data(), comenergy);
152  }
153  call_pyevnt();
154 
155  // call PYQUEN to apply parton rescattering and energy loss
156  // if doQuench=FALSE, it is pure PYTHIA
157  if (doquench_) {
159  edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN(" << abeamtarget_ << "," << cflag_ << "," << bfixed_
160  << ") ####";
161  } else {
162  edm::LogInfo("PYQUENinAction") << "##### Calling PYQUEN: QUENCHING OFF!! This is just PYTHIA !!!! ####";
163  }
164 
165  // call PYTHIA to finish the hadronization
166  pyexec_();
167 
168  // fill the HEPEVT with the PYJETS event record
169  call_pyhepc(1);
170 
171  // event information
172  HepMC::GenEvent* evt = hepevtio.read_next_event();
173 
174  evt->set_signal_process_id(pypars.msti[0]); // type of the process
175  evt->set_event_scale(pypars.pari[16]); // Q^2
176 
177  if (embedding_)
179  add_heavy_ion_rec(evt);
180 
181  HepMC::HEPEVT_Wrapper::check_hepevt_consistency();
182  //std::cout<<"Entries number: "<<HepMC::HEPEVT_Wrapper::number_entries() <<" Max. entries "<<HepMC::HEPEVT_Wrapper::max_number_entries()<<std::endl;
183 
184  event().reset(evt);
185 
186  return true;
187 }
188 
193 
194  //Proton to Nucleon fraction
195  pfrac_ = 1. / (1.98 + 0.015 * pow(abeamtarget_, 2. / 3));
196 
197  //initialize pythia
199 
200  //initilize pyquen
202 
203  return true;
204 }
205 
208 
209  // Call PYTHIA
210  call_pyinit("CMS", "p", "p", comenergy);
211 
212  return true;
213 }
214 
215 //_____________________________________________________________________
217  //Turn Hadronization Off whether or not there is quenching
218  // PYEXEC is called later anyway
219  string sHadOff("MSTP(111)=0");
220  gen::call_pygive(sHadOff);
221 
222  return true;
223 }
224 
225 //_________________________________________________________________
227  // PYQUEN initialization
228 
229  // angular emitted gluon spectrum selection
230  pyqpar.ianglu = angularspecselector_;
231 
232  // type of medium induced partonic energy loss
234  edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
235  pyqpar.ienglu = 0;
236  } else if (doradiativeenloss_) {
237  edm::LogInfo("PYQUENinRad") << "##### PYQUEN: Only RADIATIVE partonic energy loss ON ####";
238  pyqpar.ienglu = 1;
239  } else if (docollisionalenloss_) {
240  edm::LogInfo("PYQUENinColl") << "##### PYQUEN: Only COLLISIONAL partonic energy loss ON ####";
241  pyqpar.ienglu = 2;
242  } else {
243  edm::LogInfo("PYQUENinEnLoss") << "##### PYQUEN: Radiative AND Collisional partonic energy loss ON ####";
244  pyqpar.ienglu = 0;
245  }
246 
247  // number of active quark flavors in qgp
248  pyqpar.nfu = nquarkflavor_;
249  // initial temperature of QGP
250  pyqpar.T0u = qgpt0_;
251  // proper time of QGP formation
252  pyqpar.tau0u = qgptau0_;
253 
254  return true;
255 }
256 
258  int* dummy = nullptr;
259  double random = gen::pyr_(dummy);
260  const char* nuc = nullptr;
261  if (random > pfrac_)
262  nuc = "n";
263  else
264  nuc = "p";
265 
266  return nuc;
267 }
268 
270  double sinphi0 = sin(angle);
271  double cosphi0 = cos(angle);
272 
273  for (HepMC::GenEvent::vertex_iterator vt = evt->vertices_begin(); vt != evt->vertices_end(); ++vt) {
274  double x0 = (*vt)->position().x();
275  double y0 = (*vt)->position().y();
276  double z = (*vt)->position().z();
277  double t = (*vt)->position().t();
278 
279  double x = x0 * cosphi0 - y0 * sinphi0;
280  double y = y0 * cosphi0 + x0 * sinphi0;
281 
282  (*vt)->set_position(HepMC::FourVector(x, y, z, t));
283  }
284 
285  for (HepMC::GenEvent::particle_iterator vt = evt->particles_begin(); vt != evt->particles_end(); ++vt) {
286  double x0 = (*vt)->momentum().x();
287  double y0 = (*vt)->momentum().y();
288  double z = (*vt)->momentum().z();
289  double t = (*vt)->momentum().t();
290 
291  double x = x0 * cosphi0 - y0 * sinphi0;
292  double y = y0 * cosphi0 + x0 * sinphi0;
293 
294  (*vt)->set_momentum(HepMC::FourVector(x, y, z, t));
295  }
296 }
297 
298 bool PyquenHadronizer::declareStableParticles(const std::vector<int>& _pdg) {
299  std::vector<int> pdg = _pdg;
300  for (size_t i = 0; i < pdg.size(); i++) {
301  int pyCode = pycomp_(pdg[i]);
302  std::ostringstream pyCard;
303  pyCard << "MDCY(" << pyCode << ",1)=0";
304  std::cout << pyCard.str() << std::endl;
305  call_pygive(pyCard.str());
306  }
307 
308  return true;
309 }
310 
311 //____________________________________________________________________
312 
313 bool PyquenHadronizer::hadronize() { return false; }
314 
315 bool PyquenHadronizer::decay() { return true; }
316 
317 bool PyquenHadronizer::residualDecay() { return true; }
318 
320 
322 
323 const char* PyquenHadronizer::classname() const { return "gen::PyquenHadronizer"; }
gen::PyquenHadronizer::maxEventsToPrint_
unsigned int maxEventsToPrint_
Definition: PyquenHadronizer.h:88
plfpar
#define plfpar
Definition: PyquenWrapper.h:27
gen::PyquenHadronizer::nucleon
const char * nucleon()
Definition: PyquenHadronizer.cc:257
gen::PyquenHadronizer::decay
bool decay()
Definition: PyquenHadronizer.cc:315
gen::PyquenHadronizer::qgptau0_
double qgptau0_
Definition: PyquenHadronizer.h:86
gen::PyquenHadronizer::declareStableParticles
bool declareStableParticles(const std::vector< int > &)
Definition: PyquenHadronizer.cc:298
electrons_cff.bool
bool
Definition: electrons_cff.py:372
mps_fire.i
i
Definition: mps_fire.py:355
input
static const std::string input
Definition: EdmProvDump.cc:48
gen::PyquenHadronizer::classname
const char * classname() const
Definition: PyquenHadronizer.cc:323
MessageLogger.h
funct::false
false
Definition: Factorize.h:34
gen::PyquenHadronizer::embedding_
bool embedding_
Definition: PyquenHadronizer.h:78
PyquenWrapper.h
gen::PyquenHadronizer::evtPlane_
double evtPlane_
Definition: PyquenHadronizer.h:79
gen::pycomp_
int pycomp_(int &)
gen::PyquenHadronizer::filterType_
std::string filterType_
Definition: PyquenHadronizer.h:95
gen::PyquenHadronizer::pyquen_init
bool pyquen_init(const edm::ParameterSet &pset)
Definition: PyquenHadronizer.cc:226
gen::PyquenHadronizer::PyquenHadronizer
PyquenHadronizer(const edm::ParameterSet &)
Definition: PyquenHadronizer.cc:38
edm
HLT enums.
Definition: AlignableModifier.h:19
gen::pyr_
double pyr_(int *idummy)
Definition: Pythia6Service.cc:60
gather_cfg.cout
cout
Definition: gather_cfg.py:144
edm::LogInfo
Definition: MessageLogger.h:254
pypars
#define pypars
Definition: ExhumeHadronizer.cc:45
gen::PyquenHadronizer::rotateEvtPlane
void rotateEvtPlane(HepMC::GenEvent *evt, double angle)
Definition: PyquenHadronizer.cc:269
gen::FortranInstance::kFortranInstance
static const std::string kFortranInstance
Definition: FortranInstance.h:88
Pythia6Declarations.h
gen::PyquenHadronizer::abeamtarget_
double abeamtarget_
Definition: PyquenHadronizer.h:63
gen::PyquenHadronizer::finalizeEvent
void finalizeEvent()
Definition: PyquenHadronizer.cc:319
gen::PyquenHadronizer::qgpt0_
double qgpt0_
Definition: PyquenHadronizer.h:84
gen::Pythia6Service::setRandomEngine
void setRandomEngine(CLHEP::HepRandomEngine *v)
Definition: Pythia6Service.h:46
gen::PyquenHadronizer::doradiativeenloss_
bool doradiativeenloss_
if true perform quenching (default = true)
Definition: PyquenHadronizer.h:74
edm::Handle
Definition: AssociativeIterator.h:50
gen::BaseHadronizer
Definition: BaseHadronizer.h:46
PYQVER
#define PYQVER
Definition: PyquenWrapper.h:22
gen::PyquenHadronizer::docollisionalenloss_
bool docollisionalenloss_
DEFAULT = true.
Definition: PyquenHadronizer.h:75
gen::Pythia6Service::setGeneralParams
void setGeneralParams()
Definition: Pythia6Service.cc:157
gen::PyquenHadronizer::hadronize
bool hadronize()
Definition: PyquenHadronizer.cc:313
gen::PyquenHadronizer::pythiaPylistVerbosity_
unsigned int pythiaPylistVerbosity_
HepMC verbosity flag.
Definition: PyquenHadronizer.h:90
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
gen::pyexec_
void pyexec_()
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
HiGenEvtSelectorFactory::get
static BaseHiGenEvtSelector * get(std::string, const edm::ParameterSet &)
Definition: HiGenEvtSelectorFactory.cc:9
gen::PyquenHadronizer::doIsospin_
bool doIsospin_
DEFAULT = true.
Definition: PyquenHadronizer.h:76
PyquenHadronizer.h
gen::PyquenHadronizer::cflag_
int cflag_
fixed impact param (fm); valid only if cflag_=0
Definition: PyquenHadronizer.h:71
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
gen::Pythia6Service::setCSAParams
void setCSAParams()
Definition: Pythia6Service.cc:168
gen::FortranInstance::InstanceWrapper
Definition: FortranInstance.h:54
gen::BaseHadronizer::getEDMEvent
edm::Event & getEDMEvent() const
Definition: BaseHadronizer.h:69
gen::PyquenHadronizer::generatePartonsAndHadronize
bool generatePartonsAndHadronize()
Definition: PyquenHadronizer.cc:118
SharedResourceNames.h
gen
Definition: PythiaDecays.h:13
OrderedSet.t
t
Definition: OrderedSet.py:90
gen::PyquenHadronizer::initializeForInternalPartons
bool initializeForInternalPartons()
Definition: PyquenHadronizer.cc:206
FortranInstance.h
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::LogWarning
Definition: MessageLogger.h:141
gen::PyquenHadronizer::bmax_
double bmax_
min impact param (fm); valid only if cflag_!=0
Definition: PyquenHadronizer.h:69
gen::PyquenHadronizer::add_heavy_ion_rec
void add_heavy_ion_rec(HepMC::GenEvent *evt)
Definition: PyquenHadronizer.cc:96
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:670
edm::ParameterSet
Definition: ParameterSet.h:36
gen::PyquenHadronizer::pyqpythia_init
bool pyqpythia_init(const edm::ParameterSet &pset)
Definition: PyquenHadronizer.cc:216
gen::PyquenHadronizer::doquench_
bool doquench_
collision energy
Definition: PyquenHadronizer.h:73
GenEventInfoProduct.h
Event.h
gen::PyquenHadronizer::statistics
void statistics()
Definition: PyquenHadronizer.cc:321
gen::v
double v[5][pyjets_maxn]
Definition: Cascade2Hadronizer.cc:74
createfilelist.int
int
Definition: createfilelist.py:10
gen::PyquenHadronizer::angularspecselector_
unsigned int angularspecselector_
beam/target atomic mass number
Definition: PyquenHadronizer.h:64
gen::PyquenHadronizer::theSharedResources
static const std::vector< std::string > theSharedResources
Definition: PyquenHadronizer.h:53
gen::PyquenHadronizer::src_
edm::InputTag src_
Pythia PYLIST Verbosity flag.
Definition: PyquenHadronizer.h:93
gen::PyquenHadronizer::residualDecay
bool residualDecay()
Definition: PyquenHadronizer.cc:317
gen::PyquenHadronizer::pset_
edm::ParameterSet pset_
Definition: PyquenHadronizer.h:62
gen::PyquenHadronizer::bfixed_
double bfixed_
max impact param (fm); valid only if cflag_!=0
Definition: PyquenHadronizer.h:70
gen::PyquenHadronizer::pythia6Service_
Pythia6Service * pythia6Service_
Definition: PyquenHadronizer.h:94
hi
Definition: HiEvtPlaneList.h:38
gen::PyquenHadronizer::nquarkflavor_
unsigned int nquarkflavor_
Proton fraction in the nucleus.
Definition: PyquenHadronizer.h:82
std
Definition: JetResolutionObject.h:76
gen::BaseHadronizer::event
std::unique_ptr< HepMC::GenEvent > & event()
Definition: BaseHadronizer.h:86
gen::PyquenHadronizer::comenergy
double comenergy
centrality flag =0 fixed impact param, <>0 minbias
Definition: PyquenHadronizer.h:72
gen::PyquenHadronizer::pfrac_
double pfrac_
Definition: PyquenHadronizer.h:80
HiGenEvtSelectorFactory.h
edm::SharedResourceNames::kPythia6
static const std::string kPythia6
Definition: SharedResourceNames.h:26
PYQUEN
#define PYQUEN
Definition: PyquenWrapper.h:17
pyqpar
#define pyqpar
Definition: HydjetWrapper.h:112
Pythia6Service.h
angle
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
pdg
Definition: pdg_functions.h:28
gen::PyquenHadronizer::~PyquenHadronizer
~PyquenHadronizer() override
Definition: PyquenHadronizer.cc:85
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
gen::PyquenHadronizer::readSettings
bool readSettings(int)
Definition: PyquenHadronizer.cc:189
gen::Pythia6Service
Definition: Pythia6Service.h:24
gen::PyquenHadronizer::selector_
BaseHiGenEvtSelector * selector_
Definition: PyquenHadronizer.h:96
dummy
Definition: DummySelector.h:38
HepMCProduct.h
gen::PyquenHadronizer::pythiaHepMCVerbosity_
bool pythiaHepMCVerbosity_
Events to print if verbosity.
Definition: PyquenHadronizer.h:89
gen::PyquenHadronizer::bmin_
double bmin_
Definition: PyquenHadronizer.h:68
edm::Event
Definition: Event.h:73
hepevtio
HepMC::IO_HEPEVT hepevtio
Definition: PyquenHadronizer.cc:33
gen::call_pygive
bool call_pygive(const std::string &line)
Definition: ExhumeHadronizer.cc:64
edm::InputTag
Definition: InputTag.h:15
gen::PyquenHadronizer::doSetRandomEngine
void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
Definition: PyquenHadronizer.cc:93
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
gen::PyquenHadronizer::protonSide_
int protonSide_
Run n&p with proper ratios; if false, only p+p collisions.
Definition: PyquenHadronizer.h:77
reco::vertex_iterator
VertexRefVector::iterator vertex_iterator
iterator over a vector of references to Vertex objects in the same collection
Definition: VertexFwd.h:19