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 G4PrimaryParticle *p) const
 
bool particlePassesPrimaryCuts (const math::XYZTLorentzVector &mom, const double zimp) const
 
void setGenId (G4PrimaryParticle *p, int id) const
 

Private Attributes

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

Detailed Description

Definition at line 21 of file Generator.h.

Constructor & Destructor Documentation

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

Definition at line 22 of file Generator.cc.

References edm::ParameterSet::exists(), create_public_lumi_plots::exp, fEtaCuts, edm::ParameterSet::getParameter(), cuy::ii, LogDebug, pdgFilter, pdgFilterSel, theEtaCutForHector, theMaxEtaCut, theMinEtaCut, theRDecLenCut, Z_hector, Z_lmax, and Z_lmin.

22  :
23  fPCuts(p.getParameter<bool>("ApplyPCuts")),
24  fEtaCuts(p.getParameter<bool>("ApplyEtaCuts")),
25  fPhiCuts(p.getParameter<bool>("ApplyPhiCuts")),
26  theMinPhiCut(p.getParameter<double>("MinPhiCut")), // now operates in radians (CMS standard)
27  theMaxPhiCut(p.getParameter<double>("MaxPhiCut")),
28  theMinEtaCut(p.getParameter<double>("MinEtaCut")),
29  theMaxEtaCut(p.getParameter<double>("MaxEtaCut")),
30  theMinPCut(p.getParameter<double>("MinPCut")), // now operates in GeV (CMS standard)
31  theMaxPCut(p.getParameter<double>("MaxPCut")),
32  theRDecLenCut(p.getParameter<double>("RDecLenCut")*cm),
33  theEtaCutForHector(p.getParameter<double>("EtaCutForHector")),
34  verbose(p.getUntrackedParameter<int>("Verbosity",0)),
35  evt_(0),
36  vtx_(0),
37  weight_(0),
38  Z_lmin(0),
39  Z_lmax(0),
40  Z_hector(0),
41  pdgFilterSel(false)
42 {
43 
44  pdgFilter.resize(0);
45  if ( p.exists("PDGselection") ) {
46 
47  pdgFilterSel = (p.getParameter<edm::ParameterSet>("PDGselection")).getParameter<bool>("PDGfilterSel");
48  pdgFilter = (p.getParameter<edm::ParameterSet>("PDGselection")).getParameter<std::vector< int > >("PDGfilter");
49  for ( unsigned int ii = 0; ii < pdgFilter.size() ; ii++ ) {
50  if (pdgFilterSel) {
51  edm::LogWarning("SimG4CoreGenerator") << " *** Selecting only PDG ID = " << pdgFilter[ii];
52  } else {
53  edm::LogWarning("SimG4CoreGenerator") << " *** Filtering out PDG ID = " << pdgFilter[ii];
54  }
55  }
56 
57  }
58 
59  if(fEtaCuts){
60  Z_lmax = theRDecLenCut*( ( 1 - exp(-2*theMaxEtaCut) ) / ( 2*exp(-theMaxEtaCut) ) );
61  Z_lmin = theRDecLenCut*( ( 1 - exp(-2*theMinEtaCut) ) / ( 2*exp(-theMinEtaCut) ) );
62  }
63 
65 
66  if ( verbose > 0 ) LogDebug("SimG4CoreGenerator") << "Z_min = " << Z_lmin << " Z_max = " << Z_lmax << " Z_hector = " << Z_hector;
67 
68 }
#define LogDebug(id)
bool fPhiCuts
Definition: Generator.h:43
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
double theEtaCutForHector
Definition: Generator.h:51
double Z_lmax
Definition: Generator.h:56
double theMaxPhiCut
Definition: Generator.h:45
bool exists(std::string const &parameterName) const
checks if a parameter exists
double theMinEtaCut
Definition: Generator.h:46
bool pdgFilterSel
Definition: Generator.h:58
int ii
Definition: cuy.py:588
double theMinPCut
Definition: Generator.h:48
double theMinPhiCut
Definition: Generator.h:44
double weight_
Definition: Generator.h:55
double theMaxEtaCut
Definition: Generator.h:47
double theRDecLenCut
Definition: Generator.h:50
double theMaxPCut
Definition: Generator.h:49
math::XYZTLorentzVector * vtx_
Definition: Generator.h:54
HepMC::GenEvent * evt_
Definition: Generator.h:53
double Z_lmin
Definition: Generator.h:56
int verbose
Definition: Generator.h:52
std::vector< int > pdgFilter
Definition: Generator.h:57
bool fEtaCuts
Definition: Generator.h:42
double Z_hector
Definition: Generator.h:56
bool fPCuts
Definition: Generator.h:41
Generator::~Generator ( )
virtual

