BOSS 7.1.1
BESIII Offline Software System
Loading...
Searching...
No Matches
G4HepMCInterface Class Reference

#include <G4HepMCInterface.h>

+ Inheritance diagram for G4HepMCInterface:

Public Member Functions

void HepMC2G4 (HepMC::GenEvent *hepmcevt, G4Event *g4event)
 
void Print (const HepMC::GenEvent *hepmcevt)
 
G4int CheckType (const HepMC::GenEvent *hepmcevt)
 
void Boost (HepMC::GenEvent *hepmcevt)
 
void PrintHPlist ()
 
G4int GetLogLevel ()
 
void SetLogLevel (G4int level)
 
 G4HepMCInterface ()
 
virtual ~G4HepMCInterface ()
 
HepMC::GenEvent * GetHepMCGenEvent () const
 
virtual void GeneratePrimaryVertex (G4Event *anEvent)
 

Public Attributes

std::vector< G4HepMCParticle * > HPlist
 

Protected Member Functions

virtual HepMC::GenEvent * GenerateHepMCEvent ()
 

Protected Attributes

HepMC::GenEvent * hepmcEvent
 
G4int m_logLevel
 

Detailed Description

Definition at line 45 of file G4HepMCInterface.h.

Constructor & Destructor Documentation

◆ G4HepMCInterface()

G4HepMCInterface::G4HepMCInterface ( )

Definition at line 66 of file G4HepMCInterface.cpp.

67 : hepmcEvent(0), m_geantinoPdgId(100), m_chargedgeantinoPdgId(101)
68////////////////////////////////////
69{
70
71 IRealizationSvc *tmpReal;
72 ISvcLocator* svcLocator = Gaudi::svcLocator();
73 StatusCode status = svcLocator->service("RealizationSvc",tmpReal);
74 if (!status.isSuccess())
75 {
76 std::cout << " Could not initialize Realization Service in G4HepMCInterface" << std::endl;
77 } else {
78 m_RealizationSvc=dynamic_cast<RealizationSvc*>(tmpReal);
79 }
80
81 int pdlId1=-1, pdlId2=-1;
82 auto algMan = svcLocator->as<IAlgManager>();
83 if ( algMan ) {
84 if(algMan->existsAlgorithm("EvtDecay")){ // BesEvtGen used
85 EvtId evtId1 = EvtPDL::getId("geantino");
86 pdlId1 = evtId1.getId();
87 if(pdlId1!=-1){ // denote geantino has been defined in pdt.table
88 m_geantinoPdgId = EvtPDL::getStdHep(evtId1);
89 }
90 EvtId evtId2 = EvtPDL::getId("chargedgeantino");
91 pdlId2 = evtId2.getId();
92 if(pdlId2!=-1){ // denote chargedgeantino has been defined in pdt.table
93 m_chargedgeantinoPdgId = EvtPDL::getStdHep(evtId2);
94 }
95 }
96 }
97
98 if(pdlId1==-1||pdlId2==-1){ // generator is not BesEvtGen, defined in PDGTABLE.MeV; if not, will use default
99 IPartPropSvc* p_PartPropSvc;
100 StatusCode PartPropStatus = svcLocator->service("PartPropSvc", p_PartPropSvc);
101 if (!PartPropStatus.isSuccess() || 0 == p_PartPropSvc) {
102 std::cout << " Could not initialize Particle Properties Service in G4HepMCInterface" << std::endl;
103 }
104 else{
105 HepPDT::ParticleDataTable* pdt = p_PartPropSvc->PDT();
106 if(pdlId1==-1){
107 const HepPDT::ParticleData* pd = pdt->particle("geantino");
108 if(pd) m_geantinoPdgId = pd->pid();
109 }
110 if(pdlId2==-1){
111 const HepPDT::ParticleData* pd = pdt->particle("chargedgeantino");
112 if(pd) m_chargedgeantinoPdgId = pd->pid();
113 }
114 }
115 }
116
117}
Definition EvtId.hh:27
int getId() const
Definition EvtId.hh:41
static int getStdHep(EvtId id)
Definition EvtPDL.hh:56
static EvtId getId(const std::string &name)
Definition EvtPDL.cc:287
HepMC::GenEvent * hepmcEvent

