CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
Generator Class Reference

#include <Generator.h>

Public Member Functions

virtual const double eventWeight () const
 
 Generator (const edm::ParameterSet &p)
 
virtual const HepMC::GenEvent * genEvent () const
 
virtual const
math::XYZTLorentzVector
genVertex () const
 
void HepMC2G4 (const HepMC::GenEvent *g, G4Event *e)
 
void nonBeamEvent2G4 (const HepMC::GenEvent *g, G4Event *e)
 
void setGenEvent (const HepMC::GenEvent *inpevt)
 
virtual ~Generator ()
 

Private Member Functions

void particleAssignDaughters (G4PrimaryParticle *p, HepMC::GenParticle *hp, double length)
 
bool particlePassesPrimaryCuts (const G4ThreeVector &p) const
 
void setGenId (G4PrimaryParticle *p, int id) const
 

Private Attributes

HepMC::GenEvent * evt_
 
bool fEtaCuts
 
LumiMonitorFilterfLumiFilter
 
bool fPCuts
 
bool fPDGFilter
 
bool fPhiCuts
 
bool fPtransCut
 
std::vector< int > pdgFilter
 
bool pdgFilterSel
 
double theDecLenCut
 
double theDecRCut2
 
double theEtaCutForHector
 
double theMaxEtaCut
 
double theMaxPCut
 
double theMaxPhiCut
 
double theMinEtaCut
 
double theMinPCut
 
double theMinPhiCut
 
double theMinPtCut2
 
int verbose
 
math::XYZTLorentzVectorvtx_
 
double weight_
 
double Z_hector
 
double Z_lmax
 
double Z_lmin
 

Detailed Description

Definition at line 19 of file Generator.h.

Constructor & Destructor Documentation

Generator::Generator ( const edm::ParameterSet p)

Definition at line 22 of file Generator.cc.

References LumiMonitorFilter::Describe(), edm::ParameterSet::exists(), create_public_lumi_plots::exp, fEtaCuts, fLumiFilter, fPCuts, fPDGFilter, fPhiCuts, fPtransCut, edm::ParameterSet::getParameter(), cuy::ii, fjr2json::lumi, pdgFilter, pdgFilterSel, theDecLenCut, theDecRCut2, theEtaCutForHector, theMaxEtaCut, theMinEtaCut, theMinPCut, theMinPtCut2, Z_hector, Z_lmax, and Z_lmin.