Definition at line 70 of file Generator.cc.

71 {
72 }

Member Function Documentation

virtual const double Generator::eventWeight ( ) const
inlinevirtual

Definition at line 32 of file Generator.h.

References weight_.

Referenced by RunManager::produce().

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

Definition at line 30 of file Generator.h.

References evt_.

Referenced by RunManager::produce().

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

Definition at line 31 of file Generator.h.

References vtx_.

Referenced by RunManager::produce().

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

Definition at line 74 of file Generator.cc.

References spr::find(), configurableAnalysis::GenParticle, LogDebug, AlCaHLTBitMon_ParallelJobs::p, particleAssignDaughters(), particlePassesPrimaryCuts(), pdgFilter, pdgFilterSel, hitfit::return, setGenId(), mathSSE::sqrt(), theRDecLenCut, vtx_, weight_, and Z_hector.

Referenced by RunManager::generateEvent().

75 {
76 
77  if ( *(evt_orig->vertices_begin()) == 0 )
78  {
79  throw SimG4Exception( "SimG4CoreGenerator: Corrupted Event - GenEvent with no vertex" ) ;
80  }
81 
82 
83  HepMC::GenEvent* evt = new HepMC::GenEvent(*evt_orig) ;
84 
85  if ( evt->weights().size() > 0 )
86  {
87  weight_ = evt->weights()[0] ;
88  for ( unsigned int iw=1; iw<evt->weights().size(); iw++ )
89  {
90  // terminate if the versot of weights contains a zero-weight
91  if ( evt->weights()[iw] <= 0 ) break;
92  weight_ *= evt->weights()[iw] ;
93  }
94  }
95 
96  if (vtx_ != 0) delete vtx_;
97  vtx_ = new math::XYZTLorentzVector((*(evt->vertices_begin()))->position().x(),
98  (*(evt->vertices_begin()))->position().y(),
99  (*(evt->vertices_begin()))->position().z(),
100  (*(evt->vertices_begin()))->position().t());
101 
102  if(verbose >0){
103  evt->print();
104  LogDebug("SimG4CoreGenerator") << "Primary Vertex = (" << vtx_->x() << ","
105  << vtx_->y() << ","
106  << vtx_->z() << ")";
107  }
108 
109  // double x0 = vtx_->x();
110  // double y0 = vtx_->y();
111 
112  unsigned int ng4vtx = 0;
113 
114  for(HepMC::GenEvent::vertex_const_iterator vitr= evt->vertices_begin();
115  vitr != evt->vertices_end(); ++vitr ) {
116  // loop for vertex ...
117  // real vertex?
118  G4bool qvtx=false;
119 
120  for (HepMC::GenVertex::particle_iterator pitr= (*vitr)->particles_begin(HepMC::children);
121  pitr != (*vitr)->particles_end(HepMC::children); ++pitr) {
122  // Admit also status=1 && end_vertex for long vertex special decay treatment
123  if ((*pitr)->status()==1) {
124  qvtx=true;
125  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "GenVertex barcode = " << (*vitr)->barcode()
126  << " selected for GenParticle barcode = " << (*pitr)->barcode() << std::endl;
127  break;
128  }
129  // The selection is made considering if the partcile with status = 2 have the end_vertex
130  // with a radius (R) greater then the theRDecLenCut that means: the end_vertex is outside
131  // the beampipe cilinder (no requirement on the Z of the vertex is applyed).
132  else if ( (*pitr)->status()== 2 ) {
133  if ( (*pitr)->end_vertex() != 0 ) {
134  //double xx = x0-(*pitr)->end_vertex()->position().x();
135  //double yy = y0-(*pitr)->end_vertex()->position().y();
136  double xx = (*pitr)->end_vertex()->position().x();
137  double yy = (*pitr)->end_vertex()->position().y();
138  double r_dd=std::sqrt(xx*xx+yy*yy);
139  if (r_dd>theRDecLenCut){
140  qvtx=true;
141  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "GenVertex barcode = " << (*vitr)->barcode()
142  << " selected for GenParticle barcode = " << (*pitr)->barcode() << " radius = " << r_dd << std::endl;
143  break;
144  }
145  }
146  }
147  }
148 
149 
150  if (!qvtx) {
151  continue;
152  }
153 
154  double x1 = (*vitr)->position().x();
155  double y1 = (*vitr)->position().y();
156  double z1 = (*vitr)->position().z();
157  double t1 = (*vitr)->position().t();
158  G4PrimaryVertex* g4vtx=
159  new G4PrimaryVertex(x1*mm, y1*mm, z1*mm, t1*mm/c_light);
160 
161  for (HepMC::GenVertex::particle_iterator vpitr= (*vitr)->particles_begin(HepMC::children);
162  vpitr != (*vitr)->particles_end(HepMC::children); ++vpitr){
163 
164  // Special cases:
165  // 1) import in Geant4 a full decay chain (e.g. also particles with status == 2)
166  // from the generator in case the decay radius is larger than theRDecLenCut
167  // In this case no cuts will be applied to select the particles hat has to be
168  // processed by geant
169  // 2) import in Geant4 particles with status == 1 but a final end vertex.
170  // The time of the vertex is used as the time of flight to be forced for the particle
171 
172  double r_decay_length=-1;
173  double decay_length=-1;
174 
175  if ( (*vpitr)->status() == 1 || (*vpitr)->status() == 2 ) {
176  // this particle has decayed
177  if ( (*vpitr)->end_vertex() != 0 ) {
178  // needed some particles have status 2 and no end_vertex
179  // Which are the particles with status 2 and not end_vertex, what are suppose to di such kind of particles
180  // when propagated to geant?
181 
182  double x2 = (*vpitr)->end_vertex()->position().x();
183  double y2 = (*vpitr)->end_vertex()->position().y();
184  double z2 = (*vpitr)->end_vertex()->position().z();
185  r_decay_length=std::sqrt(x2*x2+y2*y2);
186  decay_length=std::sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2)+(z1-z2)*(z1-z2));
187 
188  }
189  }
190  // end modification
191 
192  bool toBeAdded = false;
193  math::XYZTLorentzVector p((*vpitr)->momentum().px(),
194  (*vpitr)->momentum().py(),
195  (*vpitr)->momentum().pz(),
196  (*vpitr)->momentum().e());
197 
198  // protection against numerical problems for extremely low momenta
199  const double minTan = 1.e-20;
200  double zimpact = Z_hector+1.;
201  if ( fabs(z1) < Z_hector && fabs(p.pt()/p.pz()) >= minTan ) {
202  // write tan(p.Theta()) as p.Pt()/p.Pz()
203  zimpact = (theRDecLenCut-sqrt(x1*x1+y1*y1))*(p.pz()/p.pt())+z1;
204  }
205 
206  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "Processing GenParticle barcode = " << (*vpitr)->barcode()
207  << " status = " << (*vpitr)->status()
208  << " zimpact = " << zimpact;
209 
210  // Filter on allowed particle species if required
211 
212  if ( pdgFilter.size() > 0 ) {
213  std::vector<int>::iterator it = find(pdgFilter.begin(),pdgFilter.end(),(*vpitr)->pdg_id());
214  if ( (it != pdgFilter.end() && !pdgFilterSel) || ( it == pdgFilter.end() && pdgFilterSel ) ) {
215  toBeAdded = false;
216  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "Skip GenParticle barcode = " << (*vpitr)->barcode()
217  << " PDG Id = " << (*vpitr)->pdg_id() << std::endl;
218  continue;
219  }
220  }
221 
222 
223  // Standard case: particles not decayed by the generator
224  if( (*vpitr)->status() == 1 && fabs(zimpact) < Z_hector ) {
225  if ( !particlePassesPrimaryCuts( p, zimpact ) ) {
226  continue ;
227  }
228  toBeAdded = true;
229  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "GenParticle barcode = " << (*vpitr)->barcode() << " passed case 1" << std::endl;
230  }
231 
232  // Decay chain entering exiting the fiducial cylinder defined by theRDecLenCut
233  else if((*vpitr)->status() == 2 && r_decay_length > theRDecLenCut &&
234  fabs(zimpact) < Z_hector ) {
235  toBeAdded=true;
236  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "GenParticle barcode = " << (*vpitr)->barcode() << " passed case 2" << std::endl;
237  }
238 
239  // Particles trasnported along the beam pipe for forward detectors (HECTOR)
240  // Always pass to Geant4 without cuts (to be checked)
241  else if( (*vpitr)->status() == 1 && fabs(zimpact) >= Z_hector && fabs(z1) >= Z_hector) {
242  toBeAdded = true;
243  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "GenParticle barcode = " << (*vpitr)->barcode() << " passed case 3" << std::endl;
244  }
245 
246  if(toBeAdded){
247 
248  G4int pdgcode= (*vpitr)-> pdg_id();
249  G4PrimaryParticle* g4prim=
250  new G4PrimaryParticle(pdgcode, p.Px()*GeV, p.Py()*GeV, p.Pz()*GeV);
251 
252  if ( g4prim->GetG4code() != 0 ){
253  g4prim->SetMass( g4prim->GetG4code()->GetPDGMass() ) ;
254  g4prim->SetCharge( g4prim->GetG4code()->GetPDGCharge() ) ;
255  }
256 
257  // V.I. do not use SetWeight but the same code
258  // value of the code compute inside TrackWithHistory
259  //g4prim->SetWeight( 10000*(*vpitr)->barcode() ) ;
260  setGenId( g4prim, (*vpitr)->barcode() ) ;
261 
262  if ( (*vpitr)->status() == 2 )
263  particleAssignDaughters(g4prim,(HepMC::GenParticle *) *vpitr, decay_length);
264  if ( verbose > 1 ) g4prim->Print();
265  g4vtx->SetPrimary(g4prim);
266  // impose also proper time for status=1 and available end_vertex
267  if ( (*vpitr)->status()==1 && (*vpitr)->end_vertex()!=0) {
268  double proper_time=decay_length/(p.Beta()*p.Gamma()*c_light);
269  if ( verbose > 1 ) LogDebug("SimG4CoreGenerator") <<"Setting proper time for beta="<<p.Beta()<<" gamma="<<p.Gamma()<<" Proper time=" <<proper_time<<" ns" ;
270  g4prim->SetProperTime(proper_time*ns);
271  }
272  }
273  }
274 
275  if (verbose > 1 ) g4vtx->Print();
276  g4evt->AddPrimaryVertex(g4vtx);
277  ng4vtx++;
278  }
279 
280  // Add a protection for completely empty events (produced by LHCTransport): add a dummy vertex with no particle attached to it
281  if ( ng4vtx == 0 ) {
282  double x1 = 0.;
283  double y1 = 0.;
284  double z1 = 0.;
285  double t1 = 0.;
286  G4PrimaryVertex* g4vtx=
287  new G4PrimaryVertex(x1*mm, y1*mm, z1*mm, t1*mm/c_light);
288  if (verbose > 1 ) g4vtx->Print();
289  g4evt->AddPrimaryVertex(g4vtx);
290  }
291 
292  delete evt ;
293 
294  return ;
295 }
#define LogDebug(id)
void setGenId(G4PrimaryParticle *p, int id) const
Definition: Generator.h:37
bool pdgFilterSel
Definition: Generator.h:58
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:297
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
T sqrt(T t)
Definition: SSEVec.h:48
double weight_
Definition: Generator.h:55
double theRDecLenCut
Definition: Generator.h:50
math::XYZTLorentzVector * vtx_
Definition: Generator.h:54
std::vector< int > pdgFilter
Definition: Generator.h:57
bool particlePassesPrimaryCuts(const G4PrimaryParticle *p) const
Definition: Generator.cc:377
double Z_hector
Definition: Generator.h:56
void Generator::nonBeamEvent2G4 ( const HepMC::GenEvent *  g,
G4Event *  e 
)

