BOSS 7.0.3
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)
 
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 ()
 
virtual HepMC::GenEvent * GenerateHepMCEvent ()
 

Protected Attributes

HepMC::GenEvent * hepmcEvent
 
G4int m_logLevel
 

Detailed Description

Constructor & Destructor Documentation

◆ G4HepMCInterface() [1/2]

G4HepMCInterface::G4HepMCInterface ( )

Definition at line 55 of file G4HepMCInterface.cpp.

56 : hepmcEvent(0)
57////////////////////////////////////
58{
59
60 IRealizationSvc *tmpReal;
61 ISvcLocator* svcLocator = Gaudi::svcLocator();
62 StatusCode status = svcLocator->service("RealizationSvc",tmpReal);
63 if (!status.isSuccess())
64 {
65 std::cout << " Could not initialize Realization Service in G4HepMCInterface" << std::endl;
66 } else {
67 m_RealizationSvc=dynamic_cast<RealizationSvc*>(tmpReal);
68 }
69
70}

◆ ~G4HepMCInterface() [1/2]

G4HepMCInterface::~G4HepMCInterface ( )
virtual

Definition at line 73 of file G4HepMCInterface.cpp.

75{
76 delete hepmcEvent;
77}

◆ G4HepMCInterface() [2/2]

G4HepMCInterface::G4HepMCInterface ( )

◆ ~G4HepMCInterface() [2/2]

virtual G4HepMCInterface::~G4HepMCInterface ( )
virtual

Member Function Documentation

◆ Boost() [1/2]

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

Definition at line 154 of file G4HepMCInterface.cpp.

156{
157 //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:
158 // v2 = vx*vx + vy*vy + vz*vz;
159 // gamma = 1.0 / sqrt(1.0 - v2);
160 // bp = vx*x + vy*y + vz*z;
161 // X=x + gamma2*bp*vx + gamma*vx*t;
162 // Y=y + gamma2*bp*vy + gamma*vy*t;
163 // Z=z + gamma2*bp*vz + gamma*vz*t;
164 // T=gamma*(t + bp);
165 // the function of these formulae is the same as xvtx.boost(ThreeVector v)
166 //For the E_cms,need to loop and splice the momentum of all the final state particles
167 G4LorentzVector ptot=0;
168 double E_cms,p_Mag,e_Mass2,M_cms,theta;
169 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
170 vitr != hepmcevt->vertices_end(); ++vitr ) { // loop for vertex ...
171
172 for(HepMC::GenVertex::particle_iterator
173 vpitr = (*vitr)->particles_begin(HepMC::children);
174 vpitr != (*vitr)->particles_end(HepMC::children); ++vpitr) {
175
176 if((*vpitr)->status() !=1)continue;// Not final state particle
177
178 G4LorentzVector p;
179 p.setX( (*vpitr)-> momentum().x());
180 p.setY( (*vpitr)-> momentum().y());
181 p.setZ( (*vpitr)-> momentum().z());
182 p.setT( (*vpitr)-> momentum().t());
183 ptot = p + ptot;
184 }
185 }
186 E_cms=ptot.e()*GeV;
187
188 //get G4Svc, allow user to access the beam momentum shift in JobOptions
189 ISvcLocator* svcLocator = Gaudi::svcLocator();
190 IG4Svc* tmpSvc;
191 G4Svc* m_G4Svc;
192 StatusCode sc=svcLocator->service("G4Svc", tmpSvc);
193 m_G4Svc=dynamic_cast<G4Svc *>(tmpSvc);
194
195
196 G4ThreeVector v=0; //velocity of cms as seen from lab
197 //for cms velocity
198 theta=11*mrad;
199 //theta=(m_G4Svc->GetBeamAngle())*mrad;
200 //theta is half of the angle between the two beams,namely,is the angle between the beam and Z axis.
201 M_cms=E_cms; //for p_cms=0 in the CMS
202 e_Mass2=electron_mass_c2*electron_mass_c2; //square of electron mass
203 p_Mag=sqrt((E_cms*E_cms-4*e_Mass2)/(2*(1-cos(pi-2*theta))));
204 G4ThreeVector p_Positron(p_Mag*sin(theta),0,p_Mag*cos(theta));
205 G4ThreeVector p_Electron(p_Mag*sin(pi-theta),0,p_Mag*cos(pi-theta));
206 G4ThreeVector vee=p_Positron+p_Electron;
207 G4ThreeVector beamshift(m_G4Svc->GetBeamShiftPx(),m_G4Svc->GetBeamShiftPy(),m_G4Svc->GetBeamShiftPz());
208 if(m_G4Svc-> GetSetBeamShift()==true && fabs(M_cms-3686)<50.0) {vee = beamshift;}
209
210 double pee=vee.r();
211 M_cms = sqrt(M_cms*M_cms + pee*pee);
212 v=vee/M_cms;
213
214
215 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
216 vitr != hepmcevt->vertices_end(); ++vitr ) { // loop for vertex ...
217
218 // Boost vertex position from cms to lab
219 G4LorentzVector xvtx;
220 xvtx.setX( (*vitr)-> position().x());
221 xvtx.setY( (*vitr)-> position().y());
222 xvtx.setZ( (*vitr)-> position().z());
223 xvtx.setT( (*vitr)-> position().t());
224 xvtx.boost(v);
225 (*vitr)->set_position(xvtx);
226
227 // Boost the particles four momentum from cms to lab
228 // the transform formulae are similary to the position boost
229 for (HepMC::GenVertex::particle_iterator
230 vpitr = (*vitr)->particles_begin(HepMC::children);
231 vpitr != (*vitr)->particles_end(HepMC::children); ++vpitr) {
232 G4LorentzVector p;
233 p.setX( (*vpitr)-> momentum().x());
234 p.setY( (*vpitr)-> momentum().y());
235 p.setZ( (*vpitr)-> momentum().z());
236 p.setT( (*vpitr)-> momentum().t());
237 p=p.boost(v);
238 (*vpitr)->set_momentum(p);
239 }
240 }
241}
**********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]
double sin(const BesAngle a)
double cos(const BesAngle a)
**********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
int t()
Definition: t.c:1