22  :
23  fPCuts(p.getParameter<bool>("ApplyPCuts")),
24  fPtransCut(p.getParameter<bool>("ApplyPtransCut")),
25  fEtaCuts(p.getParameter<bool>("ApplyEtaCuts")),
26  fPhiCuts(p.getParameter<bool>("ApplyPhiCuts")),
27  theMinPhiCut(p.getParameter<double>("MinPhiCut")), // in radians (CMS standard)
28  theMaxPhiCut(p.getParameter<double>("MaxPhiCut")),
29  theMinEtaCut(p.getParameter<double>("MinEtaCut")),
30  theMaxEtaCut(p.getParameter<double>("MaxEtaCut")),
31  theMinPCut(p.getParameter<double>("MinPCut")), // in GeV (CMS standard)
32  theMaxPCut(p.getParameter<double>("MaxPCut")),
33  theEtaCutForHector(p.getParameter<double>("EtaCutForHector")),
34  verbose(p.getUntrackedParameter<int>("Verbosity",0)),
35  fLumiFilter(0),
36  evt_(0),
37  vtx_(0),
38  weight_(0),
39  Z_lmin(0),
40  Z_lmax(0),
41  Z_hector(0),
42  pdgFilterSel(false),
43  fPDGFilter(false)
44 {
45  bool lumi = p.getParameter<bool>("ApplyLumiMonitorCuts");
46  if(lumi) { fLumiFilter = new LumiMonitorFilter(); }
47 
48  double theRDecLenCut = p.getParameter<double>("RDecLenCut")*cm;
49  theDecRCut2 = theRDecLenCut*theRDecLenCut;
50 
52 
53  double theDecLenCut = p.getParameter<double>("LDecLenCut")*cm;
54 
55  pdgFilter.resize(0);
56  if ( p.exists("PDGselection") ) {
57  pdgFilterSel = (p.getParameter<edm::ParameterSet>("PDGselection")).
58  getParameter<bool>("PDGfilterSel");
59  pdgFilter = (p.getParameter<edm::ParameterSet>("PDGselection")).
60  getParameter<std::vector< int > >("PDGfilter");
61  if(0 < pdgFilter.size()) {
62  fPDGFilter = true;
63  for ( unsigned int ii = 0; ii < pdgFilter.size(); ++ii) {
64  if (pdgFilterSel) {
65  edm::LogWarning("SimG4CoreGenerator")
66  << " *** Selecting only PDG ID = " << pdgFilter[ii];
67  } else {
68  edm::LogWarning("SimG4CoreGenerator")
69  << " *** Filtering out PDG ID = " << pdgFilter[ii];
70  }
71  }
72  }
73  }
74 
75  if(fEtaCuts){
76  Z_lmax = theRDecLenCut*((1 - exp(-2*theMaxEtaCut) )/( 2*exp(-theMaxEtaCut)));
77  Z_lmin = theRDecLenCut*((1 - exp(-2*theMinEtaCut) )/( 2*exp(-theMinEtaCut)));
78  }
79 
80  Z_hector = theRDecLenCut*((1 - exp(-2*theEtaCutForHector))
81  / ( 2*exp(-theEtaCutForHector) ) );
82 
83  edm::LogInfo("SimG4CoreGenerator")
84  << " Rdecaycut= " << theRDecLenCut/cm
85  << " cm; Zdecaycut= " << theDecLenCut/cm
86  << "Z_min= " << Z_lmin/cm << " cm; Z_max= " << Z_lmax
87  << " cm; Z_hector = " << Z_hector << " cm\n"
88  << "ApplyPCuts: " << fPCuts << " ApplyPtransCut: " << fPtransCut
89  << " ApplyEtaCuts: " << fEtaCuts
90  << " ApplyPhiCuts: " << fPhiCuts
91  << " ApplyLumiMonitorCuts: " << lumi;
92  if(lumi) { fLumiFilter->Describe(); }
93 }
double theDecRCut2
Definition: Generator.h:53
bool fPhiCuts
Definition: Generator.h:45
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
double theMinPtCut2
Definition: Generator.h:51
double theEtaCutForHector
Definition: Generator.h:54
double Z_lmax
Definition: Generator.h:61
tuple lumi
Definition: fjr2json.py:35
double theMaxPhiCut
Definition: Generator.h:47
bool fPtransCut
Definition: Generator.h:43
bool exists(std::string const &parameterName) const
checks if a parameter exists
double theDecLenCut
Definition: Generator.h:55
double theMinEtaCut
Definition: Generator.h:48
bool pdgFilterSel
Definition: Generator.h:63
int ii
Definition: cuy.py:588
double theMinPCut
Definition: Generator.h:50
LumiMonitorFilter * fLumiFilter
Definition: Generator.h:57
bool fPDGFilter
Definition: Generator.h:64
double theMinPhiCut
Definition: Generator.h:46
double weight_
Definition: Generator.h:60
double theMaxEtaCut
Definition: Generator.h:49
double theMaxPCut
Definition: Generator.h:52
void Describe() const
math::XYZTLorentzVector * vtx_
Definition: Generator.h:59
HepMC::GenEvent * evt_
Definition: Generator.h:58
double Z_lmin
Definition: Generator.h:61
int verbose
Definition: Generator.h:56
std::vector< int > pdgFilter
Definition: Generator.h:62
bool fEtaCuts
Definition: Generator.h:44
double Z_hector
Definition: Generator.h:61
bool fPCuts
Definition: Generator.h:42
Generator::~Generator ( )
virtual

Definition at line 95 of file Generator.cc.

References fLumiFilter.

96 {
97  delete fLumiFilter;
98 }
LumiMonitorFilter * fLumiFilter
Definition: Generator.h:57

Member Function Documentation

virtual const double Generator::eventWeight ( ) const
inlinevirtual

Definition at line 31 of file Generator.h.

References weight_.

Referenced by RunManagerMTWorker::produce(), and RunManager::produce().

31 { return weight_; }
double weight_
Definition: Generator.h:60
virtual const HepMC::GenEvent* Generator::genEvent ( ) const
inlinevirtual

Definition at line 29 of file Generator.h.

References evt_.

Referenced by RunManagerMTWorker::produce(), and RunManager::produce().