Definition at line 401 of file Generator.cc.

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

Referenced by RunManager::generateEvent().

402 {
403  int i = 0;
404  for(HepMC::GenEvent::particle_const_iterator it = evt->particles_begin();
405  it != evt->particles_end(); ++it )
406  {
407  i++;
408  HepMC::GenParticle * g = (*it);
409  int g_status = g->status();
410  // storing only particle with status == 1
411  if (g_status == 1)
412  {
413  int g_id = g->pdg_id();
414  G4PrimaryParticle * g4p =
415  new G4PrimaryParticle(g_id,g->momentum().px()*GeV,g->momentum().py()*GeV,g->momentum().pz()*GeV);
416  if (g4p->GetG4code() != 0)
417  {
418  g4p->SetMass(g4p->GetG4code()->GetPDGMass());
419  g4p->SetCharge(g4p->GetG4code()->GetPDGCharge()) ;
420  }
421  // V.I. do not use SetWeight but the same code
422  // value of the code compute inside TrackWithHistory
423  //g4p->SetWeight(i*10000);
424  setGenId(g4p,i);
425  if (particlePassesPrimaryCuts(g4p))
426  {
427  G4PrimaryVertex * v = new
428  G4PrimaryVertex(g->production_vertex()->position().x()*mm,
429  g->production_vertex()->position().y()*mm,
430  g->production_vertex()->position().z()*mm,
431  g->production_vertex()->position().t()*mm/c_light);
432  v->SetPrimary(g4p);
433  g4evt->AddPrimaryVertex(v);
434  if(verbose >0) {
435  v->Print();
436  }
437  }
438  }
439  } // end loop on HepMC particles
440 }
int i
Definition: DBlmapReader.cc:9
void setGenId(G4PrimaryParticle *p, int id) const
Definition: Generator.h:37
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
bool particlePassesPrimaryCuts(const G4PrimaryParticle *p) const
Definition: Generator.cc:377
void Generator::particleAssignDaughters ( G4PrimaryParticle *  p,
HepMC::GenParticle *  hp,
double  length 
)
private

