CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Cascade2Hadronizer.cc
Go to the documentation of this file.
2 
4 
5 
6 #include "HepMC/GenEvent.h"
7 #include "HepMC/PdfInfo.h"
8 #include "HepMC/PythiaWrapper6_4.h"
10 
11 #include "HepMC/HEPEVT_Wrapper.h"
12 #include "HepMC/IO_HEPEVT.h"
13 #include "HepMC/IO_GenEvent.h"
14 
17 
18 HepMC::IO_HEPEVT hepevtio;
19 
20 #include "HepPID/ParticleIDTranslations.hh"
21 
24 #include "CLHEP/Random/RandFlat.h"
26 
27 //-- Pythia6 routines and functionalities to pass around Pythia6 params
28 
31 
33 
34 using namespace edm;
35 using namespace std;
36 
37 #define debug 0
38 
39 CLHEP::RandFlat* fFlat_extern;
40 
41 extern "C" {
42 
43  double dcasrn_(int *idummy) {
44 
45  static int call = 0;
46 
47  double rdm_nb = fFlat_extern->fire();
48 
49  if(debug && ++call < 100) cout<<"dcasrn from c++, call: "<<call<<" random number: "<<rdm_nb<<endl;
50 
51  return rdm_nb;
52  }
53 }
54 
55 
56 namespace gen {
57 
59 
60  public:
62 
63  private:
64  void upInit() {
65  FortranCallback::getInstance()->fillHeader();
66  }
67 
68  void upEvnt() {
69  FortranCallback::getInstance()->fillEvent();
70  }
71 
72  bool upVeto() {
73  bool veto = false;
74  if (!hepeup_.nup) veto = true; //-- LHE Common Blocks
75  return(veto);
76  }
77 
78  };
79 
80  static struct {
81  int n, npad, k[5][pyjets_maxn];
82  double p[5][pyjets_maxn], v[5][pyjets_maxn];
83  } pyjets_local;
84 
85  Cascade2Hadronizer::Cascade2Hadronizer(edm::ParameterSet const& pset)
86  : BaseHadronizer(pset),
87  fPy6Service(new Pythia6ServiceWithCallback(pset)), //-- this will store py6 parameters for further settings
88 
89  //-- fRandomEngine(&getEngineReference()),
90 
91  //-- defined in GeneratorInterface/Core/src/RNDMEngineAccess.cc
92  //-- CLHEP::HepRandomEngine& gen::getEngineReference()
93  //-- { edm::Service<edm::RandomNumberGenerator> rng;
94  //-- return rng->getEngine(); }
95 
96  fComEnergy(pset.getParameter<double>("comEnergy")),
97  fextCrossSection(pset.getUntrackedParameter<double>("crossSection",-1.)),
98  fextCrossSectionError(pset.getUntrackedParameter<double>("crossSectionError",-1.)),
99  fFilterEfficiency(pset.getUntrackedParameter<double>("filterEfficiency",-1.)),
100 
101  fMaxEventsToPrint(pset.getUntrackedParameter<int>("maxEventsToPrint",0)),
102  fHepMCVerbosity(pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false)),
103  fPythiaListVerbosity(pset.getUntrackedParameter<int>("pythiaPylistVerbosity",0)),
104 
105  fDisplayPythiaBanner(pset.getUntrackedParameter<bool>("displayPythiaBanner",false)),
106  fDisplayPythiaCards(pset.getUntrackedParameter<bool>("displayPythiaCards",false)){
107 
108  fParameters = pset.getParameter<ParameterSet>("Cascade2Parameters");
109 
111  if(debug) cout<<"seed: "<<rng->mySeed()<<endl;
112  CLHEP::HepRandomEngine& engine = rng->getEngine();
113  fFlat = new CLHEP::RandFlat(engine);
115 
116  fConvertToPDG = false;
117  if(pset.exists("doPDGConvert"))
118  fConvertToPDG = pset.getParameter<bool>("doPDGConvert");
119 
120  //-- silence Pythia6 banner printout, unless display requested
121 
123  if(!call_pygive("MSTU(12)=12345")){
124  throw edm::Exception(edm::errors::Configuration,"PythiaError")
125  <<" Pythia did not accept MSTU(12)=12345";
126  }
127  }
128 
129  //-- silence printouts from PYGIVE, unless display requested
130 
131  if(!fDisplayPythiaCards){
132  if(!call_pygive("MSTU(13)=0")){
133  throw edm::Exception(edm::errors::Configuration,"PythiaError")
134  <<" Pythia did not accept MSTU(13)=0";
135  }
136  }
137 
138  //-- tmp stuff to deal with EvtGen corrupting pyjets
139  flushTmpStorage();
140  }
141 
143  if(fPy6Service != 0) delete fPy6Service;
144  }
145 
147 
148  pyjets_local.n = 0 ;
149  pyjets_local.npad = 0 ;
150 
151  for(int ip=0; ip<pyjets_maxn; ip++){
152  for(int i=0; i<5; i++){
153  pyjets_local.k[i][ip] = 0;
154  pyjets_local.p[i][ip] = 0.;
155  pyjets_local.v[i][ip] = 0.;
156  }
157  }
158  return;
159  }
160 
162 
163  pyjets_local.n = pyjets.n;
164  pyjets_local.npad = pyjets.npad;
165 
166  for(int ip=0; ip<pyjets_maxn; ip++){
167  for(int i=0; i<5; i++){
168  pyjets_local.k[i][ip] = pyjets.k[i][ip];
169  pyjets_local.p[i][ip] = pyjets.p[i][ip];
170  pyjets_local.v[i][ip] = pyjets.v[i][ip];
171  }
172  }
173  return ;
174  }
175 
177 
178  HepMC::PdfInfo pdf;
179 
180  //-- filling factorization "Q scale" now! pthat moved to binningValues()
181 
182  if(event()->signal_process_id() <= 0) event()->set_signal_process_id(pypars.msti[0]);
183  if(event()->event_scale() <=0 ) event()->set_event_scale(pypars.pari[22]);
184  if(event()->alphaQED() <= 0) event()->set_alphaQED(pyint1.vint[56]);
185  if(event()->alphaQCD() <= 0) event()->set_alphaQCD(pyint1.vint[57]);
186 
187  //-- get pdf info directly from Pythia6 and set it up into HepMC::GenEvent
188  //-- S. Mrenna: Prefer vint block
189 
190  if(pdf.id1() <= 0) pdf.set_id1(pyint1.mint[14] == 21 ? 0 : pyint1.mint[14]);
191  if(pdf.id2() <= 0) pdf.set_id2(pyint1.mint[15] == 21 ? 0 : pyint1.mint[15]);
192  if(pdf.x1() <= 0) pdf.set_x1(pyint1.vint[40]);
193  if(pdf.x2() <= 0) pdf.set_x2(pyint1.vint[41]);
194  if(pdf.pdf1() <= 0) pdf.set_pdf1(pyint1.vint[38] / pyint1.vint[40]);
195  if(pdf.pdf2() <= 0) pdf.set_pdf2(pyint1.vint[39] / pyint1.vint[41]);
196  if(pdf.scalePDF() <= 0) pdf.set_scalePDF(pyint1.vint[50]);
197 
198  event()->set_pdf_info(pdf) ;
199 
200  if(debug) {
201  cout<<"standard Py6 event weight: pyint1.vint[96]: "<<pyint1.vint[96]<<endl;
202  cout<<"event weight returned by PYEVWT: 1./(pyint1.vint[98]): "<<1./(pyint1.vint[98])<<endl;
203  }
204 
205  //-- this is "standard" Py6 event weight (corresponds to PYINT1/VINT(97))
206  // event()->weights().push_back(pyint1.vint[96]);
207 
208  //-- this is event weight as 1./VINT(99) (PYINT1/VINT(99) is returned by the PYEVWT)
209  // event()->weights().push_back(1./(pyint1.vint[98]));
210 
211  //-- all cascade events have weight = 1
212  event()->weights().push_back(1.);
213 
214  //-- now create the GenEventInfo product from the GenEvent and fill the missing pieces
215 
216  eventInfo().reset(new GenEventInfoProduct(event().get()));
217 
218  //-- in Pythia6 pthat is used to subdivide samples into different bins
219  //-- in LHE mode the binning is done by the external ME generator
220  //-- which is likely not pthat, so only filling it for Py6 internal mode
221 
222  eventInfo()->setBinningValues(vector<double>(1,pypars.pari[16]));
223 
224  //-- here we treat long-lived particles
225 
226  if (pydat1.mstj[21]==3 || pydat1.mstj[21]==4) imposeProperTime();
227 
228  //-- convert particle IDs Py6 -> PDG (if requested)
229 
230  if(fConvertToPDG){
231  for(HepMC::GenEvent::particle_iterator part = event()->particles_begin();
232  part != event()->particles_end(); ++part){
233  (*part)->set_pdg_id(HepPID::translatePythiatoPDT((*part)->pdg_id()));
234  }
235  }
236 
237  //-- service printouts, if requested
238  if(fMaxEventsToPrint > 0){
241  if(fHepMCVerbosity){
242  cout <<"Event process = "<<pypars.msti[0]<<endl<<"----------------------"<<endl;
243  event()->print();
244  }
245  }
246 
247  return;
248  }
249 
251  return false;
252  }
253 
255 
256  //-- grab Py6 instance
258 
260 
261  //-- skip event if py6 considers it bad
262  if(pyint1.mint[50] != 0 ){
263  event().reset();
264  return false;
265  }
266 
267  //-- generation of the event with CASCADE
268  call_event();
269 
270  //-- pythia pyhepc routine converts common PYJETS in common HEPEVT
271  call_pyhepc(1);
272 
273  //-- delete the created event from memory
274  event().reset(hepevtio.read_next_event());
275 
276  //-- this is to deal with post-gen tools & residualDecay() that may reuse PYJETS
277  flushTmpStorage();
278  fillTmpStorage();
279 
280  return true;
281  }
282 
284  return true;
285  }
286 
288 
289  //-- grab Py6 instance
291 
292  int NPartsBeforeDecays = pyjets_local.n ;
293  int NPartsAfterDecays = event().get()->particles_size();
294  int barcode = NPartsAfterDecays;
295 
296  // JVY: well, in principle, it's not a 100% fair to go up to NPartsBeforeDecays,
297  // because Photos will attach gamma's to existing vertices, i.e. in the middle
298  // of the event rather than at the end; but this will only shift pointers down,
299  // so we'll be going again over a few "original" particle...
300  // in the alternative, we may go all the way up to the beginning of the event
301  // and re-check if anything remains to decay, that's fine even if it'll take
302  // some extra CPU...
303 
304  for(int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- ){
305 
306  HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
307  int status = part->status();
308  if ( status != 1 ) continue; // check only "stable" particles,
309  // as some undecayed may still be marked as such
310  int pdgid = part->pdg_id();
311  int py6ptr = pycomp_( pdgid );
312 
313  if(pydat3.mdcy[0][py6ptr-1] != 1) continue; // particle is not expected to decay
314 
315  int py6id = HepPID::translatePDTtoPythia(pdgid);
316 
317  // first, will need to zero out, then fill up PYJETS
318  // I better do it directly (by hands) rather than via py1ent
319  // - to avoid (additional) mass smearing
320 
321  if(part->momentum().t() <= part->generated_mass()){
322  continue ; // e == m -> 0-momentum, nothing to decay...
323  }
324 
325  else{
326  pyjets.n = 0;
327 
328  for(int i=0; i<5; i++){
329  pyjets.k[i][0] = 0;
330  pyjets.p[i][0] = 0.;
331  pyjets.v[i][0] = 0.;
332  }
333 
334  pyjets.k[0][0] = 1;
335  pyjets.k[1][0] = py6id;
336  pyjets.p[4][0] = part->generated_mass();
337  pyjets.p[3][0] = part->momentum().t();
338  pyjets.p[0][0] = part->momentum().x();
339  pyjets.p[1][0] = part->momentum().y();
340  pyjets.p[2][0] = part->momentum().z();
341  HepMC::GenVertex* prod_vtx = part->production_vertex();
342  if ( !prod_vtx ) continue; // in principle, should never happen but...
343  pyjets.v[0][0] = prod_vtx->position().x();
344  pyjets.v[1][0] = prod_vtx->position().y();
345  pyjets.v[2][0] = prod_vtx->position().z();
346  pyjets.v[3][0] = prod_vtx->position().t();
347  pyjets.v[4][0] = 0.;
348  pyjets.n = 1;
349  pyjets.npad = pyjets_local.npad;
350  }
351 
352  //-- now call Py6 decay routine
353 
354  int parent = 1; // since we always pass to Py6 a single particle
355  pydecy_(parent);
356 
357  //-- now attach decay products to mother
358 
359  for(int iprt1=1; iprt1<pyjets.n; iprt1++){
360 
361  part->set_status( 2 );
362 
363  HepMC::GenVertex* DecVtx = new HepMC::GenVertex( HepMC::FourVector(pyjets.v[0][iprt1],
364  pyjets.v[1][iprt1],
365  pyjets.v[2][iprt1],
366  pyjets.v[3][iprt1]) );
367  DecVtx->add_particle_in( part ); // this will cleanup end_vertex if exists, replace with the new one
368  // I presume (vtx) barcode will be given automatically
369 
370  HepMC::FourVector pmom(pyjets.p[0][iprt1],pyjets.p[1][iprt1],
371  pyjets.p[2][iprt1],pyjets.p[3][iprt1] );
372 
373  int dstatus = 0;
374  if(pyjets.k[0][iprt1] >= 1 && pyjets.k[0][iprt1] <= 10) dstatus = 1;
375  else if(pyjets.k[0][iprt1] >= 11 && pyjets.k[0][iprt1] <= 20) dstatus = 2;
376  else if(pyjets.k[0][iprt1] >= 21 && pyjets.k[0][iprt1] <= 30) dstatus = 3;
377  else if(pyjets.k[0][iprt1] >= 31 && pyjets.k[0][iprt1] <= 100) dstatus = pyjets.k[0][iprt1];
378 
379  HepMC::GenParticle* daughter = new HepMC::GenParticle(pmom,
380  HepPID::translatePythiatoPDT(pyjets.k[1][iprt1]),
381  dstatus);
382  barcode++;
383  daughter->suggest_barcode( barcode );
384  DecVtx->add_particle_out( daughter );
385 
386  int iprt2;
387  for(iprt2=iprt1+1; iprt2<pyjets.n; iprt2++){ // the pointer is shifted by -1, c++ style
388 
389  if(pyjets.k[2][iprt2] != parent){
390  parent = pyjets.k[2][iprt2];
391  break; // another parent particle; reset & break the loop
392  }
393 
394  HepMC::FourVector pmomN(pyjets.p[0][iprt2],pyjets.p[1][iprt2],
395  pyjets.p[2][iprt2],pyjets.p[3][iprt2] );
396 
397  dstatus = 0;
398  if(pyjets.k[0][iprt2] >= 1 && pyjets.k[0][iprt2] <= 10) dstatus = 1;
399  else if(pyjets.k[0][iprt2] >= 11 && pyjets.k[0][iprt2] <= 20) dstatus = 2;
400  else if(pyjets.k[0][iprt2] >= 21 && pyjets.k[0][iprt2] <= 30) dstatus = 3;
401  else if(pyjets.k[0][iprt2] >= 31 && pyjets.k[0][iprt2] <= 100) dstatus = pyjets.k[0][iprt2];
402 
403  HepMC::GenParticle* daughterN = new HepMC::GenParticle(pmomN,
404  HepPID::translatePythiatoPDT( pyjets.k[1][iprt2] ),
405  dstatus);
406  barcode++;
407  daughterN->suggest_barcode( barcode );
408  DecVtx->add_particle_out( daughterN );
409  } //-- end iprt2 loop
410 
411  iprt1 = iprt2-1; // reset counter such that it doesn't go over the same child more than once
412  // don't forget to offset back into c++ counting, as it's already +1 forward
413 
414  event().get()->add_vertex( DecVtx );
415  } //-- end iprt1 loop
416  } //-- end loop over decay products
417 
418  // now restore the very original Py6 event record
419 
420  if(pyjets_local.n != pyjets.n){
421 
422  // restore pyjets to its state as it was before external decays -
423  // might have been jammed by action above or by py1ent calls in EvtGen
424 
425  pyjets.n = pyjets_local.n;
426  pyjets.npad = pyjets_local.npad;
427  for(int ip=0; ip<pyjets_local.n; ip++){
428  for(int i=0; i<5; i++){
429  pyjets.k[i][ip] = pyjets_local.k[i][ip];
430  pyjets.p[i][ip] = pyjets_local.p[i][ip];
431  pyjets.v[i][ip] = pyjets_local.v[i][ip];
432  }
433  }
434  }
435 
436  return true;
437  }
438 
440 
441  //-- grab Py6 instance
443 
445 
446  if (key == 0) fPy6Service->setCSAParams();
448 
449  return true;
450  }
451 
453  return false;
454  }
455 
457 
458  //-- grab Py6 instance
460 
461  //-- change standard parameters of JETSET/PYTHIA - replace call_pytcha()
462 
464 
468 
470 
471  //-- mstu(8) is set to NMXHEP in this dummy call (version >=6.404)
472  call_pyhepc(1);
473 
474  //-- initialise random number generator: has been changed to be CMSSW compliant
475  //-- dcasrn overloaded: call to rluxgo not needed anymore
476  //-- call_rluxgo(4,314159265,0,0);
477 
478  //-- initialise CASCADE parameters (default values)
479  call_casini();
480 
481  //-- Read the parameters and pass them to the common blocks
482  //-- call_steer();
483 
484  //-- initialise CASCADE parameters (user values)
485 
486  //-- retrieve all the different sets
487  vector<string> AllSets = fParameters.getParameter<vector<string> >("parameterSets");
488 
489  //-- loop over the different sets
490  for(unsigned i=0; i<AllSets.size();++i) {
491 
492  string Set = AllSets[i];
493  vector<string> Para_Set = fParameters.getParameter<vector<string> >(Set);
494 
495  //-- loop over all the parameters and stop in case of mistake
496  for(vector<string>::const_iterator itPara = Para_Set.begin(); itPara != Para_Set.end(); ++itPara) {
497 
498  if(!cascadeReadParameters(*itPara)) {
499  throw edm::Exception(edm::errors::Configuration,"CascadeError")
500  <<" Cascade did not accept the following parameter: \""<<*itPara<<"\""<<endl;
501  }
502  } //-- end loop over all the parameters
503  } //-- end loop over the different sets
504 
505  cainpu.plepin = -fComEnergy/2;
506  cainpu.ppin = fComEnergy/2;
507 
509 
510  //-- change standard parameters of CASCADE
511  call_cascha();
512 
513  //-- change standard parameters of JETSET/PYTHIA
514  //-- call_pytcha();
515 
516  //-- set up for running CASCADE (integration of the cross-section)
517  call_cascade();
518 
519  //-- print cross-section result from integration
520  call_caend(1);
521 
523 
525 
526  return true;
527  }
528 
530 
531  for(size_t i=0; i<pdg.size(); i++){
532  int PyID = HepPID::translatePDTtoPythia(pdg[i]);
533  int pyCode = pycomp_( PyID );
534 
535  if(pyCode > 0){
536  ostringstream pyCard ;
537  pyCard << "MDCY(" << pyCode << ",1)=0";
538  //-- cout << "pdg= " << pdg[i] << " " << pyCard.str() << endl;
539 
540  call_pygive(pyCard.str());
541  }
542  }
543 
544  return true;
545  }
546 
548 
549  // this is practically a copy/paste of the original code by J.Alcaraz,
550  // taken directly from PythiaSource
551 
552  int dumm=0;
553  HepMC::GenEvent::vertex_const_iterator vbegin = event()->vertices_begin();
554  HepMC::GenEvent::vertex_const_iterator vend = event()->vertices_end();
555  HepMC::GenEvent::vertex_const_iterator vitr = vbegin;
556 
557  for(; vitr != vend; ++vitr){
558 
559  HepMC::GenVertex::particle_iterator pbegin = (*vitr)->particles_begin(HepMC::children);
560  HepMC::GenVertex::particle_iterator pend = (*vitr)->particles_end(HepMC::children);
561  HepMC::GenVertex::particle_iterator pitr = pbegin;
562 
563  for(; pitr != pend; ++pitr){
564 
565  if((*pitr)->end_vertex()) continue;
566  if((*pitr)->status()!=1) continue;
567 
568  int pdgcode= abs((*pitr)->pdg_id());
569  // Do nothing if the particle is not expected to decay
570  if ( pydat3.mdcy[0][pycomp_(pdgcode)-1] !=1 ) continue;
571 
572  double ctau = pydat2.pmas[3][pycomp_(pdgcode)-1];
573  HepMC::FourVector mom = (*pitr)->momentum();
574  HepMC::FourVector vin = (*vitr)->position();
575  double x = 0.;
576  double y = 0.;
577  double z = 0.;
578  double t = 0.;
579  bool decayInRange = false;
580  while(!decayInRange){
581 
582  double unif_rand = fPy6Service->call(pyr_, &dumm);
583  // Value of 0 is excluded, so following line is OK
584  double proper_length = - ctau * log(unif_rand);
585  double factor = proper_length/mom.m();
586  x = vin.x() + factor * mom.px();
587  y = vin.y() + factor * mom.py();
588  z = vin.z() + factor * mom.pz();
589  t = vin.t() + factor * mom.e();
590 
591  // Decay must be happen outside a cylindrical region
592  if(pydat1.mstj[21]==4){
593  if (sqrt(x*x+y*y)>pydat1.parj[72] || fabs(z)>pydat1.parj[73]) decayInRange = true;
594  // Decay must be happen outside a given sphere
595  }
596  else if (pydat1.mstj[21]==3){
597  if (sqrt(x*x+y*y+z*z)>pydat1.parj[71]) decayInRange = true;
598  }
599  // Decay is always OK otherwise
600  else {
601  decayInRange = true;
602  }
603  }
604 
605  HepMC::GenVertex* vdec = new HepMC::GenVertex(HepMC::FourVector(x,y,z,t));
606  event()->add_vertex(vdec);
607  vdec->add_particle_in((*pitr));
608  }
609  }
610 
611  return;
612  }
613 
615 
617 
619 
620  if(!runInfo().internalXSec()){
621  double sigma_CMS = 0.001 * caeffic.avgi;
622  double error_CMS = 0.001 * caeffic.sd;
623  runInfo().setInternalXSec(GenRunInfoProduct::XSec(sigma_CMS,error_CMS));
624  }
625 
626  //-- print out generated event summary
627  call_caend(2);
628 
629  //-- write out some information from Pythia to the screen
630  call_pystat(1);
631 
632  return;
633  }
634 
635  const char* Cascade2Hadronizer::classname() const {
636  return "gen::Cascade2Hadronizer";
637  }
638 
639  //-- Read the parameters and pass them to the common blocks
640 
641  bool Cascade2Hadronizer::cascadeReadParameters(const string& ParameterString) {
642 
643  bool accepted = 1;
644 
645  if(!strncmp(ParameterString.c_str(),"KE",2))
646  caluco.ke = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
647 
648  else if(!strncmp(ParameterString.c_str(),"IRES(1)",7))
649  capar6.ires[0] = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
650 
651  else if(!strncmp(ParameterString.c_str(),"KP",2))
652  caluco.kp = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
653 
654  else if(!strncmp(ParameterString.c_str(),"IRES(2)",7))
655  capar6.ires[1] = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
656 
657  else if(!strncmp(ParameterString.c_str(),"NFRAG",5))
658  cainpu.nfrag = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
659 
660  else if(!strncmp(ParameterString.c_str(),"IPST",4))
661  cashower.ipst = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
662 
663  else if(!strncmp(ParameterString.c_str(),"IPSIPOL",7))
664  jpsi.ipsipol = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
665 
666  //-- from version 2.2.03 on
667  // else if(!strncmp(ParameterString.c_str(),"I23S",4))
668  // jpsi.i23s = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
669 
670  else if(!strncmp(ParameterString.c_str(),"IFPS",4))
671  cainpu.ifps = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
672 
673  else if(!strncmp(ParameterString.c_str(),"ITIMSHR",7))
674  casshwr.itimshr = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
675 
676  else if(!strncmp(ParameterString.c_str(),"IRUNAEM",7))
677  capar1.irunaem = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
678 
679  else if(!strncmp(ParameterString.c_str(),"IRUNA",5))
680  capar1.iruna = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
681 
682  else if(!strncmp(ParameterString.c_str(),"IQ2",3))
683  capar1.iq2 = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
684 
685  else if(!strncmp(ParameterString.c_str(),"IPRO",4))
686  capar1.ipro = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
687 
688  else if(!strncmp(ParameterString.c_str(),"NFLAV",5))
689  caluco.nflav = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
690 
691  else if(!strncmp(ParameterString.c_str(),"INTER",5))
692  cainpu.inter = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
693 
694  else if(!strncmp(ParameterString.c_str(),"IHFLA",5))
695  cahflav.ihfla = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
696 
697  else if(!strncmp(ParameterString.c_str(),"IRPA",4))
698  cascol.irpa = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
699 
700  else if(!strncmp(ParameterString.c_str(),"IRPB",4))
701  cascol.irpb = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
702 
703  else if(!strncmp(ParameterString.c_str(),"IRPC",4))
704  cascol.irpc = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
705 
706  else if(!strncmp(ParameterString.c_str(),"ICCFM",5))
707  casshwr.iccfm = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
708 
709  else if(!strncmp(ParameterString.c_str(),"IFINAL",6))
710  cainpu.ifinal = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
711 
712  else if(!strncmp(ParameterString.c_str(),"IGLU",4))
713  cagluon.iglu = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
714 
715  else if(!strncmp(ParameterString.c_str(),"IRspl",5))
716  casprre.irspl = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
717 
718  else if(!strncmp(ParameterString.c_str(),"PT2CUT",6))
719  captcut.pt2cut[capar1.ipro-1] = atof(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
720 
721  else if(!strncmp(ParameterString.c_str(),"NCB",3))
722  integr.ncb = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
723 
724  else if(!strncmp(ParameterString.c_str(),"ACC1",4))
725  integr.acc1 = atof(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
726 
727  else if(!strncmp(ParameterString.c_str(),"ACC2",4))
728  integr.acc2 = atof(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
729 
730  else if(!strncmp(ParameterString.c_str(),"SCALFAF",7))
731  scalf.scalfaf = atof(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
732 
733  else if(!strncmp(ParameterString.c_str(),"SCALFA",6))
734  scalf.scalfa = atof(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
735 
736  else accepted = 0;
737 
738  return accepted;
739  }
740 
742 
743  cout<<""<<endl;
744  cout<<"----------------------------------"<<endl;
745  cout<<"---- Cascade Parameters ----"<<endl;
746  cout<<"----------------------------------"<<endl;
747  cout<<""<<endl;
748 
749  cout<<"computation switch: "<<endl;
750  cout<<"number of calls per iteration for bases: "<<integr.ncb<<endl;
751  cout<<"relative precision for grid optimisation: "<<integr.acc1<<" %"<<endl;
752  cout<<"relative precision for integration: "<<integr.acc2<<" %"<<endl;
753  cout<<""<<endl;
754 
755  cout<<"kinematics: "<<endl;
756  cout<<"flavour code of beam 1: "<<caluco.ke<<endl;
757  cout<<"direct or resolved particle 1: "<<capar6.ires[0]<<endl;
758  cout<<"pz of incoming beam 1: "<<cainpu.plepin<<" GeV"<<endl;
759  cout<<"flavour code of beam 2: "<<caluco.kp<<endl;
760  cout<<"direct or resolved particle 2: "<<capar6.ires[1]<<endl;
761  cout<<"pz of incoming beam 2: "<<cainpu.ppin<<" GeV"<<endl;
762  cout<<"number of active flavors: "<<caluco.nflav<<endl;
763  cout<<""<<endl;
764 
765  cout<<"hard subprocess selection: "<<endl;
766  cout<<"hard subprocess number (IPRO): "<<capar1.ipro<<endl;
767  cout<<"IPRO = 10, switch to select QCD process g* g* -> q qbar: "<<cascol.irpa<<endl;
768  cout<<"IPRO = 10, switch to select QCD process g* g -> g g: "<<cascol.irpb<<endl;
769  cout<<"IPRO = 10, switch to select QCD process g* q -> g q: "<<cascol.irpc<<endl;
770  cout<<"flavor of heavy quark produced (IPRO = 11, 504, 514): "<<cahflav.ihfla<<endl;
771  // cout<<"select vector meson state: "<<jpsi.i23s<<endl;
772  cout<<"use matrix element including J/psi (Upsilon) polarisation: "<<jpsi.ipsipol<<endl;
773  cout<<"pt2 cut in matrix element for massless partons: "<<captcut.pt2cut[capar1.ipro-1]<<" GeV^2"<<endl;
774  cout<<""<<endl;
775 
776  cout<<"parton shower and fragmentation: "<<endl;
777  cout<<"switch for fragmentation: "<<cainpu.nfrag<<endl;
778  cout<<"switch for parton shower (0 off - 1 initial state - 2 final state - 3 initial & final state): "<<cainpu.ifps<<endl;
779  cout<<"switch for time like parton shower in initial state: "<<casshwr.itimshr<<endl;
780  cout<<"select CCFM (1) or DGLAP (0) evolution: "<<casshwr.iccfm<<endl;
781  cout<<"scale switch for final state parton shower: "<<cainpu.ifinal<<endl;
782  cout<<"scale factor for final state parton shower: "<<scalf.scalfaf<<endl;
783  cout<<"switch for proton remnant treatment: "<<casprre.irspl<<endl;
784  cout<<"keep track of intermediate state in parton shower: "<<cashower.ipst<<endl;
785  cout<<"mode of interaction for e p: "<<cainpu.inter<<endl;
786  cout<<""<<endl;
787 
788  cout<<"pdfs, couplings and scales: "<<endl;
789  cout<<"switch for running alphaem: "<<capar1.irunaem<<endl;
790  cout<<"switch for running alphas: "<<capar1.iruna<<endl;
791  cout<<"scale in alphas: "<<capar1.iq2<<endl;
792  cout<<"scale factor for scale in alphas: "<<scalf.scalfa<<endl;
793  cout<<"unintegrated pdf: "<<cagluon.iglu<<endl;
794  cout<<"path where updf are stored: "<<caspdf.pdfpath<<endl;
795  cout<<""<<endl;
796 
797  cout<<"center of mass energy, cross section, filter efficiency: "<<endl;
798  cout<<"center of mass energy: "<<fComEnergy<<" GeV"<<endl;
799  cout<<"external LO cross section: "<<fextCrossSection<<" +/- "<<fextCrossSectionError<<" pb"<<endl;
800  cout<<"filter efficiency: "<<fFilterEfficiency<<endl;
801  cout<<""<<endl;
802  }
803 
805 
806  cout<<""<<endl;
807  cout<<"----------------------------------"<<endl;
808  cout<<"---- Pythia6 Parameters ----"<<endl;
809  cout<<"----------------------------------"<<endl;
810  cout<<""<<endl;
811 
812  cout<<"charm mass: "<<pydat2.pmas[0][3]<<" GeV (default value = 1.5 GeV)"<<endl;
813  cout<<"bottom mass: "<<pydat2.pmas[0][4]<<" GeV (default value = 4.8 GeV)"<<endl;
814  cout<<"Higgs mass: "<<pydat2.pmas[0][24]<<" GeV (default value = 115 GeV)"<<endl;
815 
816  cout<<"lambda QCD: "<<pydat1.paru[111]<<" GeV (default value = 0.25 GeV)"<<endl;
817 
818  cout<<"alphas behaviour: "<<pydat1.mstu[110]<<" (default value = 1)"<<endl;
819  cout<<"nr of flavours wrt lambda QCD: "<<pydat1.mstu[111]<<" (default value = 5)"<<endl;
820  cout<<"min nr of flavours for alphas: "<<pydat1.mstu[112]<<" (default value = 3)"<<endl;
821  cout<<"max nr of flavours for alphas: "<<pydat1.mstu[113]<<" (default value = 5)"<<endl;
822 
823  cout<<"maximum angle settings: "<<pydat1.mstj[47]<<" (default value = 0)"<<endl;
824 
825  cout<<""<<endl;
826  }
827 
828 } //-- namespace gen
#define cainpu
#define pydat1
static struct gen::@232 pyjets_local
T getParameter(std::string const &) const
double dcasrn_(int *idummy)
int i
Definition: DBlmapReader.cc:9
void call(void(&fn)())
list parent
Definition: dbtoconf.py:74
#define captcut
Pythia6Service * fPy6Service
#define cahflav
#define caluco
struct HEPEUP_ hepeup_
bool call_pygive(const std::string &line)
bool exists(std::string const &parameterName) const
checks if a parameter exists
void call_casini()
#define abs(x)
Definition: mlp_lapack.h:159
#define cagluon
#define integr
bool declareStableParticles(const std::vector< int >)
void setInternalXSec(const XSec &xsec)
#define vend()
Definition: vmac.h:42
double double double z
std::auto_ptr< HepMC::GenEvent > & event()
#define cascol
#define pypars
GenRunInfoProduct & runInfo()
#define caspdf
void call_pylist(int mode)
CLHEP::RandFlat * fFlat_extern
Pythia6ServiceWithCallback(const edm::ParameterSet &pset)
edm::ParameterSet fParameters
void call_cascade()
T sqrt(T t)
Definition: SSEVec.h:46
void call_caend(int mode)
#define vbegin()
Definition: vmac.h:35
#define scalf
#define pyint1
void pydecy_(int &ip)
void setPYUPDAParams(bool afterPyinit)
const char * classname() const
#define capar6
std::auto_ptr< GenEventInfoProduct > & eventInfo()
bool cascadeReadParameters(const std::string &ParameterString)
int k[5][pyjets_maxn]
#define capar1
static FortranCallback * getInstance()
int pycomp_(int &)
part
Definition: HCALResponse.h:21
#define cashower
void call_event()
#define casshwr
void call_cascha()
#define jpsi
list key
Definition: combine.py:13
#define casprre
tuple cout
Definition: gather_cfg.py:121
Definition: DDAxes.h:10
void setExternalXSecNLO(const XSec &xsec)
tuple status
Definition: ntuplemaker.py:245
#define debug
Definition: MEtoEDMFormat.h:34
double pyr_(int *idummy)
void setExternalXSecLO(const XSec &xsec)
#define caeffic
mathSSE::Vec4< T > v
HepMC::IO_HEPEVT hepevtio