29 { return evt_; }
HepMC::GenEvent * evt_
Definition: Generator.h:58
virtual const math::XYZTLorentzVector* Generator::genVertex ( ) const
inlinevirtual

Definition at line 30 of file Generator.h.

References vtx_.

Referenced by RunManagerMTWorker::produce(), and RunManager::produce().

30 { return vtx_; }
math::XYZTLorentzVector * vtx_
Definition: Generator.h:59
void Generator::HepMC2G4 ( const HepMC::GenEvent *  g,
G4Event *  e 
)

Definition at line 100 of file Generator.cc.

References RecoTauCleanerPlugins::charge, fEtaCuts, spr::find(), fLumiFilter, fPCuts, fPDGFilter, fPhiCuts, fPtransCut, GenParticle::GenParticle, GeV, LumiMonitorFilter::isGoodForLumiMonitor(), LogDebug, AlCaHLTBitMon_ParallelJobs::p, particleAssignDaughters(), pdgFilter, pdgFilterSel, phi, setGenId(), mathSSE::sqrt(), ntuplemaker::status, theDecRCut2, theMaxPCut, theMaxPhiCut, theMinPtCut2, vtx_, weight_, Z_hector, Z_lmax, and Z_lmin.

Referenced by RunManagerMTWorker::generateEvent(), and RunManager::generateEvent().