Definition at line 297 of file Generator.cc.

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

Referenced by HepMC2G4().

298 {
299 
300  if ( !(vp->end_vertex()) ) return ;
301 
302  if ( verbose > 1 )
303  LogDebug("SimG4CoreGenerator") << "Special case of long decay length \n"
304  << "Assign daughters with to mother with decaylength=" << decaylength << "mm";
305  math::XYZTLorentzVector p(vp->momentum().px(), vp->momentum().py(), vp->momentum().pz(), vp->momentum().e());
306  double proper_time=decaylength/(p.Beta()*p.Gamma()*c_light);
307  if ( verbose > 2 ) {
308  LogDebug("SimG4CoreGenerator") <<" px="<<vp->momentum().px()
309  <<" py="<<vp->momentum().py()
310  <<" pz="<<vp->momentum().pz()
311  <<" e="<<vp->momentum().e()
312  <<" beta="<<p.Beta()
313  <<" gamma="<<p.Gamma()
314  <<" Proper time=" <<proper_time<<" ns" ;
315  }
316  g4p->SetProperTime(proper_time*ns); // the particle will decay after the same length if it has not interacted before
317  double x1 = vp->end_vertex()->position().x();
318  double y1 = vp->end_vertex()->position().y();
319  double z1 = vp->end_vertex()->position().z();
320  for (HepMC::GenVertex::particle_iterator
321  vpdec= vp->end_vertex()->particles_begin(HepMC::children);
322  vpdec != vp->end_vertex()->particles_end(HepMC::children); ++vpdec) {
323 
324  //transform decay products such that in the rest frame of mother
325  math::XYZTLorentzVector pdec((*vpdec)->momentum().px(),
326  (*vpdec)->momentum().py(),
327  (*vpdec)->momentum().pz(),
328  (*vpdec)->momentum().e());
329  // children should only be taken into account once
330  G4PrimaryParticle * g4daught=
331  new G4PrimaryParticle((*vpdec)->pdg_id(), pdec.x()*GeV, pdec.y()*GeV, pdec.z()*GeV);
332  if ( g4daught->GetG4code() != 0 )
333  {
334  g4daught->SetMass( g4daught->GetG4code()->GetPDGMass() ) ;
335  g4daught->SetCharge( g4daught->GetG4code()->GetPDGCharge() ) ;
336  }
337  // V.I. do not use SetWeight but the same code
338  // value of the code compute inside TrackWithHistory
339  //g4daught->SetWeight( 10000*(*vpdec)->barcode() ) ;
340  setGenId( g4daught, (*vpdec)->barcode() ) ;
341  if ( verbose > 1 ) LogDebug("SimG4CoreGenerator") <<"Assigning a "<<(*vpdec)->pdg_id()
342  <<" as daughter of a " <<vp->pdg_id() ;
343  if ( (*vpdec)->status() == 2 && (*vpdec)->end_vertex() != 0 )
344  {
345  double x2 = (*vpdec)->end_vertex()->position().x();
346  double y2 = (*vpdec)->end_vertex()->position().y();
347  double z2 = (*vpdec)->end_vertex()->position().z();
348  double dd = std::sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2)+(z1-z2)*(z1-z2));
349  particleAssignDaughters(g4daught,*vpdec,dd);
350  }
351  (*vpdec)->set_status(1000+(*vpdec)->status());
352  g4p->SetDaughter(g4daught);
353  }
354  return;
355 }
#define LogDebug(id)
void setGenId(G4PrimaryParticle *p, int id) const
Definition: Generator.h:37
void particleAssignDaughters(G4PrimaryParticle *p, HepMC::GenParticle *hp, double length)
Definition: Generator.cc:297
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
T sqrt(T t)
Definition: SSEVec.h:48
bool Generator::particlePassesPrimaryCuts ( const G4PrimaryParticle *  p) const
private