Referenced by HepMC2G4().

◆ Boost() [2/2]

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

◆ CheckType() [1/2]

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

Definition at line 129 of file G4HepMCInterface.cpp.

130{
131 G4int flag =0;
132 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
133 vitr != hepmcevt->vertices_end(); ++vitr ) { // loop for vertex ...
134 for (HepMC::GenVertex::particle_iterator
135 pitr= (*vitr)->particles_begin(HepMC::children);
136 pitr != (*vitr)->particles_end(HepMC::children); ++pitr) {
137 if((*pitr)->end_vertex())
138 {flag = 1; break;}
139 }
140 if(flag) break;
141 }
142 if(m_logLevel <= 4)
143 {
144 if(flag == 0)
145 G4cout<<G4endl<<"generator is GENBES!"<<G4endl;
146 else
147 G4cout<<G4endl<<"generator is not GENBES!"<<G4endl;
148 }
149 return flag;
150}

Referenced by HepMC2G4().

◆ CheckType() [2/2]

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

◆ GenerateHepMCEvent() [1/2]

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

Reimplemented in BesHepMCInterface, and BesHepMCInterface.

Definition at line 465 of file G4HepMCInterface.cpp.

467{
468 HepMC::GenEvent* aevent= new HepMC::GenEvent();
469 return aevent;
470}

Referenced by GeneratePrimaryVertex().

◆ GenerateHepMCEvent() [2/2]

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

Reimplemented in BesHepMCInterface, and BesHepMCInterface.

◆ GeneratePrimaryVertex() [1/2]

void G4HepMCInterface::GeneratePrimaryVertex ( G4Event *  anEvent)
virtual

Definition at line 473 of file G4HepMCInterface.cpp.