101 {
102 
103  HepMC::GenEvent *evt=new HepMC::GenEvent(*evt_orig);
104 
105  if ( *(evt->vertices_begin()) == 0 ) {
106  throw SimG4Exception("SimG4CoreGenerator: Corrupted Event - GenEvent with no vertex");
107  }
108 
109  if (evt->weights().size() > 0) {
110 
111  weight_ = evt->weights()[0] ;
112  for (unsigned int iw=1; iw<evt->weights().size(); ++iw) {
113 
114  // terminate if the vector of weights contains a zero-weight
115  if ( evt->weights()[iw] <= 0 ) break;
116  weight_ *= evt->weights()[iw] ;
117  }
118  }
119 
120  if (vtx_ != 0) { delete vtx_; }
121  vtx_ = new math::XYZTLorentzVector((*(evt->vertices_begin()))->position().x(),
122  (*(evt->vertices_begin()))->position().y(),
123  (*(evt->vertices_begin()))->position().z(),
124  (*(evt->vertices_begin()))->position().t());
125 
126  if(verbose > 0) {
127  evt->print(G4cout);
128  LogDebug("SimG4CoreGenerator") << "Primary Vertex = ("
129  << vtx_->x() << ","
130  << vtx_->y() << ","
131  << vtx_->z() << ")";
132  }
133 
134  unsigned int ng4vtx = 0;
135 
136  for(HepMC::GenEvent::vertex_const_iterator vitr= evt->vertices_begin();
137  vitr != evt->vertices_end(); ++vitr ) {
138 
139  // loop for vertex, is it a real vertex?
140  G4bool qvtx=false;
141  HepMC::GenVertex::particle_iterator pitr;
142  for (pitr= (*vitr)->particles_begin(HepMC::children);
143  pitr != (*vitr)->particles_end(HepMC::children); ++pitr) {
144 
145  // Particles which are not decayed by generator
146  if (1 == (*pitr)->status()) {
147  qvtx = true;
148  if (verbose > 2) LogDebug("SimG4CoreGenerator")
149  << "GenVertex barcode = " << (*vitr)->barcode()
150  << " selected for GenParticle barcode = " << (*pitr)->barcode();
151  break;
152  }
153  // The selection is made considering if the partcile with status = 2
154  // have the end_vertex with a radius greater than the radius of beampipe
155  // cilinder (no requirement on the Z of the vertex is applyed).
156  // Status 104 (metastable R-hadrons in Pythia8) is treated the same way.
157  else if (2 == (*pitr)->status() || 104 == (*pitr)->status()) {
158 
159  if ( (*pitr)->end_vertex() != 0 ) {
160  double xx = (*pitr)->end_vertex()->position().x();
161  double yy = (*pitr)->end_vertex()->position().y();
162  double r_dd = xx*xx+yy*yy;
163  if (r_dd > theDecRCut2) {
164  qvtx = true;
165  if (verbose > 2) LogDebug("SimG4CoreGenerator")
166  << "GenVertex barcode = " << (*vitr)->barcode()
167  << " selected for GenParticle barcode = "
168  << (*pitr)->barcode() << " radius = " << std::sqrt(r_dd);
169  break;
170  }
171  } else {
172  // particles with status 2 without end_vertex are
173  // equivalent to stable
174  qvtx = true;
175  }
176  }
177  }
178 
179  // if this vertex is inside fiductial volume inside the beam pipe
180  // and has no long-lived secondary the vertex is not saved
181  if (!qvtx) {
182  continue;
183  }
184 
185  double x1 = (*vitr)->position().x()*mm;
186  double y1 = (*vitr)->position().y()*mm;
187  double z1 = (*vitr)->position().z()*mm;
188  double t1 = (*vitr)->position().t()*mm/c_light;
189 
190  G4PrimaryVertex* g4vtx = new G4PrimaryVertex(x1, y1, z1, t1);
191 
192  for (pitr= (*vitr)->particles_begin(HepMC::children);
193  pitr != (*vitr)->particles_end(HepMC::children); ++pitr){
194 
195  // Filter on allowed particle species if required
196  if (fPDGFilter) {
197  std::vector<int>::iterator it = find(pdgFilter.begin(),pdgFilter.end(),
198  (*pitr)->pdg_id());
199  if ( (it != pdgFilter.end() && !pdgFilterSel)
200  || ( it == pdgFilter.end() && pdgFilterSel ) ) {
201  if (verbose > 2) LogDebug("SimG4CoreGenerator")
202  << "Skip GenParticle barcode = " << (*pitr)->barcode()
203  << " PDG Id = " << (*pitr)->pdg_id();
204  continue;
205  }
206  }
207 
208  double x2 = 0.0;
209  double y2 = 0.0;
210  double z2 = 0.0;
211  double decay_length = 0.0;
212  int status = (*pitr)->status();
213 
214  // treat status 104 (metastable R-hadrons in Pythia8) the same as status 2
215  if (104 == status)
216  status = 2;
217 
218  // check the status, 2 has end point with decay defined by generator
219  if (1 == status || 2 == status) {
220 
221  // this particle has decayed but have no vertex
222  // it is an exotic case
223  if ( !(*pitr)->end_vertex() ) {
224  status = 1;
225  } else {
226  x2 = (*pitr)->end_vertex()->position().x();
227  y2 = (*pitr)->end_vertex()->position().y();
228  z2 = (*pitr)->end_vertex()->position().z();
229  decay_length =
230  std::sqrt((x2-x1)*(x2-x1)+(y2-y1)*(y2-y1)+(z2-z1)*(z2-z1));
231  }
232  }
233 
234  bool toBeAdded = false;
235  double px = (*pitr)->momentum().px();
236  double py = (*pitr)->momentum().py();
237  double pz = (*pitr)->momentum().pz();
238  double ptot = std::sqrt(px*px + py*py + pz*pz);
239  math::XYZTLorentzVector p(px, py, pz, (*pitr)->momentum().e());
240 
241  double ximpact = x1;
242  double yimpact = y1;
243  double zimpact = z1;
244 
245  // protection against numerical problems for extremely low momenta
246  // compute impact point at transition to Hector
247  const double minTan = 1.e-20;
248  if(fabs(z1) < Z_hector && fabs(pz) >= minTan*ptot) {
249  if(pz > 0.0) { zimpact = Z_hector; }
250  else { zimpact = -Z_hector; }
251  double del = (zimpact - z1)/pz;
252  ximpact += del*px;
253  yimpact += del*py;
254  }
255  double rimpact2 = ximpact*ximpact + yimpact*yimpact;
256 
257  if (verbose > 2) LogDebug("SimG4CoreGenerator")
258  << "Processing GenParticle barcode= " << (*pitr)->barcode()
259  << " pdg= " << (*pitr)-> pdg_id()
260  << " status= " << (*pitr)->status()
261  << " st= " << status
262  << " rimpact(cm)= " << std::sqrt(rimpact2)/cm
263  << " zimpact(cm)= " << zimpact/cm
264  << " ptot(GeV)= " << ptot
265  << " pz(GeV)= " << pz;
266 
267  // Particles of status 1 trasnported along the beam pipe for forward
268  // detectors (HECTOR) always pass to Geant4 without cuts
269  if( 1 == status &&
270  fabs(zimpact) >= Z_hector && rimpact2 <= theDecRCut2) {
271  toBeAdded = true;
272  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator")
273  << "GenParticle barcode = " << (*pitr)->barcode()
274  << " passed case 3";
275  } else {
276 
277  // Standard case: particles not decayed by the generator
278  if(1 == status &&
279  (fabs(zimpact) < Z_hector || rimpact2 > theDecRCut2)) {
280 
281  // Ptot cut for all particles
282  if (fPCuts && (ptot < theMinPCut || ptot > theMaxPCut)) {
283  continue;
284  }
285  // phi cut is applied mainly for particle gun
286  if (fPhiCuts) {
287  double phi = p.phi();
288  if(phi < theMinPhiCut || phi > theMaxPhiCut) {
289  continue;
290  }
291  }
292  // eta cut is applied if position of the decay
293  // is within vacuum chamber and limited in Z
294  if(fEtaCuts) {
295 
296  // eta cut
297  double xi = x1;
298  double yi = y1;
299  double zi = z1;
300 
301  // can be propagated along Z
302  if(fabs(pz) >= minTan*ptot) {
303  if((zi >= Z_lmax) & (pz < 0.0)) {
304  zi = Z_lmax;
305  } else if((zi <= Z_lmin) & (pz > 0.0)) {
306  zi = Z_lmin;
307  } else {
308  if(pz > 0) { zi = Z_lmax; }
309  else { zi = Z_lmin; }
310  }
311  double del = (zi - z1)/pz;
312  xi += del*px;
313  yi += del*py;
314  }
315  // check eta cut
316  if((zi >= Z_lmin) & (zi <= Z_lmax)
317  & (xi*xi + yi*yi < theDecRCut2)) {
318  continue;
319  }
320  }
322  { continue; }
323  toBeAdded = true;
324  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator")
325  << "GenParticle barcode = " << (*pitr)->barcode()
326  << " passed case 1";
327 
328  // Decay chain outside the fiducial cylinder defined by theRDecLenCut
329  // are used for Geant4 tracking with predefined decay channel
330  // In the case of decay in vacuum particle is not tracked by Geant4
331  } else if(2 == status &&
332  x2*x2 + y2*y2 >= theDecRCut2 && fabs(z2) < Z_hector) {
333 
334  toBeAdded = true;
335  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator")
336  << "GenParticle barcode = " << (*pitr)->barcode()
337  << " passed case 2"
338  << " decay_length(cm)= " << decay_length/cm;
339  }
340  }
341  if(toBeAdded){
342  G4int pdgcode= (*pitr)-> pdg_id();
343  G4PrimaryParticle* g4prim=
344  new G4PrimaryParticle(pdgcode, px*GeV, py*GeV, pz*GeV);
345 
346  if ( g4prim->GetG4code() != 0 ){
347  g4prim->SetMass( g4prim->GetG4code()->GetPDGMass() );
348  double charge = g4prim->GetG4code()->GetPDGCharge();
349 
350  // apply Pt cut
351  if (fPtransCut &&
352  0.0 != charge && px*px + py*py < theMinPtCut2) {
353  delete g4prim;
354  continue;
355  }
356  g4prim->SetCharge(charge);
357  }
358 
359  // V.I. do not use SetWeight but the same code
360  // value of the code compute inside TrackWithHistory
361  //g4prim->SetWeight( 10000*(*vpitr)->barcode() ) ;
362  setGenId( g4prim, (*pitr)->barcode() );
363 
364  if (2 == status) {
366  (HepMC::GenParticle *) *pitr, decay_length);
367  }
368  if ( verbose > 1 ) g4prim->Print();
369  g4vtx->SetPrimary(g4prim);
370  }
371  }
372 
373  if ( verbose > 1 ) g4vtx->Print();
374  g4evt->AddPrimaryVertex(g4vtx);
375  ++ng4vtx;
376  }
377 
378  // Add a protection for completely empty events (produced by LHCTransport):
379  // add a dummy vertex with no particle attached to it
380  if ( ng4vtx == 0 ) {
381  G4PrimaryVertex* g4vtx = new G4PrimaryVertex(0.0, 0.0, 0.0, 0.0);
382  if ( verbose > 1 ) g4vtx->Print();
383  g4evt->AddPrimaryVertex(g4vtx);
384  }
385 
386  delete evt;
387 }
#define LogDebug(id)
double theDecRCut2
Definition: Generator.h:53
bool fPhiCuts
Definition: Generator.h:45
void setGenId(G4PrimaryParticle *p, int id) const
Definition: Generator.h:38
const double GeV
Definition: MathUtil.h:16
double theMinPtCut2
Definition: Generator.h:51
double Z_lmax
Definition: Generator.h:61
double theMaxPhiCut
Definition: Generator.h:47
bool fPtransCut
Definition: Generator.h:43
bool pdgFilterSel
Definition: Generator.h:63
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
void particleAssignDaughters(G4PrimaryParticle *p, HepMC::GenParticle *hp, double length)
Definition: Generator.cc:389
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
LumiMonitorFilter * fLumiFilter
Definition: Generator.h:57
T sqrt(T t)
Definition: SSEVec.h:48
bool fPDGFilter
Definition: Generator.h:64
double weight_
Definition: Generator.h:60
double theMaxPCut
Definition: Generator.h:52
bool isGoodForLumiMonitor(const HepMC::GenParticle *) const
math::XYZTLorentzVector * vtx_
Definition: Generator.h:59
double Z_lmin
Definition: Generator.h:61
std::vector< int > pdgFilter
Definition: Generator.h:62
tuple status
Definition: ntuplemaker.py:245
bool fEtaCuts
Definition: Generator.h:44
double Z_hector
Definition: Generator.h:61
bool fPCuts
Definition: Generator.h:42
void Generator::nonBeamEvent2G4 ( const HepMC::GenEvent *  g,
G4Event *  e 
)