◆ ~G4HepMCInterface()

G4HepMCInterface::~G4HepMCInterface ( )
virtual

Definition at line 120 of file G4HepMCInterface.cpp.

122{
123 delete hepmcEvent;
124}

Member Function Documentation

◆ Boost()

void G4HepMCInterface::Boost ( HepMC::GenEvent * hepmcevt)

Definition at line 201 of file G4HepMCInterface.cpp.

203{
204 //suppose relavtive vectoy is v(vx,vy,vz),position are (x,y,z,t)in CMS and(X,Y,Z,T)in Lab,the formulae for boost vertex position from Cms to Lab in natural units are following:
205 // v2 = vx*vx + vy*vy + vz*vz;
206 // gamma = 1.0 / sqrt(1.0 - v2);
207 // bp = vx*x + vy*y + vz*z;
208 // X=x + gamma2*bp*vx + gamma*vx*t;
209 // Y=y + gamma2*bp*vy + gamma*vy*t;
210 // Z=z + gamma2*bp*vz + gamma*vz*t;
211 // T=gamma*(t + bp);
212 // the function of these formulae is the same as xvtx.boost(ThreeVector v)
213 //For the E_cms,need to loop and splice the momentum of all the final state particles
214 G4LorentzVector ptot(0,0,0,0);
215 double E_cms,p_Mag,e_Mass2,M_cms,theta;
216 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
217 vitr != hepmcevt->vertices_end(); ++vitr ) { // loop for vertex ...
218
219 for(HepMC::GenVertex::particle_iterator
220 vpitr = (*vitr)->particles_begin(HepMC::children);
221 vpitr != (*vitr)->particles_end(HepMC::children); ++vpitr) {
222
223 if((*vpitr)->status() !=1)continue;// Not final state particle
224
225 G4LorentzVector p;
226 p.setX( (*vpitr)-> momentum().x());
227 p.setY( (*vpitr)-> momentum().y());
228 p.setZ( (*vpitr)-> momentum().z());
229 p.setT( (*vpitr)-> momentum().t());
230 ptot = p + ptot;
231 }
232 }
233 E_cms=ptot.e()*GeV;
234
235 //get G4Svc, allow user to access the beam momentum shift in JobOptions
236 ISvcLocator* svcLocator = Gaudi::svcLocator();
237 IG4Svc* tmpSvc;
238 G4Svc* m_G4Svc;
239 StatusCode sc=svcLocator->service("G4Svc", tmpSvc);
240 m_G4Svc=dynamic_cast<G4Svc *>(tmpSvc);
241
242
243 G4ThreeVector v(0,0,0); //velocity of cms as seen from lab
244 //for cms velocity
245 theta=11*mrad;
246 //theta=(m_G4Svc->GetBeamAngle())*mrad;
247 //theta is half of the angle between the two beams,namely,is the angle between the beam and Z axis.
248 M_cms=E_cms; //for p_cms=0 in the CMS
249 e_Mass2=electron_mass_c2*electron_mass_c2; //square of electron mass
250 p_Mag=sqrt((E_cms*E_cms-4*e_Mass2)/(2*(1-cos(pi-2*theta))));
251 G4ThreeVector p_Positron(p_Mag*sin(theta),0,p_Mag*cos(theta));
252 G4ThreeVector p_Electron(p_Mag*sin(pi-theta),0,p_Mag*cos(pi-theta));
253 G4ThreeVector vee=p_Positron+p_Electron;
254 G4ThreeVector beamshift(m_G4Svc->GetBeamShiftPx(),m_G4Svc->GetBeamShiftPy(),m_G4Svc->GetBeamShiftPz());
255 if(m_G4Svc-> GetSetBeamShift()==true && fabs(M_cms-3686)<50.0) {vee = beamshift;}
256
257 double pee=vee.r();
258 M_cms = sqrt(M_cms*M_cms + pee*pee);
259 v=vee/M_cms;
260
261
262 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
263 vitr != hepmcevt->vertices_end(); ++vitr ) { // loop for vertex ...
264
265 // Boost vertex position from cms to lab
266 G4LorentzVector xvtx;
267 xvtx.setX( (*vitr)-> position().x());
268 xvtx.setY( (*vitr)-> position().y());
269 xvtx.setZ( (*vitr)-> position().z());
270 xvtx.setT( (*vitr)-> position().t());
271 xvtx.boost(v);
272 (*vitr)->set_position(xvtx);
273
274 // Boost the particles four momentum from cms to lab
275 // the transform formulae are similary to the position boost
276 for (HepMC::GenVertex::particle_iterator
277 vpitr = (*vitr)->particles_begin(HepMC::children);
278 vpitr != (*vitr)->particles_end(HepMC::children); ++vpitr) {
279 G4LorentzVector p;
280 p.setX( (*vpitr)-> momentum().x());
281 p.setY( (*vpitr)-> momentum().y());
282 p.setZ( (*vpitr)-> momentum().z());
283 p.setT( (*vpitr)-> momentum().t());
284 p=p.boost(v);
285 (*vpitr)->set_momentum(p);
286 }
287 }
288}
double sin(const BesAngle a)
Definition BesAngle.h:210
double cos(const BesAngle a)
Definition BesAngle.h:213
**********INTEGER nmxhep !maximum number of particles DOUBLE PRECISION vhep INTEGER jdahep COMMON hepevt $ !serial number $ !number of particles $ !status code $ !particle ident KF $ !parent particles $ !childreen particles $ !four momentum
Double_t x[10]
**********Class see also m_nmax DOUBLE PRECISION m_amel DOUBLE PRECISION m_x2 DOUBLE PRECISION m_alfinv DOUBLE PRECISION m_Xenph INTEGER m_KeyWtm INTEGER m_idyfs DOUBLE PRECISION m_zini DOUBLE PRECISION m_q2 DOUBLE PRECISION m_Wt_KF DOUBLE PRECISION m_WtCut INTEGER m_KFfin *COMMON c_KarLud $ !Input CMS energy[GeV] $ !CMS energy after beam spread beam strahlung[GeV] $ !Beam energy spread[GeV] $ !z boost due to beam spread $ !electron beam mass *ff pair spectrum $ !minimum v
Definition KarLud.h:35
TTree * t
Definition binning.cxx:23
Definition G4Svc.h:33
double GetBeamShiftPz()
Definition G4Svc.h:86
double GetBeamShiftPx()
Definition G4Svc.h:84
double GetBeamShiftPy()
Definition G4Svc.h:85
double y[1000]
const float pi
Definition vector3.h:133

