Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PrimaryTransformer Class Reference

#include <G4PrimaryTransformer.hh>

Public Member Functions

 G4PrimaryTransformer ()
 
virtual ~G4PrimaryTransformer ()
 
G4TrackVectorGimmePrimaries (G4Event *anEvent, G4int trackIDCounter=0)
 
void CheckUnknown ()
 
void SetVerboseLevel (G4int vl)
 
void SetUnknnownParticleDefined (G4bool vl)
 
G4bool GetUnknownParticleDefined () const
 

Protected Member Functions

void GenerateTracks (G4PrimaryVertex *primaryVertex)
 
void GenerateSingleTrack (G4PrimaryParticle *primaryParticle, G4double x0, G4double y0, G4double z0, G4double t0, G4double wv)
 
void SetDecayProducts (G4PrimaryParticle *mother, G4DynamicParticle *motherDP)
 
G4bool CheckDynamicParticle (G4DynamicParticle *DP)
 
virtual G4ParticleDefinitionGetDefinition (G4PrimaryParticle *pp)
 
virtual G4bool IsGoodForTrack (G4ParticleDefinition *pd)
 

Protected Attributes

G4TrackVector TV
 
G4ParticleTableparticleTable
 
G4int verboseLevel
 
G4int trackID
 
G4ParticleDefinitionunknown
 
G4bool unknownParticleDefined
 

Detailed Description

Definition at line 47 of file G4PrimaryTransformer.hh.

Constructor & Destructor Documentation

◆ G4PrimaryTransformer()

G4PrimaryTransformer::G4PrimaryTransformer ( )

◆ ~G4PrimaryTransformer()

G4PrimaryTransformer::~G4PrimaryTransformer ( )
virtual

Definition at line 50 of file G4PrimaryTransformer.cc.

51{;}

Member Function Documentation

◆ CheckDynamicParticle()

G4bool G4PrimaryTransformer::CheckDynamicParticle ( G4DynamicParticle DP)
protected

Definition at line 301 of file G4PrimaryTransformer.cc.

302{
303 if(IsGoodForTrack(DP->GetDefinition())) return true;
305 if(decayProducts && decayProducts->entries()>0) return true;
306 G4cerr << G4endl
307 << "G4PrimaryTransformer: a shortlived primary particle is found" << G4endl
308 << " without any valid decay table nor pre-assigned decay mode." << G4endl;
309 G4Exception("G4PrimaryTransformer","InvalidPrimary",JustWarning,
310 "This primary particle will be ignored.");
311 return false;
312}
@ JustWarning
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4int entries() const
const G4DecayProducts * GetPreAssignedDecayProducts() const
G4ParticleDefinition * GetDefinition() const
virtual G4bool IsGoodForTrack(G4ParticleDefinition *pd)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by GenerateSingleTrack(), and SetDecayProducts().

◆ CheckUnknown()

void G4PrimaryTransformer::CheckUnknown ( )

Definition at line 53 of file G4PrimaryTransformer.cc.

54{
56 if(unknown)
57 { unknownParticleDefined = true; }
58 else
59 { unknownParticleDefined = false; }
60}
G4ParticleDefinition * FindParticle(G4int PDGEncoding)

Referenced by G4PrimaryTransformer(), and G4RunManagerKernel::RunInitialization().

◆ GenerateSingleTrack()

void G4PrimaryTransformer::GenerateSingleTrack ( G4PrimaryParticle primaryParticle,
G4double  x0,
G4double  y0,
G4double  z0,
G4double  t0,
G4double  wv 
)
protected

Definition at line 107 of file G4PrimaryTransformer.cc.