Definition at line 500 of file Generator.cc.

References GenParticle::GenParticle, GeV, i, particlePassesPrimaryCuts(), setGenId(), and findQualityFiles::v.

Referenced by RunManagerMTWorker::generateEvent(), and RunManager::generateEvent().

501 {
502  int i = 0;
503  for(HepMC::GenEvent::particle_const_iterator it = evt->particles_begin();
504  it != evt->particles_end(); ++it ) {
505  ++i;
506  HepMC::GenParticle * gp = (*it);
507  int g_status = gp->status();
508  // storing only particle with status == 1
509  if (g_status == 1) {
510  int g_id = gp->pdg_id();
511  G4PrimaryParticle * g4p =
512  new G4PrimaryParticle(g_id,gp->momentum().px()*GeV,
513  gp->momentum().py()*GeV,
514  gp->momentum().pz()*GeV);
515  if (g4p->GetG4code() != 0) {
516  g4p->SetMass(g4p->GetG4code()->GetPDGMass());
517  g4p->SetCharge(g4p->GetG4code()->GetPDGCharge());
518  }
519  // V.I. do not use SetWeight but the same code
520  // value of the code compute inside TrackWithHistory
521  //g4p->SetWeight(i*10000);
522  setGenId(g4p,i);
523  if (particlePassesPrimaryCuts(g4p->GetMomentum())) {
524  G4PrimaryVertex * v =
525  new G4PrimaryVertex(gp->production_vertex()->position().x()*mm,
526  gp->production_vertex()->position().y()*mm,
527  gp->production_vertex()->position().z()*mm,
528  gp->production_vertex()->position().t()*mm/c_light);
529  v->SetPrimary(g4p);
530  g4evt->AddPrimaryVertex(v);
531  if(verbose > 0) { v->Print(); }
532 
533  } else {
534  delete g4p;
535  }
536  }
537  } // end loop on HepMC particles
538 }
int i
Definition: DBlmapReader.cc:9
void setGenId(G4PrimaryParticle *p, int id) const
Definition: Generator.h:38
const double GeV
Definition: MathUtil.h:16
bool particlePassesPrimaryCuts(const G4ThreeVector &p) const
Definition: Generator.cc:467
void Generator::particleAssignDaughters ( G4PrimaryParticle *  p,
HepMC::GenParticle *  hp,
double  length 
)
private