Referenced by HepMC2G4().

◆ CheckType()

G4int G4HepMCInterface::CheckType ( const HepMC::GenEvent * hepmcevt)

Definition at line 176 of file G4HepMCInterface.cpp.

177{
178 G4int flag =0;
179 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
180 vitr != hepmcevt->vertices_end(); ++vitr ) { // loop for vertex ...
181 for (HepMC::GenVertex::particle_iterator
182 pitr= (*vitr)->particles_begin(HepMC::children);
183 pitr != (*vitr)->particles_end(HepMC::children); ++pitr) {
184 if((*pitr)->end_vertex())
185 {flag = 1; break;}
186 }
187 if(flag) break;
188 }
189 if(m_logLevel <= 4)
190 {
191 if(flag == 0)
192 G4cout<<G4endl<<"generator is GENBES!"<<G4endl;
193 else
194 G4cout<<G4endl<<"generator is not GENBES!"<<G4endl;
195 }
196 return flag;
197}

Referenced by HepMC2G4().

◆ GenerateHepMCEvent()

HepMC::GenEvent * G4HepMCInterface::GenerateHepMCEvent ( )
protectedvirtual

Reimplemented in BesHepMCInterface.

Definition at line 558 of file G4HepMCInterface.cpp.

560{
561 HepMC::GenEvent* aevent= new HepMC::GenEvent();
562 return aevent;
563}

Referenced by GeneratePrimaryVertex().