Definition at line 377 of file Generator.cc.

References eta(), fEtaCuts, fPCuts, fPhiCuts, create_public_lumi_plots::log, LogDebug, phi, mathSSE::sqrt(), funct::tan(), theMaxEtaCut, theMaxPCut, and theMaxPhiCut.

Referenced by HepMC2G4(), and nonBeamEvent2G4().

378 {
379 
380  G4ThreeVector mom = p->GetMomentum();
381  double phi = mom.phi() ;
382  double nrg = sqrt(p->GetPx()*p->GetPx() + p->GetPy()*p->GetPy() + p->GetPz()*p->GetPz());
383  nrg /= GeV ; // need to convert, since Geant4 operates in MeV
384  double eta = -log(tan(mom.theta()/2));
385  bool flag = true;
386 
387  if ( (fEtaCuts) && (eta < theMinEtaCut || eta > theMaxEtaCut) ) {
388  flag = false;
389  } else if ( (fPCuts) && (nrg < theMinPCut || nrg > theMaxPCut) ) {
390  flag = false;
391  } else if ( (fPhiCuts) && (phi < theMinPhiCut || phi > theMaxPhiCut) ) {
392  flag = false;
393  }
394 
395  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "Generator p = " << nrg << " eta = " << eta << " theta = " << mom.theta() << " phi = " << phi << " Flag = " << flag;
396 
397  return flag;
398 
399 }
#define LogDebug(id)
bool fPhiCuts
Definition: Generator.h:43
long int flag
Definition: mlp_lapack.h:47
double theMaxPhiCut
Definition: Generator.h:45
T eta() const
T sqrt(T t)
Definition: SSEVec.h:48
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
double theMaxEtaCut
Definition: Generator.h:47
double theMaxPCut
Definition: Generator.h:49
bool fEtaCuts
Definition: Generator.h:42
Definition: DDAxes.h:10
bool fPCuts
Definition: Generator.h:41
bool Generator::particlePassesPrimaryCuts ( const math::XYZTLorentzVector mom,
const double  zimp 
) const
private

Definition at line 357 of file Generator.cc.

References fEtaCuts, fPCuts, fPhiCuts, LogDebug, phi, theMaxPCut, theMaxPhiCut, and Z_lmax.

358 {
359 
360  double phi = mom.Phi() ;
361  double nrg = mom.P() ;
362  bool flag = true;
363 
364  if ( (fEtaCuts) && (zimp < Z_lmin || zimp > Z_lmax) ) {
365  flag = false;
366  } else if ( (fPCuts) && (nrg < theMinPCut || nrg > theMaxPCut) ) {
367  flag = false;
368  } else if ( (fPhiCuts) && (phi < theMinPhiCut || phi > theMaxPhiCut) ) {
369  flag = false;
370  }
371 
372  if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "Generator p = " << nrg << " z_imp = " << zimp << " phi = " << mom.Phi() << " Flag = " << flag;
373 
374  return flag;
375 }
#define LogDebug(id)
bool fPhiCuts
Definition: Generator.h:43
long int flag
Definition: mlp_lapack.h:47
double Z_lmax
Definition: Generator.h:56
double theMaxPhiCut
Definition: Generator.h:45
double theMaxPCut
Definition: Generator.h:49
bool fEtaCuts
Definition: Generator.h:42
Definition: DDAxes.h:10
bool fPCuts
Definition: Generator.h:41
void Generator::setGenEvent ( const HepMC::GenEvent *  inpevt)
inline