110{
111 static G4ParticleDefinition* optPhoton = 0;
112 static G4int nWarn = 0;
113 if(!optPhoton) optPhoton = particleTable->FindParticle("opticalphoton");
114
115 G4ParticleDefinition* partDef = GetDefinition(primaryParticle);
116 if(!IsGoodForTrack(partDef))
117 // The particle cannot be converted to G4Track, check daughters
118 {
119#ifdef G4VERBOSE
120 if(verboseLevel>2)
121 {
122 G4cout << "Primary particle (PDGcode " << primaryParticle->GetPDGcode()
123 << ") --- Ignored" << G4endl;
124 }
125#endif
126 G4PrimaryParticle* daughter = primaryParticle->GetDaughter();
127 while(daughter)
128 {
129 GenerateSingleTrack(daughter,x0,y0,z0,t0,wv);
130 daughter = daughter->GetNext();
131 }
132 }
133
134 // The particle is defined in GEANT4
135 else
136 {
137 // Create G4DynamicParticle object
138#ifdef G4VERBOSE
139 if(verboseLevel>1)
140 {
141 G4cout << "Primary particle (" << partDef->GetParticleName()
142 << ") --- Transfered with momentum " << primaryParticle->GetMomentum()
143 << G4endl;
144 }
145#endif
146 G4DynamicParticle* DP =
147 new G4DynamicParticle(partDef,
148 primaryParticle->GetMomentumDirection(),
149 primaryParticle->GetKineticEnergy());
150 if(partDef==optPhoton && primaryParticle->GetPolarization().mag2()==0.)
151 {
152 if(nWarn<10)
153 {
154 G4Exception("G4PrimaryTransformer::GenerateSingleTrack","ZeroPolarization",JustWarning,
155 "Polarization of the optical photon is null. Random polarization is assumed.");
156 G4cerr << "This warning message is issued up to 10 times." << G4endl;
157 nWarn++;
158 }
159
160 G4double angle = G4UniformRand() * 360.0*deg;
161 G4ThreeVector normal (1., 0., 0.);
162 G4ThreeVector kphoton = DP->GetMomentumDirection();
163 G4ThreeVector product = normal.cross(kphoton);
164 G4double modul2 = product*product;
165
166 G4ThreeVector e_perpend (0., 0., 1.);
167 if (modul2 > 0.) e_perpend = (1./std::sqrt(modul2))*product;
168 G4ThreeVector e_paralle = e_perpend.cross(kphoton);
169
170 G4ThreeVector polar = std::cos(angle)*e_paralle + std::sin(angle)*e_perpend;
171 DP->SetPolarization(polar.x(),polar.y(),polar.z());
172 }
173 else
174 {
175 DP->SetPolarization(primaryParticle->GetPolX(),
176 primaryParticle->GetPolY(),
177 primaryParticle->GetPolZ());
178 }
179 if(primaryParticle->GetProperTime()>0.0)
180 { DP->SetPreAssignedDecayProperTime(primaryParticle->GetProperTime()); }
181
182 // Set Mass if it is specified
183 G4double pmas = primaryParticle->GetMass();
184 if(pmas>=0.)
185 { DP->SetMass(pmas); }
186
187 // Set Charge if it is specified
188 if (primaryParticle->GetCharge()<DBL_MAX) {
189 if (partDef->GetAtomicNumber() <0) {
190 DP->SetCharge(primaryParticle->GetCharge());
191 } else {
192 // ions
193 G4int iz = partDef->GetAtomicNumber();
194 G4int iq = static_cast<int>(primaryParticle->GetCharge()/eplus);
195 G4int n_e = iz - iq;
196 if (n_e>0) DP->AddElectron(0,n_e);
197 }
198 }
199 // Set decay products to the DynamicParticle
200 SetDecayProducts( primaryParticle, DP );
201 // Set primary particle
202 DP->SetPrimaryParticle(primaryParticle);
203 // Set PDG code if it is different from G4ParticleDefinition
204 if(partDef->GetPDGEncoding()==0 && primaryParticle->GetPDGcode()!=0)
205 {
206 DP->SetPDGcode(primaryParticle->GetPDGcode());
207 }
208 // Check the particle is properly constructed
209 if(!CheckDynamicParticle(DP))
210 {
211 delete DP;
212 return;
213 }
214 // Create G4Track object
215 G4Track* track = new G4Track(DP,t0,G4ThreeVector(x0,y0,z0));
216 // Set trackID and let primary particle know it
217 trackID++;
218 track->SetTrackID(trackID);
219 primaryParticle->SetTrackID(trackID);
220 // Set parentID to 0 as a primary particle
221 track->SetParentID(0);
222 // Set weight ( vertex weight * particle weight )
223 track->SetWeight(wv*(primaryParticle->GetWeight()));
224 // Store it to G4TrackVector
225 TV.push_back( track );
226
227 }
228}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
double z() const
double x() const
double mag2() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
void SetCharge(G4double charge)
void SetPDGcode(G4int c)
const G4ThreeVector & GetMomentumDirection() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
void SetMass(G4double mass)
void SetPrimaryParticle(G4PrimaryParticle *p)
void SetPreAssignedDecayProperTime(G4double)
void AddElectron(G4int orbit, G4int number=1)
G4int GetAtomicNumber() const
const G4String & GetParticleName() const
G4double GetWeight() const
G4double GetCharge() const
G4double GetKineticEnergy() const
G4double GetProperTime() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetPolY() const
void SetTrackID(G4int id)
G4ThreeVector GetPolarization() const
G4PrimaryParticle * GetNext() const
G4double GetMass() const
G4ThreeVector GetMomentum() const
G4int GetPDGcode() const
G4double GetPolZ() const
G4double GetPolX() const
G4PrimaryParticle * GetDaughter() const
G4bool CheckDynamicParticle(G4DynamicParticle *DP)
virtual G4ParticleDefinition * GetDefinition(G4PrimaryParticle *pp)
void SetDecayProducts(G4PrimaryParticle *mother, G4DynamicParticle *motherDP)
void GenerateSingleTrack(G4PrimaryParticle *primaryParticle, G4double x0, G4double y0, G4double z0, G4double t0, G4double wv)
void SetWeight(G4double aValue)
void SetTrackID(const G4int aValue)
void SetParentID(const G4int aValue)
#define DBL_MAX
Definition: templates.hh:83