475{
476 // delete previous event object
477 delete hepmcEvent;
478
479 // generate next event
481 if(! hepmcEvent) {
482 G4Exception("G4HepMCInterface","GeneratePrimaryVertex",RunMustBeAborted,
483 "HepMCInterface: no generated particles. run terminated..." );
484 return;
485 }
486 HepMC2G4(hepmcEvent, anEvent);
487}
void HepMC2G4(HepMC::GenEvent *hepmcevt, G4Event *g4event)
virtual HepMC::GenEvent * GenerateHepMCEvent()

◆ GeneratePrimaryVertex() [2/2]

virtual void G4HepMCInterface::GeneratePrimaryVertex ( G4Event *  anEvent)
virtual

◆ GetHepMCGenEvent() [1/2]

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

Definition at line 94 of file InstallArea/include/G4Svc/G4Svc/G4HepMCInterface.h.

95{
96 return hepmcEvent;
97}

◆ GetHepMCGenEvent() [2/2]

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

◆ GetLogLevel() [1/2]

G4int G4HepMCInterface::GetLogLevel ( )
inline

Definition at line 55 of file InstallArea/include/G4Svc/G4Svc/G4HepMCInterface.h.

55{return m_logLevel;}

◆ GetLogLevel() [2/2]

G4int G4HepMCInterface::GetLogLevel ( )
inline

◆ HepMC2G4() [1/2]

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

Definition at line 244 of file G4HepMCInterface.cpp.

247{
248 //Print HepMC Event before boost
249 if(m_logLevel <=4 )
250 Print(hepmcevt);
251
252 //get G4Svc
253 ISvcLocator* svcLocator = Gaudi::svcLocator();
254 IG4Svc* tmpSvc;
255 G4Svc* m_G4Svc;
256 StatusCode sc=svcLocator->service("G4Svc", tmpSvc);
257 m_G4Svc=dynamic_cast<G4Svc *>(tmpSvc);
258
259 //considering 22rmad beam clossing angle ,need to boost the Cms to Lab
260 // need to boost ?
261 if(m_G4Svc->GetBoostLab()== true)
262 Boost(hepmcevt);
263
264 //Print HepMC Event after boost
265 if(m_logLevel <=2 && m_G4Svc->GetBoostLab()== true)
266 Print(hepmcevt);
267
268 //*********************send particles to HPlist*************************//
269 G4int flag = CheckType(hepmcevt);
270 if(flag)
271 {
272 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
273 vitr != hepmcevt->vertices_end(); ++vitr )
274 { // loop for vertex ...
275 G4int vtxFlag=0;
276 G4int pOut = (*vitr)->particles_out_size();
277 HepMC::GenVertex::particle_iterator pitr = (*vitr)->particles_begin(HepMC::children);
278 G4int pdgcode= (*pitr)-> pdg_id();
279 if(pOut == 1 && abs(pdgcode) == 11)
280 vtxFlag=1;
281 G4int barcodeVtx = (*vitr)->barcode();
282
283 for (HepMC::GenVertex::particle_iterator pitr= (*vitr)->particles_begin(HepMC::children);
284 pitr != (*vitr)->particles_end(HepMC::children); ++pitr)
285 {// loop for particle from this vertex
286 G4int pdgcode = (*pitr)->pdg_id();
287 if(vtxFlag==0)
288 {
289 if(pdgcode==-22) pdgcode=22;
290 G4LorentzVector p;
291 p.setX( (*pitr)-> momentum().x());
292 p.setY( (*pitr)-> momentum().y());
293 p.setZ( (*pitr)-> momentum().z());
294 p.setT( (*pitr)-> momentum().t());
295 G4PrimaryParticle* particle = new G4PrimaryParticle(pdgcode, p.x()*GeV, p.y()*GeV, p.z()*GeV);
296 //G4PrimaryParticle* particle = new G4PrimaryParticle(pdgcode);
297 //particle->Set4Momentum(p.x()*GeV,p.y()*GeV,p.z()*GeV,p.t()*GeV);
298 G4int ISTHEP = (*pitr)->status();
299
300 G4int barcodeEndVtx = 0;
301 if((*pitr)->end_vertex())
302 barcodeEndVtx = (*pitr)->end_vertex()->barcode();
303
304 G4HepMCParticle* hepmcParticle = new G4HepMCParticle(particle,ISTHEP,barcodeEndVtx);
305 HPlist.push_back(hepmcParticle);
306
307 for( size_t ii=0; ii<HPlist.size(); ii++ )
308 {
309 if(barcodeVtx == HPlist[ii]->GetBarcodeEndVtx())
310 {
311 HPlist[ii]->GetTheParticle()->SetDaughter(particle);
312 hepmcParticle->Done();
313 break;
314 }
315 }
316 }
317 }
318 }
319 }
320 else
321 {
322 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
323 vitr != hepmcevt->vertices_end(); ++vitr )
324 {
325 for (HepMC::GenVertex::particle_iterator pitr= (*vitr)->particles_begin(HepMC::children);
326 pitr != (*vitr)->particles_end(HepMC::children); ++pitr)
327 {
328 G4int ISTHEP = (*pitr)->status();
329 G4LorentzVector p;
330 p.setX( (*pitr)-> momentum().x());
331 p.setY( (*pitr)-> momentum().y());
332 p.setZ( (*pitr)-> momentum().z());
333 p.setT( (*pitr)-> momentum().t());
334
335 G4int pdgcode = (*pitr)->pdg_id();
336 G4int barcodeEndVtx = 0;
337 if((*pitr)->end_vertex())
338 barcodeEndVtx = (*pitr)->end_vertex()->barcode();
339 G4PrimaryParticle* particle = new G4PrimaryParticle(pdgcode, p.x()*GeV, p.y()*GeV, p.z()*GeV);
340 //G4PrimaryParticle* particle = new G4PrimaryParticle(pdgcode);
341 //particle->Set4Momentum(p.x()*GeV,p.y()*GeV,p.z()*GeV,p.t()*GeV);
342
343 G4HepMCParticle* hepmcParticle = new G4HepMCParticle(particle,ISTHEP,barcodeEndVtx);
344 HPlist.push_back(hepmcParticle);
345 if(ISTHEP>1)
346 hepmcParticle->Done();
347 }
348 }
349 }
350
351 //print particles in HPlist
352 if(m_logLevel<=4)
353 PrintHPlist();
354
355 //get time and position information from G4Svc
356 G4double pmPosX,pmPosY,pmPosZ,pmTime;
357 G4double tmpPosX,tmpPosY,tmpPosZ,tmpT;
358 G4double beamPosX,beamPosY,beamPosZ,beamSizeX,beamSizeY,beamSizeZ;
359
360 if(m_RealizationSvc->UseDBFlag() == false) {
361 beamPosX = m_G4Svc->GetBeamPosX();
362 beamPosY = m_G4Svc->GetBeamPosY();
363 beamPosZ = m_G4Svc->GetBeamPosZ();
364
365 beamSizeX = m_G4Svc->GetBeamSizeX();
366 beamSizeY = m_G4Svc->GetBeamSizeY();
367 beamSizeZ = m_G4Svc->GetBeamSizeZ()/sqrt(2);
368 }
369 if(m_RealizationSvc->UseDBFlag() == true && m_RealizationSvc->ifReadBunch() == true) {
370 beamPosX = m_RealizationSvc->getBunchPosX();
371 beamPosY = m_RealizationSvc->getBunchPosY();
372 beamPosZ = m_RealizationSvc->getBunchPosZ();
373
374 beamSizeX = m_RealizationSvc->getBunchSizeX();
375 beamSizeY = m_RealizationSvc->getBunchSizeY();
376 beamSizeZ = m_RealizationSvc->getBunchSizeZ();
377 }
378 if(m_RealizationSvc->UseDBFlag() == true && m_RealizationSvc->ifReadBunch() == false) {
379 beamPosX = m_G4Svc->GetBeamPosX();
380 beamPosY = m_G4Svc->GetBeamPosY();
381 beamPosZ = m_G4Svc->GetBeamPosZ();
382
383 beamSizeX = m_G4Svc->GetBeamSizeX();
384 beamSizeY = m_G4Svc->GetBeamSizeY();
385 beamSizeZ = m_G4Svc->GetBeamSizeZ()/sqrt(2);
386 }
387
388 G4double gaussX = G4RandGauss::shoot();
389 G4double gaussY = G4RandGauss::shoot();
390 G4double gaussZ = G4RandGauss::shoot();
391 G4double gaussT = G4RandGauss::shoot();
392
393 G4double beamStartTime = m_G4Svc->GetBeamStartTime() * ns;
394 G4double beamDeltaTime = m_G4Svc->GetBeamDeltaTime() * ns;
395 G4double nBunch = m_G4Svc->GetNBunch();
396 G4double bunchTimeSigma = m_G4Svc->GetBunchTimeSigma() * ns;
397
398 G4double ran=G4UniformRand();
399 G4double beamTime = bunchTimeSigma*G4RandGauss::shoot() + beamStartTime + beamDeltaTime*int(ran*nBunch);
400
401 tmpPosX = (beamPosX + beamSizeX*gaussX ) *mm;
402 tmpPosY = (beamPosY + beamSizeY*gaussY ) *mm;
403 tmpPosZ = (beamPosZ + beamSizeZ*gaussZ ) *mm;
404 tmpT = (beamSizeZ * gaussT ) * mm/c_light +beamTime;
405
406 G4LorentzVector tmpv(tmpPosX,tmpPosY,tmpPosZ,tmpT*c_light/mm);
407
408 HepMC::GenEvent::vertex_const_iterator vitr0= hepmcevt->vertices_begin();
409 G4LorentzVector xvtx0 ;
410 xvtx0.setX( (*vitr0)-> position().x());
411 xvtx0.setY( (*vitr0)-> position().y());
412 xvtx0.setZ( (*vitr0)-> position().z());
413 xvtx0.setT( (*vitr0)-> position().t());
414 pmPosX = xvtx0.x()*mm + tmpPosX;
415 pmPosY = xvtx0.y()*mm + tmpPosY;
416 pmPosZ = xvtx0.z()*mm + tmpPosZ;
417 pmTime = xvtx0.t()*mm/c_light + tmpT;
418
419 if(m_logLevel<=4)
420 {
421 G4cout<<G4endl;
422 G4cout<<xvtx0.x()<<" "<<xvtx0.y()<<" "<<xvtx0.z()<<" "
423 <<beamSizeX*gaussX<<" "
424 <<beamSizeY*gaussY<<" "
425 <<beamSizeZ*gaussZ<<" "<<G4endl;
426 G4cout<<xvtx0.t()* mm/c_light<<" "
427 <<beamSizeZ * gaussT/sqrt(2)*mm/c_light<<" "
428 <<beamTime<<G4endl;
429 }
430 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
431 vitr != hepmcevt->vertices_end(); ++vitr )
432 {
433 G4LorentzVector xvtx;
434 xvtx.setX((*vitr)-> position().x());
435 xvtx.setY((*vitr)-> position().y());
436 xvtx.setZ((*vitr)-> position().z());
437 xvtx.setT((*vitr)-> position().t());
438 (*vitr)->set_position(xvtx+tmpv);
439 }
440 m_G4Svc->SetBeamTime(pmTime);
441
442 G4PrimaryVertex* g4vtx= new G4PrimaryVertex(pmPosX,pmPosY,pmPosZ,pmTime);
443
444 // put initial particles to the vertex
445 for( size_t ii=0; ii<HPlist.size(); ii++ )
446 {
447 if( HPlist[ii]->GetISTHEP() > 0 ) // ISTHEP of daughters had been
448 // set to negative
449 {
450 G4PrimaryParticle* initialParticle = HPlist[ii]->GetTheParticle();
451 g4vtx->SetPrimary( initialParticle );
452 }
453 }
454
455 //clear G4HepMCInterface
456 for(size_t iii=0;iii<HPlist.size();iii++)
457 { delete HPlist[iii]; }
458 HPlist.clear();
459
460 g4event->AddPrimaryVertex(g4vtx);
461}
void Print(const HepMC::GenEvent *hepmcevt)
void Boost(HepMC::GenEvent *hepmcevt)
G4int CheckType(const HepMC::GenEvent *hepmcevt)
void SetBeamTime(double value)
#define ns(x)
Definition: xmltok.c:1504

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

◆ HepMC2G4() [2/2]

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

◆ Print() [1/2]

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

Definition at line 78 of file G4HepMCInterface.cpp.

79{
80 G4int i=0;
81 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
82 vitr != hepmcevt->vertices_end(); ++vitr ) { // loop for vertex ...
83 G4cout<<G4endl<<"vertex:"<<i<<" barcode:"<<(*vitr)->barcode()<<" ";
84 i++;
85 //if((*vitr)->mother())
86 // G4cout<<"mother particle: "<<(*vitr)->mother()-> pdg_id()<<" ";
87 //if((*vitr)->secondMother())
88 // G4cout<<" second mother: "<<(*vitr)->secondMother()-> pdg_id()<<G4endl;
89 G4LorentzVector xvtx;
90 xvtx.setX( (*vitr)-> position().x());
91 xvtx.setY( (*vitr)-> position().y());
92 xvtx.setZ( (*vitr)-> position().z());
93 xvtx.setT( (*vitr)-> position().t());
94 G4cout<<"x:"<<xvtx.x()<<" y:"<<xvtx.y()<<" z:"<<xvtx.z()
95 <<" t:"<<xvtx.t()*mm/c_light<<G4endl;
96
97 G4int j=0;
98 for (HepMC::GenVertex::particle_iterator
99 pitr= (*vitr)->particles_begin(HepMC::children);
100 pitr != (*vitr)->particles_end(HepMC::children); ++pitr) {
101 G4int pdgcode= (*pitr)-> pdg_id();
102 G4LorentzVector p;
103 p.setX( (*pitr)-> momentum().x());
104 p.setY( (*pitr)-> momentum().y());
105 p.setZ( (*pitr)-> momentum().z());
106 p.setT( (*pitr)-> momentum().t());
107 G4cout<<"particle:"<<j<<" pdgcode:"<<pdgcode<<" ";
108 G4cout<<"p:"<<p.x()<<" "<<p.y()<<" "<<p.z()<<" ";
109 if((*pitr)->end_vertex())
110 G4cout<<"endVtx:"<<(*pitr)->end_vertex()->barcode()<<" ";
111 G4cout<<"status:"<<(*pitr)->status()<<G4endl;
112 j++;
113 }
114 }
115}