◆ GeneratePrimaryVertex()

void G4HepMCInterface::GeneratePrimaryVertex ( G4Event * anEvent)
virtual

Definition at line 566 of file G4HepMCInterface.cpp.

568{
569 // delete previous event object
570 delete hepmcEvent;
571
572 // generate next event
574 if(! hepmcEvent) {
575 G4Exception("G4HepMCInterface","GeneratePrimaryVertex",RunMustBeAborted,
576 "HepMCInterface: no generated particles. run terminated..." );
577 return;
578 }
579 HepMC2G4(hepmcEvent, anEvent);
580}
void HepMC2G4(HepMC::GenEvent *hepmcevt, G4Event *g4event)
virtual HepMC::GenEvent * GenerateHepMCEvent()

◆ GetHepMCGenEvent()

HepMC::GenEvent * G4HepMCInterface::GetHepMCGenEvent ( ) const
inline

Definition at line 98 of file G4HepMCInterface.h.

99{
100 return hepmcEvent;
101}

◆ GetLogLevel()

G4int G4HepMCInterface::GetLogLevel ( )
inline

Definition at line 55 of file G4HepMCInterface.h.

55{return m_logLevel;}

◆ HepMC2G4()

void G4HepMCInterface::HepMC2G4 ( HepMC::GenEvent * hepmcevt,
G4Event * g4event )

Definition at line 291 of file G4HepMCInterface.cpp.

