Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
HepMCG4Interface Class Reference

A base class for primary generation via HepMC object. More...

#include <Doxymodules_eventgenerator.h>

Inheritance diagram for HepMCG4Interface:
G4VPrimaryGenerator HepMCG4AsciiReader HepMCG4PythiaInterface

Public Member Functions

 HepMCG4Interface ()
 
virtual ~HepMCG4Interface ()
 
HepMC::GenEvent * GetHepMCGenEvent () const
 
virtual void GeneratePrimaryVertex (G4Event *anEvent)
 

Protected Member Functions

virtual G4bool CheckVertexInsideWorld (const G4ThreeVector &pos) const
 
void HepMC2G4 (const HepMC::GenEvent *hepmcevt, G4Event *g4event)
 
virtual HepMC::GenEvent * GenerateHepMCEvent ()
 

Protected Attributes

HepMC::GenEvent * hepmcEvent
 

Detailed Description

A base class for primary generation via HepMC object.

This class is derived from G4VPrimaryGenerator.

Definition at line 58 of file Doxymodules_eventgenerator.h.

Constructor & Destructor Documentation

◆ HepMCG4Interface()

HepMCG4Interface::HepMCG4Interface ( )

Definition at line 43 of file HepMCG4Interface.cc.

44 : hepmcEvent(0)
45{
46}
HepMC::GenEvent * hepmcEvent

◆ ~HepMCG4Interface()

HepMCG4Interface::~HepMCG4Interface ( )
virtual

Definition at line 49 of file HepMCG4Interface.cc.

50{
51 delete hepmcEvent;
52}

Member Function Documentation

◆ CheckVertexInsideWorld()

G4bool HepMCG4Interface::CheckVertexInsideWorld ( const G4ThreeVector &  pos) const
protectedvirtual

Definition at line 55 of file HepMCG4Interface.cc.

57{
58 G4Navigator* navigator= G4TransportationManager::GetTransportationManager()
59 -> GetNavigatorForTracking();
60
61 G4VPhysicalVolume* world= navigator-> GetWorldVolume();
62 G4VSolid* solid= world-> GetLogicalVolume()-> GetSolid();
63 EInside qinside= solid-> Inside(pos);
64
65 if( qinside != kInside) return false;
66 else return true;
67}

◆ HepMC2G4()

void HepMCG4Interface::HepMC2G4 ( const HepMC::GenEvent *  hepmcevt,
G4Event g4event 
)
protected

Definition at line 70 of file HepMCG4Interface.cc.

72{
73 for(HepMC::GenEvent::vertex_const_iterator vitr= hepmcevt->vertices_begin();
74 vitr != hepmcevt->vertices_end(); ++vitr ) { // loop for vertex ...
75
76 // real vertex?
77 G4bool qvtx=false;
78 for (HepMC::GenVertex::particle_iterator
79 pitr= (*vitr)->particles_begin(HepMC::children);
80 pitr != (*vitr)->particles_end(HepMC::children); ++pitr) {
81
82 if (!(*pitr)->end_vertex() && (*pitr)->status()==1) {
83 qvtx=true;
84 break;
85 }
86 }
87 if (!qvtx) continue;
88
89 // check world boundary
90 HepMC::FourVector pos= (*vitr)-> position();
91 G4LorentzVector xvtx(pos.x(), pos.y(), pos.z(), pos.t());
92 if (! CheckVertexInsideWorld(xvtx.vect()*mm)) continue;
93
94 // create G4PrimaryVertex and associated G4PrimaryParticles
95 G4PrimaryVertex* g4vtx=
96 new G4PrimaryVertex(xvtx.x()*mm, xvtx.y()*mm, xvtx.z()*mm,
97 xvtx.t()*mm/c_light);
98
99 for (HepMC::GenVertex::particle_iterator
100 vpitr= (*vitr)->particles_begin(HepMC::children);
101 vpitr != (*vitr)->particles_end(HepMC::children); ++vpitr) {
102
103 if( (*vpitr)->status() != 1 ) continue;
104
105 G4int pdgcode= (*vpitr)-> pdg_id();
106 pos= (*vpitr)-> momentum();
107 G4LorentzVector p(pos.px(), pos.py(), pos.pz(), pos.e());
108 G4PrimaryParticle* g4prim=
109 new G4PrimaryParticle(pdgcode, p.x()*GeV, p.y()*GeV, p.z()*GeV);
110
111 g4vtx-> SetPrimary(g4prim);
112 }
113 g4event-> AddPrimaryVertex(g4vtx);
114 }
115}
virtual G4bool CheckVertexInsideWorld(const G4ThreeVector &pos) const

◆ GenerateHepMCEvent()

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

Reimplemented in HepMCG4AsciiReader, and HepMCG4PythiaInterface.

Definition at line 118 of file HepMCG4Interface.cc.

119{
120 HepMC::GenEvent* aevent= new HepMC::GenEvent();
121 return aevent;
122}

◆ GetHepMCGenEvent()

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

Definition at line 73 of file HepMCG4Interface.hh.

74{
75 return hepmcEvent;
76}

◆ GeneratePrimaryVertex()

void HepMCG4Interface::GeneratePrimaryVertex ( G4Event anEvent)
virtual

Definition at line 125 of file HepMCG4Interface.cc.

126{
127 // delete previous event object
128 delete hepmcEvent;
129
130 // generate next event
132 if(! hepmcEvent) {
133 G4cout << "HepMCInterface: no generated particles. run terminated..."
134 << G4endl;
135 G4RunManager::GetRunManager()-> AbortRun();
136 return;
137 }
138 HepMC2G4(hepmcEvent, anEvent);
139}
void HepMC2G4(const HepMC::GenEvent *hepmcevt, G4Event *g4event)
virtual HepMC::GenEvent * GenerateHepMCEvent()

Member Data Documentation

◆ hepmcEvent

HepMC::GenEvent* HepMCG4Interface::hepmcEvent
protected

Definition at line 45 of file HepMCG4Interface.hh.


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

Applications | User Support | Publications | Collaboration