Referenced by GenerateSingleTrack(), and GenerateTracks().

◆ GenerateTracks()

void G4PrimaryTransformer::GenerateTracks ( G4PrimaryVertex primaryVertex)
protected

Definition at line 79 of file G4PrimaryTransformer.cc.

80{
81 G4double X0 = primaryVertex->GetX0();
82 G4double Y0 = primaryVertex->GetY0();
83 G4double Z0 = primaryVertex->GetZ0();
84 G4double T0 = primaryVertex->GetT0();
85 G4double WV = primaryVertex->GetWeight();
86
87#ifdef G4VERBOSE
88 if(verboseLevel>2) {
89 primaryVertex->Print();
90 } else if (verboseLevel==1) {
91 G4cout << "G4PrimaryTransformer::PrimaryVertex ("
92 << X0 / mm << "(mm),"
93 << Y0 / mm << "(mm),"
94 << Z0 / mm << "(mm),"
95 << T0 / nanosecond << "(nsec))" << G4endl;
96 }
97#endif
98
99 G4PrimaryParticle* primaryParticle = primaryVertex->GetPrimary();
100 while( primaryParticle != 0 )
101 {
102 GenerateSingleTrack( primaryParticle, X0, Y0, Z0, T0, WV );
103 primaryParticle = primaryParticle->GetNext();
104 }
105}
G4double GetT0() const
void Print() const
G4double GetWeight() const
G4double GetZ0() const
G4double GetX0() const
G4double GetY0() const
G4PrimaryParticle * GetPrimary(G4int i=0) const

Referenced by GimmePrimaries().

◆ GetDefinition()

G4ParticleDefinition * G4PrimaryTransformer::GetDefinition ( G4PrimaryParticle pp)
protectedvirtual

Definition at line 314 of file G4PrimaryTransformer.cc.

315{
316 G4ParticleDefinition* partDef = pp->GetG4code();
317 if(!partDef) partDef = particleTable->FindParticle(pp->GetPDGcode());
318 if(unknownParticleDefined && ((!partDef)||partDef->IsShortLived())) partDef = unknown;
319 return partDef;
320}

Referenced by GenerateSingleTrack(), and SetDecayProducts().

◆ GetUnknownParticleDefined()

G4bool G4PrimaryTransformer::GetUnknownParticleDefined ( ) const
inline

Definition at line 74 of file G4PrimaryTransformer.hh.

75 { return unknownParticleDefined; }

◆ GimmePrimaries()

G4TrackVector * G4PrimaryTransformer::GimmePrimaries ( G4Event anEvent,
G4int  trackIDCounter = 0 
)

Definition at line 62 of file G4PrimaryTransformer.cc.

63{
64 trackID = trackIDCounter;
65
66 //TV.clearAndDestroy();
67 for( size_t ii=0; ii<TV.size();ii++)
68 { delete TV[ii]; }
69 TV.clear();
70 G4PrimaryVertex* nextVertex = anEvent->GetPrimaryVertex();
71 while(nextVertex)
72 {
73 GenerateTracks(nextVertex);
74 nextVertex = nextVertex->GetNext();
75 }
76 return &TV;
77}
G4PrimaryVertex * GetPrimaryVertex(G4int i=0) const
Definition: G4Event.hh:155
void GenerateTracks(G4PrimaryVertex *primaryVertex)
G4PrimaryVertex * GetNext() const

◆ IsGoodForTrack()

G4bool G4PrimaryTransformer::IsGoodForTrack ( G4ParticleDefinition pd)
protectedvirtual

Definition at line 322 of file G4PrimaryTransformer.cc.

323{
324 if(!pd)
325 { return false; }
326 else if(!(pd->IsShortLived()))
327 { return true; }
328// Following two lines should be removed if the user does not want to make shortlived
329// primary particle with proper decay table to be converted into a track.
330 else if(pd->GetDecayTable())
331 { return true; }
332 return false;
333}
G4DecayTable * GetDecayTable() const

Referenced by CheckDynamicParticle(), GenerateSingleTrack(), and SetDecayProducts().

◆ SetDecayProducts()

void G4PrimaryTransformer::SetDecayProducts ( G4PrimaryParticle mother,
G4DynamicParticle motherDP 
)
protected

Definition at line 230 of file G4PrimaryTransformer.cc.

232{
233 G4PrimaryParticle* daughter = mother->GetDaughter();
234 if(!daughter) return;
235 G4DecayProducts* decayProducts = (G4DecayProducts*)(motherDP->GetPreAssignedDecayProducts() );
236 if(!decayProducts)
237 {
238 decayProducts = new G4DecayProducts(*motherDP);
239 motherDP->SetPreAssignedDecayProducts(decayProducts);
240 }
241 while(daughter)
242 {
243 G4ParticleDefinition* partDef = GetDefinition(daughter);
244 if(!IsGoodForTrack(partDef))
245 {
246#ifdef G4VERBOSE
247 if(verboseLevel>2)
248 {
249 G4cout << " >> Decay product (PDGcode " << daughter->GetPDGcode()
250 << ") --- Ignored" << G4endl;
251 }
252#endif
253 SetDecayProducts(daughter,motherDP);
254 }
255 else
256 {
257#ifdef G4VERBOSE
258 if(verboseLevel>1)
259 {
260 G4cout << " >> Decay product (" << partDef->GetParticleName()
261 << ") --- Attached with momentum " << daughter->GetMomentum()
262 << G4endl;
263 }
264#endif
266 = new G4DynamicParticle(partDef,daughter->GetMomentum());
267 DP->SetPrimaryParticle(daughter);
268 // Decay proper time for daughter
269 if(daughter->GetProperTime()>0.0)
270 { DP->SetPreAssignedDecayProperTime(daughter->GetProperTime()); }
271 // Set Charge is specified
272 if (daughter->GetCharge()<DBL_MAX) {
273 DP->SetCharge(daughter->GetCharge());
274 }
275 G4double pmas = daughter->GetMass();
276 if(pmas>=0.)
277 { DP->SetMass(pmas); }
278 decayProducts->PushProducts(DP);
279 SetDecayProducts(daughter,DP);
280 // Check the particle is properly constructed
281 if(!CheckDynamicParticle(DP))
282 {
283 delete DP;
284 return;
285 }
286 }
287 daughter = daughter->GetNext();
288 }
289}
G4int PushProducts(G4DynamicParticle *aParticle)
void SetPreAssignedDecayProducts(G4DecayProducts *aDecayProducts)

Referenced by GenerateSingleTrack(), and SetDecayProducts().

◆ SetUnknnownParticleDefined()

void G4PrimaryTransformer::SetUnknnownParticleDefined ( G4bool  vl)

Definition at line 291 of file G4PrimaryTransformer.cc.

292{
295 { G4cerr << "unknownParticleDefined cannot be set true because G4UnknownParticle is not defined in the physics list."
296 << G4endl << "Command ignored." << G4endl;
298 }
299}

◆ SetVerboseLevel()

void G4PrimaryTransformer::SetVerboseLevel ( G4int  vl)
inline

Definition at line 66 of file G4PrimaryTransformer.hh.

67 { verboseLevel = vl; };

Referenced by G4EventManager::SetVerboseLevel().

Member Data Documentation

◆ particleTable

G4ParticleTable* G4PrimaryTransformer::particleTable
protected

◆ trackID

G4int G4PrimaryTransformer::trackID
protected

Definition at line 60 of file G4PrimaryTransformer.hh.

Referenced by GenerateSingleTrack(), and GimmePrimaries().

◆ TV

G4TrackVector G4PrimaryTransformer::TV
protected

Definition at line 57 of file G4PrimaryTransformer.hh.

Referenced by GenerateSingleTrack(), and GimmePrimaries().

◆ unknown

G4ParticleDefinition* G4PrimaryTransformer::unknown
protected

Definition at line 62 of file G4PrimaryTransformer.hh.

Referenced by CheckUnknown(), GetDefinition(), and SetUnknnownParticleDefined().

◆ unknownParticleDefined

G4bool G4PrimaryTransformer::unknownParticleDefined
protected

◆ verboseLevel

G4int G4PrimaryTransformer::verboseLevel
protected

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