294{
295 //Print HepMC Event before boost
296 if(m_logLevel <=4 )
297 Print(hepmcevt);
298
299 //get G4Svc
300 ISvcLocator* svcLocator = Gaudi::svcLocator();
301 IG4Svc* tmpSvc;
302 G4Svc* m_G4Svc;
303 StatusCode sc=svcLocator->service("G4Svc", tmpSvc);
304 m_G4Svc=dynamic_cast<G4Svc *>(tmpSvc);
305
306 //considering 22rmad beam clossing angle ,need to boost the Cms to Lab
307 // need to boost ?
308 if(m_G4Svc->GetBoostLab()== true)
309 Boost(hepmcevt);
310
311 //Print HepMC Event after boost
312 if(m_logLevel <=2 && m_G4Svc->GetBoostLab()== true)
313 Print(hepmcevt);
314
315 //*********************send particles to HPlist*************************//
316 G4int flag = CheckType(hepmcevt);
317 if(flag)
318 {
319 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
320 vitr != hepmcevt->vertices_end(); ++vitr )
321 { // loop for vertex ...
322 G4int vtxFlag=0;
323 G4int pOut = (*vitr)->particles_out_size();
324 HepMC::GenVertex::particle_iterator pitr = (*vitr)->particles_begin(HepMC::children);
325 G4int pdgcode= (*pitr)-> pdg_id();
326 if(pOut == 1 && abs(pdgcode) == 11)
327 vtxFlag=1;
328 G4int barcodeVtx = (*vitr)->barcode();
329
330 for (HepMC::GenVertex::particle_iterator pitr= (*vitr)->particles_begin(HepMC::children);
331 pitr != (*vitr)->particles_end(HepMC::children); ++pitr)
332 {// loop for particle from this vertex
333 G4int pdgcode = (*pitr)->pdg_id();
334 if(vtxFlag==0)
335 {
336 if(pdgcode==-22) pdgcode=22;
337 G4LorentzVector p;
338 p.setX( (*pitr)-> momentum().x());
339 p.setY( (*pitr)-> momentum().y());
340 p.setZ( (*pitr)-> momentum().z());
341 p.setT( (*pitr)-> momentum().t());
342 G4PrimaryParticle* particle = 0;
343 // geantino and chargedgeantino can't create through pdgcode;
344 // their PdgID = 0 in Geant4.10.07.p02, and if particle->SetPDGcode(100), Geant4 will not know.
345 // So reset PdgID of input HepMC::GenParticle to G4ParticleDefinition's;
346 // if not, BesSensitiveManager::UpdatePrimaryTrack(const G4Track* track) will exit(1) because MatchError.
347 if(pdgcode==m_geantinoPdgId){
348 G4ParticleDefinition* pd = G4ParticleTable::GetParticleTable()->FindParticle("geantino");
349 particle = new G4PrimaryParticle(pd, p.x()*GeV, p.y()*GeV, p.z()*GeV);
350 (*pitr)->set_pdg_id(particle->GetPDGcode());
351 }
352 else if(pdgcode==m_chargedgeantinoPdgId){
353 G4ParticleDefinition* pd = G4ParticleTable::GetParticleTable()->FindParticle("chargedgeantino");
354 particle = new G4PrimaryParticle(pd, p.x()*GeV, p.y()*GeV, p.z()*GeV);
355 (*pitr)->set_pdg_id(particle->GetPDGcode());
356 }
357 // pdgcode==0 dot define both in Geant4 and pdt.table
358 //else if(pdgcode==0){
359 //}
360 else{
361 particle = new G4PrimaryParticle(pdgcode, p.x()*GeV, p.y()*GeV, p.z()*GeV);
362 }
363 if(!particle){
364 G4cout << "ERROR! Cannot create particle PDGCode=" << pdgcode << G4endl;
365 }
366 //G4PrimaryParticle* particle = new G4PrimaryParticle(pdgcode);
367 //particle->Set4Momentum(p.x()*GeV,p.y()*GeV,p.z()*GeV,p.t()*GeV);
368 G4int ISTHEP = (*pitr)->status();
369
370 G4int barcodeEndVtx = 0;
371 if((*pitr)->end_vertex())
372 barcodeEndVtx = (*pitr)->end_vertex()->barcode();
373
374 G4HepMCParticle* hepmcParticle = new G4HepMCParticle(particle,ISTHEP,barcodeEndVtx);
375 HPlist.push_back(hepmcParticle);
376
377 for( size_t ii=0; ii<HPlist.size(); ii++ )
378 {
379 if(barcodeVtx == HPlist[ii]->GetBarcodeEndVtx())
380 {
381 HPlist[ii]->GetTheParticle()->SetDaughter(particle);
382 hepmcParticle->Done();
383 break;
384 }
385 }
386 }
387 }
388 }
389 }
390 else
391 {
392 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
393 vitr != hepmcevt->vertices_end(); ++vitr )
394 {
395 for (HepMC::GenVertex::particle_iterator pitr= (*vitr)->particles_begin(HepMC::children);
396 pitr != (*vitr)->particles_end(HepMC::children); ++pitr)
397 {
398 G4int ISTHEP = (*pitr)->status();
399 G4LorentzVector p;
400 p.setX( (*pitr)-> momentum().x());
401 p.setY( (*pitr)-> momentum().y());
402 p.setZ( (*pitr)-> momentum().z());
403 p.setT( (*pitr)-> momentum().t());
404
405 G4int pdgcode = (*pitr)->pdg_id();
406 G4int barcodeEndVtx = 0;
407 if((*pitr)->end_vertex())
408 barcodeEndVtx = (*pitr)->end_vertex()->barcode();
409 G4PrimaryParticle* particle = 0;
410 // geantino and chargedgeantino can't create through pdgcode;
411 // their PdgID = 0 in Geant4.10.07.p02, and if particle->SetPDGcode(100), Geant4 will not know.
412 // So reset PdgID of input HepMC::GenParticle to G4ParticleDefinition's;
413 // if not, BesSensitiveManager::UpdatePrimaryTrack(const G4Track* track) will exit(1) because MatchError.
414 if(pdgcode==m_geantinoPdgId){
415 G4ParticleDefinition* pd = G4ParticleTable::GetParticleTable()->FindParticle("geantino");
416 particle = new G4PrimaryParticle(pd, p.x()*GeV, p.y()*GeV, p.z()*GeV);
417 (*pitr)->set_pdg_id(particle->GetPDGcode());
418 }
419 else if(pdgcode==m_chargedgeantinoPdgId){
420 G4ParticleDefinition* pd = G4ParticleTable::GetParticleTable()->FindParticle("chargedgeantino");
421 particle = new G4PrimaryParticle(pd, p.x()*GeV, p.y()*GeV, p.z()*GeV);
422 (*pitr)->set_pdg_id(particle->GetPDGcode());
423 }
424 // pdgcode==0 dot define both in Geant4 and pdt.table
425 //else if(pdgcode==0){
426 //}
427 else{
428 particle = new G4PrimaryParticle(pdgcode, p.x()*GeV, p.y()*GeV, p.z()*GeV);
429 }
430 if(!particle){
431 G4cout << "ERROR! Cannot create particle PDGCode=" << pdgcode << G4endl;
432 }
433 //G4PrimaryParticle* particle = new G4PrimaryParticle(pdgcode);
434 //particle->Set4Momentum(p.x()*GeV,p.y()*GeV,p.z()*GeV,p.t()*GeV);
435
436 G4HepMCParticle* hepmcParticle = new G4HepMCParticle(particle,ISTHEP,barcodeEndVtx);
437 HPlist.push_back(hepmcParticle);
438 if(ISTHEP>1)
439 hepmcParticle->Done();
440 }
441 }
442 }
443
444 //print particles in HPlist
445 if(m_logLevel<=4)
446 PrintHPlist();
447
448 //get time and position information from G4Svc
449 G4double pmPosX,pmPosY,pmPosZ,pmTime;
450 G4double tmpPosX,tmpPosY,tmpPosZ,tmpT;
451 G4double beamPosX,beamPosY,beamPosZ,beamSizeX,beamSizeY,beamSizeZ;
452
453 if(m_RealizationSvc->UseDBFlag() == false) {
454 beamPosX = m_G4Svc->GetBeamPosX();
455 beamPosY = m_G4Svc->GetBeamPosY();
456 beamPosZ = m_G4Svc->GetBeamPosZ();
457
458 beamSizeX = m_G4Svc->GetBeamSizeX();
459 beamSizeY = m_G4Svc->GetBeamSizeY();
460 beamSizeZ = m_G4Svc->GetBeamSizeZ()/sqrt(2);
461 }
462 if(m_RealizationSvc->UseDBFlag() == true && m_RealizationSvc->ifReadBunch() == true) {
463 beamPosX = m_RealizationSvc->getBunchPosX();
464 beamPosY = m_RealizationSvc->getBunchPosY();
465 beamPosZ = m_RealizationSvc->getBunchPosZ();
466
467 beamSizeX = m_RealizationSvc->getBunchSizeX();
468 beamSizeY = m_RealizationSvc->getBunchSizeY();
469 beamSizeZ = m_RealizationSvc->getBunchSizeZ();
470 }
471 if(m_RealizationSvc->UseDBFlag() == true && m_RealizationSvc->ifReadBunch() == false) {
472 beamPosX = m_G4Svc->GetBeamPosX();
473 beamPosY = m_G4Svc->GetBeamPosY();
474 beamPosZ = m_G4Svc->GetBeamPosZ();
475
476 beamSizeX = m_G4Svc->GetBeamSizeX();
477 beamSizeY = m_G4Svc->GetBeamSizeY();
478 beamSizeZ = m_G4Svc->GetBeamSizeZ()/sqrt(2);
479 }
480
481 G4double gaussX = G4RandGauss::shoot();
482 G4double gaussY = G4RandGauss::shoot();
483 G4double gaussZ = G4RandGauss::shoot();
484 G4double gaussT = G4RandGauss::shoot();
485
486 G4double beamStartTime = m_G4Svc->GetBeamStartTime() * ns;
487 G4double beamDeltaTime = m_G4Svc->GetBeamDeltaTime() * ns;
488 G4double nBunch = m_G4Svc->GetNBunch();
489 G4double bunchTimeSigma = m_G4Svc->GetBunchTimeSigma() * ns;
490
491 G4double ran=G4UniformRand();
492 G4double beamTime = bunchTimeSigma*G4RandGauss::shoot() + beamStartTime + beamDeltaTime*int(ran*nBunch);
493
494 tmpPosX = (beamPosX + beamSizeX*gaussX ) *mm;
495 tmpPosY = (beamPosY + beamSizeY*gaussY ) *mm;
496 tmpPosZ = (beamPosZ + beamSizeZ*gaussZ ) *mm;
497 tmpT = (beamSizeZ * gaussT ) * mm/c_light +beamTime;
498
499 G4LorentzVector tmpv(tmpPosX,tmpPosY,tmpPosZ,tmpT*c_light/mm);
500
501 HepMC::GenEvent::vertex_const_iterator vitr0= hepmcevt->vertices_begin();
502 G4LorentzVector xvtx0 ;
503 xvtx0.setX( (*vitr0)-> position().x());
504 xvtx0.setY( (*vitr0)-> position().y());
505 xvtx0.setZ( (*vitr0)-> position().z());
506 xvtx0.setT( (*vitr0)-> position().t());
507 pmPosX = xvtx0.x()*mm + tmpPosX;
508 pmPosY = xvtx0.y()*mm + tmpPosY;
509 pmPosZ = xvtx0.z()*mm + tmpPosZ;
510 pmTime = xvtx0.t()*mm/c_light + tmpT;
511
512 if(m_logLevel<=4)
513 {
514 G4cout<<G4endl;
515 G4cout<<xvtx0.x()<<" "<<xvtx0.y()<<" "<<xvtx0.z()<<" "
516 <<beamSizeX*gaussX<<" "
517 <<beamSizeY*gaussY<<" "
518 <<beamSizeZ*gaussZ<<" "<<G4endl;
519 G4cout<<xvtx0.t()* mm/c_light<<" "
520 <<beamSizeZ * gaussT/sqrt(2)*mm/c_light<<" "
521 <<beamTime<<G4endl;
522 }
523 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
524 vitr != hepmcevt->vertices_end(); ++vitr )
525 {
526 G4LorentzVector xvtx;
527 xvtx.setX((*vitr)-> position().x());
528 xvtx.setY((*vitr)-> position().y());
529 xvtx.setZ((*vitr)-> position().z());
530 xvtx.setT((*vitr)-> position().t());
531 (*vitr)->set_position(xvtx+tmpv);
532 }
533 m_G4Svc->SetBeamTime(pmTime);
534
535 G4PrimaryVertex* g4vtx= new G4PrimaryVertex(pmPosX,pmPosY,pmPosZ,pmTime);
536
537 // put initial particles to the vertex
538 for( size_t ii=0; ii<HPlist.size(); ii++ )
539 {
540 if( HPlist[ii]->GetISTHEP() > 0 ) // ISTHEP of daughters had been
541 // set to negative
542 {
543 G4PrimaryParticle* initialParticle = HPlist[ii]->GetTheParticle();
544 g4vtx->SetPrimary( initialParticle );
545 }
546 }
547
548 //clear G4HepMCInterface
549 for(size_t iii=0;iii<HPlist.size();iii++)
550 { delete HPlist[iii]; }
551 HPlist.clear();
552
553 g4event->AddPrimaryVertex(g4vtx);
554}
void Print(const HepMC::GenEvent *hepmcevt)
void Boost(HepMC::GenEvent *hepmcevt)
std::vector< G4HepMCParticle * > HPlist
G4int CheckType(const HepMC::GenEvent *hepmcevt)
double GetBeamSizeZ()
Definition G4Svc.h:82
double GetBeamPosX()
Definition G4Svc.h:76
double GetBeamPosZ()
Definition G4Svc.h:78
bool GetBoostLab()
Definition G4Svc.h:98
double GetNBunch()
Definition G4Svc.h:90
double GetBeamDeltaTime()
Definition G4Svc.h:89
double GetBunchTimeSigma()
Definition G4Svc.h:91
double GetBeamStartTime()
Definition G4Svc.h:88
void SetBeamTime(double value)
Definition G4Svc.h:94
double GetBeamSizeX()
Definition G4Svc.h:80
double GetBeamPosY()
Definition G4Svc.h:77
double GetBeamSizeY()
Definition G4Svc.h:81
double getBunchSizeZ()
double getBunchSizeY()
double getBunchPosY()
double getBunchPosZ()
double getBunchSizeX()
double getBunchPosX()
#define ns(x)
Definition xmltok.c:1504