Definition at line 27 of file Generator.h.

References evt_, and hitfit::return.

Referenced by RunManager::generateEvent().

27 { evt_ = (HepMC::GenEvent*)inpevt; return ; }
HepMC::GenEvent * evt_
Definition: Generator.h:53
void Generator::setGenId ( G4PrimaryParticle *  p,
int  id 
) const
inlineprivate

Definition at line 37 of file Generator.h.

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

Member Data Documentation

HepMC::GenEvent* Generator::evt_
private

Definition at line 53 of file Generator.h.

Referenced by genEvent(), and setGenEvent().

bool Generator::fEtaCuts
private

Definition at line 42 of file Generator.h.

Referenced by Generator(), and particlePassesPrimaryCuts().

bool Generator::fPCuts
private

Definition at line 41 of file Generator.h.

Referenced by particlePassesPrimaryCuts().

bool Generator::fPhiCuts
private

Definition at line 43 of file Generator.h.

Referenced by particlePassesPrimaryCuts().

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

Definition at line 57 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

bool Generator::pdgFilterSel
private

Definition at line 58 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

double Generator::theEtaCutForHector
private

Definition at line 51 of file Generator.h.

Referenced by Generator().

double Generator::theMaxEtaCut
private

Definition at line 47 of file Generator.h.

Referenced by Generator(), and particlePassesPrimaryCuts().

double Generator::theMaxPCut
private

Definition at line 49 of file Generator.h.

Referenced by particlePassesPrimaryCuts().

double Generator::theMaxPhiCut
private

Definition at line 45 of file Generator.h.

Referenced by particlePassesPrimaryCuts().

double Generator::theMinEtaCut
private

Definition at line 46 of file Generator.h.

Referenced by Generator().

double Generator::theMinPCut
private

Definition at line 48 of file Generator.h.

double Generator::theMinPhiCut
private

Definition at line 44 of file Generator.h.

double Generator::theRDecLenCut
private

Definition at line 50 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

int Generator::verbose
private

Definition at line 52 of file Generator.h.

math::XYZTLorentzVector* Generator::vtx_
private

Definition at line 54 of file Generator.h.

Referenced by genVertex(), and HepMC2G4().

double Generator::weight_
private

Definition at line 55 of file Generator.h.

Referenced by eventWeight(), and HepMC2G4().

double Generator::Z_hector
private

Definition at line 56 of file Generator.h.

Referenced by Generator(), and HepMC2G4().

double Generator::Z_lmax
private

Definition at line 56 of file Generator.h.

Referenced by Generator(), and particlePassesPrimaryCuts().

double Generator::Z_lmin
private

Definition at line 56 of file Generator.h.

Referenced by Generator().