Referenced by HepMC2G4().

◆ Print() [2/2]

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

◆ PrintHPlist() [1/2]

void G4HepMCInterface::PrintHPlist ( )

Definition at line 117 of file G4HepMCInterface.cpp.

118{
119 G4cout<<G4endl;
120 G4cout<<"particles sent to Geant4: "<<G4endl;
121 for( size_t ii=0; ii<HPlist.size(); ii++ )
122 {
123 G4int pdgcode = HPlist[ii]->GetTheParticle()->GetPDGcode();
124 G4cout<<pdgcode<<" ";
125 }
126 G4cout<<G4endl;
127}

Referenced by HepMC2G4().

◆ PrintHPlist() [2/2]

void G4HepMCInterface::PrintHPlist ( )

◆ SetLogLevel() [1/2]

void G4HepMCInterface::SetLogLevel ( G4int  level)
inline

◆ SetLogLevel() [2/2]

void G4HepMCInterface::SetLogLevel ( G4int  level)
inline

Member Data Documentation

◆ hepmcEvent

HepMC::GenEvent * G4HepMCInterface::hepmcEvent
protected

◆ HPlist

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

Definition at line 84 of file InstallArea/include/G4Svc/G4Svc/G4HepMCInterface.h.

Referenced by HepMC2G4(), and PrintHPlist().

◆ m_logLevel

G4int G4HepMCInterface::m_logLevel
protected

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