Referenced by G4SvcRunManager::GenerateEvent(), and GeneratePrimaryVertex().

◆ Print()

void G4HepMCInterface::Print ( const HepMC::GenEvent * hepmcevt)

Definition at line 125 of file G4HepMCInterface.cpp.

126{
127 G4int i=0;
128 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
129 vitr != hepmcevt->vertices_end(); ++vitr ) { // loop for vertex ...
130 G4cout<<G4endl<<"vertex:"<<i<<" barcode:"<<(*vitr)->barcode()<<" ";
131 i++;
132 //if((*vitr)->mother())
133 // G4cout<<"mother particle: "<<(*vitr)->mother()-> pdg_id()<<" ";
134 //if((*vitr)->secondMother())
135 // G4cout<<" second mother: "<<(*vitr)->secondMother()-> pdg_id()<<G4endl;
136 G4LorentzVector xvtx;
137 xvtx.setX( (*vitr)-> position().x());
138 xvtx.setY( (*vitr)-> position().y());
139 xvtx.setZ( (*vitr)-> position().z());
140 xvtx.setT( (*vitr)-> position().t());
141 G4cout<<"x:"<<xvtx.x()<<" y:"<<xvtx.y()<<" z:"<<xvtx.z()
142 <<" t:"<<xvtx.t()*mm/c_light<<G4endl;
143
144 G4int j=0;
145 for (HepMC::GenVertex::particle_iterator
146 pitr= (*vitr)->particles_begin(HepMC::children);
147 pitr != (*vitr)->particles_end(HepMC::children); ++pitr) {
148 G4int pdgcode= (*pitr)-> pdg_id();
149 G4LorentzVector p;
150 p.setX( (*pitr)-> momentum().x());
151 p.setY( (*pitr)-> momentum().y());
152 p.setZ( (*pitr)-> momentum().z());
153 p.setT( (*pitr)-> momentum().t());
154 G4cout<<"particle:"<<j<<" pdgcode:"<<pdgcode<<" ";
155 G4cout<<"p:"<<p.x()<<" "<<p.y()<<" "<<p.z()<<" ";
156 if((*pitr)->end_vertex())
157 G4cout<<"endVtx:"<<(*pitr)->end_vertex()->barcode()<<" ";
158 G4cout<<"status:"<<(*pitr)->status()<<G4endl;
159 j++;
160 }
161 }
162}