Definition at line 389 of file Generator.cc.

References createTree::dd, GeV, LogDebug, AlCaHLTBitMon_ParallelJobs::p, setGenId(), and mathSSE::sqrt().

Referenced by HepMC2G4().

392 {
393  if ( verbose > 1 ) {
394  LogDebug("SimG4CoreGenerator")
395  << "Special case of long decay length \n"
396  << "Assign daughters with to mother with decaylength="
397  << decaylength/cm << " cm";
398  }
399  math::XYZTLorentzVector p(vp->momentum().px(), vp->momentum().py(),
400  vp->momentum().pz(), vp->momentum().e());
401 
402  // defined preassigned decay time
403  double proper_time = decaylength/(p.Beta()*p.Gamma()*c_light);
404  g4p->SetProperTime(proper_time);
405 
406  if( verbose > 2 ) {
407  LogDebug("SimG4CoreGenerator") <<" px= "<< p.px()
408  <<" py= "<< p.py()
409  <<" pz= "<< p.pz()
410  <<" e= " << p.e()
411  <<" beta= "<< p.Beta()
412  <<" gamma= " << p.Gamma()
413  <<" Proper time= " <<proper_time/ns <<" ns";
414  }
415 
416  // the particle will decay after the same length if it
417  // has not interacted before
418  double x1 = vp->end_vertex()->position().x();
419  double y1 = vp->end_vertex()->position().y();
420  double z1 = vp->end_vertex()->position().z();
421 
422  for (HepMC::GenVertex::particle_iterator
423  vpdec= vp->end_vertex()->particles_begin(HepMC::children);
424  vpdec != vp->end_vertex()->particles_end(HepMC::children); ++vpdec) {
425 
426  //transform decay products such that in the rest frame of mother
427  math::XYZTLorentzVector pdec((*vpdec)->momentum().px(),
428  (*vpdec)->momentum().py(),
429  (*vpdec)->momentum().pz(),
430  (*vpdec)->momentum().e());
431 
432  // children should only be taken into account once
433  G4PrimaryParticle * g4daught=
434  new G4PrimaryParticle((*vpdec)->pdg_id(), pdec.x()*GeV,
435  pdec.y()*GeV, pdec.z()*GeV);
436 
437  if ( g4daught->GetG4code() != 0 )
438  {
439  g4daught->SetMass( g4daught->GetG4code()->GetPDGMass() ) ;
440  g4daught->SetCharge( g4daught->GetG4code()->GetPDGCharge() ) ;
441  }
442 
443  // V.I. do not use SetWeight but the same code
444  // value of the code compute inside TrackWithHistory
445  //g4daught->SetWeight( 10000*(*vpdec)->barcode() ) ;
446  setGenId( g4daught, (*vpdec)->barcode() );
447 
448  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator")
449  <<"Assigning a "<< (*vpdec)->pdg_id()
450  <<" as daughter of a " << vp->pdg_id();
451  if ( (*vpdec)->status() == 2 && (*vpdec)->end_vertex() != 0 )
452  {
453  double x2 = (*vpdec)->end_vertex()->position().x();
454  double y2 = (*vpdec)->end_vertex()->position().y();
455  double z2 = (*vpdec)->end_vertex()->position().z();
456  double dd = std::sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2)+(z1-z2)*(z1-z2));
457  particleAssignDaughters(g4daught,*vpdec,dd);
458  }
459  (*vpdec)->set_status(1000+(*vpdec)->status());
460  g4p->SetDaughter(g4daught);
461 
462  if ( verbose > 1 ) g4daught->Print();
463  }
464 }
#define LogDebug(id)
void setGenId(G4PrimaryParticle *p, int id) const
Definition: Generator.h:38
const double GeV
Definition: MathUtil.h:16
void particleAssignDaughters(G4PrimaryParticle *p, HepMC::GenParticle *hp, double length)
Definition: Generator.cc:389
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
T sqrt(T t)
Definition: SSEVec.h:48
bool Generator::particlePassesPrimaryCuts ( const G4ThreeVector &  p) const
private

Definition at line 467 of file Generator.cc.

References alignCSCRings::e, eta, fEtaCuts, fPCuts, fPhiCuts, GeV, LogDebug, phi, theMaxEtaCut, theMaxPCut, and theMaxPhiCut.

Referenced by nonBeamEvent2G4().

468 {
469  bool flag = true;
470  double ptot = p.mag();
471  if (fPCuts && (ptot < theMinPCut*GeV || ptot > theMaxPCut*GeV)) {
472  flag = false;
473  }
474  if (fEtaCuts && flag) {
475  double pz = p.z();
476  if(ptot < pz + 1.e-10) {
477  flag = false;
478 
479  } else {
480  double eta = 0.5*G4Log((ptot + pz)/(ptot - pz));
481  if(eta < theMinEtaCut || eta > theMaxEtaCut) {
482  flag = false;
483  }
484  }
485  }
486  if (fPhiCuts && flag) {
487  double phi = p.phi();
488  if(phi < theMinPhiCut || phi > theMaxPhiCut) {
489  flag = false;
490  }
491  }
492 
493  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator")
494  << "Generator ptot(GeV)= " << ptot/GeV << " eta= " << p.eta()
495  << " phi= " << p.phi() << " Flag= " << flag;
496 
497  return flag;
498 }
#define LogDebug(id)
bool fPhiCuts
Definition: Generator.h:45
const double GeV
Definition: MathUtil.h:16
double theMaxPhiCut
Definition: Generator.h:47
double theMaxEtaCut
Definition: Generator.h:49
double theMaxPCut
Definition: Generator.h:52
bool fEtaCuts
Definition: Generator.h:44
bool fPCuts
Definition: Generator.h:42
void Generator::setGenEvent ( const HepMC::GenEvent *  inpevt)
inline

Definition at line 25 of file Generator.h.

References evt_, and mathSSE::return().

Referenced by RunManagerMTWorker::generateEvent(), and RunManager::generateEvent().

26  { evt_ = (HepMC::GenEvent*)inpevt; return ; }
return((rh^lh)&mask)
HepMC::GenEvent * evt_
Definition: Generator.h:58
void Generator::setGenId ( G4PrimaryParticle *  p,
int  id 
) const
inlineprivate

Definition at line 38 of file Generator.h.

Referenced by HepMC2G4(), nonBeamEvent2G4(), and particleAssignDaughters().

39  { p->SetUserInformation(new GenParticleInfo(id));}

Member Data Documentation

HepMC::GenEvent* Generator::evt_
private

Definition at line 58 of file Generator.h.

Referenced by genEvent(), and setGenEvent().

bool Generator::fEtaCuts
private

Definition at line 44 of file Generator.h.

Referenced by Generator(), HepMC2G4(), and particlePassesPrimaryCuts().

LumiMonitorFilter* Generator::fLumiFilter
private

Definition at line 57 of file Generator.h.

Referenced by Generator(), HepMC2G4(), and ~Generator().

bool Generator::fPCuts
private

Definition at line 42 of file Generator.h.

Referenced by Generator(), HepMC2G4(), and particlePassesPrimaryCuts().

bool Generator::fPDGFilter
private

Definition at line 64 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

bool Generator::fPhiCuts
private

Definition at line 45 of file Generator.h.

Referenced by Generator(), HepMC2G4(), and particlePassesPrimaryCuts().

bool Generator::fPtransCut
private

Definition at line 43 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

std::vector<int> Generator::pdgFilter
private

Definition at line 62 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

bool Generator::pdgFilterSel
private

Definition at line 63 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

double Generator::theDecLenCut
private

Definition at line 55 of file Generator.h.

Referenced by Generator().

double Generator::theDecRCut2
private

Definition at line 53 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

double Generator::theEtaCutForHector
private

Definition at line 54 of file Generator.h.

Referenced by Generator().

double Generator::theMaxEtaCut
private

Definition at line 49 of file Generator.h.

Referenced by Generator(), and particlePassesPrimaryCuts().

double Generator::theMaxPCut
private

Definition at line 52 of file Generator.h.

Referenced by HepMC2G4(), and particlePassesPrimaryCuts().

double Generator::theMaxPhiCut
private

Definition at line 47 of file Generator.h.

Referenced by HepMC2G4(), and particlePassesPrimaryCuts().

double Generator::theMinEtaCut
private

Definition at line 48 of file Generator.h.

Referenced by Generator().

double Generator::theMinPCut
private

Definition at line 50 of file Generator.h.

Referenced by Generator().

double Generator::theMinPhiCut
private

Definition at line 46 of file Generator.h.

double Generator::theMinPtCut2
private

Definition at line 51 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

int Generator::verbose
private

Definition at line 56 of file Generator.h.

math::XYZTLorentzVector* Generator::vtx_
private

Definition at line 59 of file Generator.h.

Referenced by genVertex(), and HepMC2G4().

double Generator::weight_
private

Definition at line 60 of file Generator.h.

Referenced by eventWeight(), and HepMC2G4().

double Generator::Z_hector
private

Definition at line 61 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

double Generator::Z_lmax
private

Definition at line 61 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

double Generator::Z_lmin
private

Definition at line 61 of file Generator.h.

Referenced by Generator(), and HepMC2G4().