Referenced by HepMC2G4().

◆ PrintHPlist()

void G4HepMCInterface::PrintHPlist ( )

Definition at line 164 of file G4HepMCInterface.cpp.

165{
166 G4cout<<G4endl;
167 G4cout<<"particles sent to Geant4: "<<G4endl;
168 for( size_t ii=0; ii<HPlist.size(); ii++ )
169 {
170 G4int pdgcode = HPlist[ii]->GetTheParticle()->GetPDGcode();
171 G4cout<<pdgcode<<" ";
172 }
173 G4cout<<G4endl;
174}

Referenced by HepMC2G4().

◆ SetLogLevel()

void G4HepMCInterface::SetLogLevel ( G4int level)
inline

Definition at line 56 of file G4HepMCInterface.h.

56{m_logLevel = level;}

Referenced by G4SvcRunManager::GenerateEvent().

Member Data Documentation

◆ hepmcEvent

HepMC::GenEvent* G4HepMCInterface::hepmcEvent
protected

Definition at line 66 of file G4HepMCInterface.h.

Referenced by GeneratePrimaryVertex(), GetHepMCGenEvent(), and ~G4HepMCInterface().

◆ HPlist

std::vector<G4HepMCParticle*> G4HepMCInterface::HPlist

Definition at line 84 of file G4HepMCInterface.h.

Referenced by HepMC2G4(), and PrintHPlist().

◆ m_logLevel

G4int G4HepMCInterface::m_logLevel
protected

Definition at line 71 of file G4HepMCInterface.h.

Referenced by CheckType(), GetLogLevel(), HepMC2G4(), and SetLogLevel().


The documentation for this class was generated from the following files: