Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PrimaryTransformer.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4PrimaryTransformer class implementation
27//
28// Author: Makoto Asai, 1999
29// --------------------------------------------------------------------
30
32#include "G4SystemOfUnits.hh"
33#include "G4Event.hh"
34#include "G4PrimaryVertex.hh"
36#include "G4DynamicParticle.hh"
37#include "G4Track.hh"
38#include "G4ThreeVector.hh"
39#include "G4DecayProducts.hh"
40#include "G4UnitsTable.hh"
41#include "G4ios.hh"
42#include "Randomize.hh"
43
49
57
60{
61 trackID = trackIDCounter;
62
63 for(auto tr : TV) delete tr;
64 TV.clear();
65
66 // Loop over vertices
67 //
68 G4PrimaryVertex* nextVertex = anEvent->GetPrimaryVertex();
69 while(nextVertex != nullptr) // Loop checking 12.28.2015 M.Asai
70 {
71 GenerateTracks(nextVertex);
72 nextVertex = nextVertex->GetNext();
73 }
74 return &TV;
75}
76
78{
79 G4double X0 = primaryVertex->GetX0();
80 G4double Y0 = primaryVertex->GetY0();
81 G4double Z0 = primaryVertex->GetZ0();
82 G4double T0 = primaryVertex->GetT0();
83 G4double WV = primaryVertex->GetWeight();
84
85#ifdef G4VERBOSE
86 if(verboseLevel>2)
87 {
88 primaryVertex->Print();
89 }
90 else if (verboseLevel==1)
91 {
92 G4cout << "G4PrimaryTransformer::PrimaryVertex ("
93 << X0 / mm << "(mm),"
94 << Y0 / mm << "(mm),"
95 << Z0 / mm << "(mm),"
96 << T0 / nanosecond << "(nsec))" << G4endl;
97 }
98#endif
99
100 G4PrimaryParticle* primaryParticle = primaryVertex->GetPrimary();
101 while( primaryParticle != nullptr ) // Loop checking 12.28.2015 M.Asai
102 {
103 GenerateSingleTrack( primaryParticle, X0, Y0, Z0, T0, WV );
104 primaryParticle = primaryParticle->GetNext();
105 }
106}
107
109GenerateSingleTrack( G4PrimaryParticle* primaryParticle,
110 G4double x0, G4double y0, G4double z0,
111 G4double t0, G4double wv)
112{
113 G4ParticleDefinition* partDef = GetDefinition(primaryParticle);
114 if(!IsGoodForTrack(partDef))
115 { // The particle cannot be converted to G4Track, check daughters
116#ifdef G4VERBOSE
117 if(verboseLevel>2)
118 {
119 G4cout << "Primary particle (PDGcode " << primaryParticle->GetPDGcode()
120 << ") --- Ignored" << G4endl;
121 }
122#endif
123 G4PrimaryParticle* daughter = primaryParticle->GetDaughter();
124 while(daughter != nullptr) // Loop checking 12.28.2015 M.Asai
125 {
126 GenerateSingleTrack(daughter,x0,y0,z0,t0,wv);
127 daughter = daughter->GetNext();
128 }
129 }
130 else // The particle is defined in GEANT4
131 {
132 // Create G4DynamicParticle object
133#ifdef G4VERBOSE
134 if(verboseLevel>1)
135 {
136 G4cout << "Primary particle (" << partDef->GetParticleName()
137 << ") --- Transferred with momentum "
138 << primaryParticle->GetMomentum()
139 << G4endl;
140 }
141#endif
142 auto* DP =
143 new G4DynamicParticle(partDef,
144 primaryParticle->GetMomentumDirection(),
145 primaryParticle->GetKineticEnergy());
147 && primaryParticle->GetPolarization().mag2()==0.)
148 {
149 if(nWarn<10)
150 {
151 G4Exception("G4PrimaryTransformer::GenerateSingleTrack",
152 "ZeroPolarization", JustWarning,
153 "Polarization of the optical photon is null.\
154 Random polarization is assumed.");
155 G4cerr << "This warning message is issued up to 10 times." << G4endl;
156 ++nWarn;
157 }
158
159 G4double angle = G4UniformRand() * 360.0*deg;
160 G4ThreeVector normal (1., 0., 0.);
161 G4ThreeVector kphoton = DP->GetMomentumDirection();
162 G4ThreeVector product = normal.cross(kphoton);
163 G4double modul2 = product*product;
164
165 G4ThreeVector e_perpend (0., 0., 1.);
166 if (modul2 > 0.) e_perpend = (1./std::sqrt(modul2))*product;
167 G4ThreeVector e_paralle = e_perpend.cross(kphoton);
168
169 G4ThreeVector polar = std::cos(angle)*e_paralle
170 + 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 {
181 DP->SetPreAssignedDecayProperTime(primaryParticle->GetProperTime());
182 }
183
184 // Set Mass if it is specified
185 //
186 G4double pmas = primaryParticle->GetMass();
187 if(pmas>=0.) { DP->SetMass(pmas); }
188
189 // Set Charge if it is specified
190 //
191 if (primaryParticle->GetCharge()<DBL_MAX)
192 {
193 if (partDef->GetAtomicNumber() <0)
194 {
195 DP->SetCharge(primaryParticle->GetCharge());
196 }
197 else // ions
198 {
199 G4int iz = partDef->GetAtomicNumber();
200 auto iq = static_cast<G4int>(primaryParticle->GetCharge()/eplus);
201 G4int n_e = iz - iq;
202 if (n_e>0) DP->AddElectron(0,n_e);
203 }
204 }
205 // Set decay products to the DynamicParticle
206 //
207 SetDecayProducts( primaryParticle, DP );
208
209 // Set primary particle
210 //
211 DP->SetPrimaryParticle(primaryParticle);
212
213 // Set PDG code if it is different from G4ParticleDefinition
214 //
215 if(partDef->GetPDGEncoding()==0 && primaryParticle->GetPDGcode()!=0)
216 {
217 DP->SetPDGcode(primaryParticle->GetPDGcode());
218 }
219
220 // Check the particle is properly constructed
221 //
222 if(!CheckDynamicParticle(DP))
223 {
224 delete DP;
225 return;
226 }
227
228 // Create G4Track object
229 //
230 auto track = new G4Track(DP,t0,G4ThreeVector(x0,y0,z0));
231
232 // Set trackID and let primary particle know it
233 //
234 ++trackID;
235 track->SetTrackID(trackID);
236 primaryParticle->SetTrackID(trackID);
237
238 // Set parentID to 0 as a primary particle
239 //
240 track->SetParentID(0);
241
242 // Set weight ( vertex weight * particle weight )
243 //
244 track->SetWeight(wv*(primaryParticle->GetWeight()));
245
246 // Store it to G4TrackVector
247 //
248 TV.push_back( track );
249 }
250}
251
254{
255 G4PrimaryParticle* daughter = mother->GetDaughter();
256 if(daughter == nullptr) return;
257 auto* decayProducts
259 if(decayProducts == nullptr)
260 {
261 decayProducts = new G4DecayProducts(*motherDP);
262 motherDP->SetPreAssignedDecayProducts(decayProducts);
263 }
264 while(daughter != nullptr)
265 {
266 G4ParticleDefinition* partDef = GetDefinition(daughter);
267 if(!IsGoodForTrack(partDef))
268 {
269#ifdef G4VERBOSE
270 if(verboseLevel>2)
271 {
272 G4cout << " >> Decay product (PDGcode " << daughter->GetPDGcode()
273 << ") --- Ignored" << G4endl;
274 }
275#endif
276 SetDecayProducts(daughter,motherDP);
277 }
278 else
279 {
280#ifdef G4VERBOSE
281 if(verboseLevel>1)
282 {
283 G4cout << " >> Decay product (" << partDef->GetParticleName()
284 << ") --- Attached with momentum " << daughter->GetMomentum()
285 << G4endl;
286 }
287#endif
288 auto* DP
289 = new G4DynamicParticle(partDef,daughter->GetMomentum());
290 DP->SetPrimaryParticle(daughter);
291
292 // Decay proper time for daughter
293 //
294 if(daughter->GetProperTime()>=0.0)
295 {
296 DP->SetPreAssignedDecayProperTime(daughter->GetProperTime());
297 }
298
299 // Set Charge and Mass is specified
300 //
301 if (daughter->GetCharge()<DBL_MAX)
302 {
303 DP->SetCharge(daughter->GetCharge());
304 }
305 G4double pmas = daughter->GetMass();
306 if(pmas>=0.)
307 {
308 DP->SetMass(pmas);
309 }
310
311 // Set Polarization
312 //
313 DP->SetPolarization(daughter->GetPolX(),
314 daughter->GetPolY(),
315 daughter->GetPolZ());
316 decayProducts->PushProducts(DP);
317 SetDecayProducts(daughter,DP);
318
319 // Check the particle is properly constructed
320 //
321 if(!CheckDynamicParticle(DP))
322 {
323 delete DP;
324 return;
325 }
326 }
327 daughter = daughter->GetNext();
328 }
329}
330
332{
334 if(unknownParticleDefined && (unknown == nullptr))
335 {
336 G4cerr << "unknownParticleDefined cannot be set true because" << G4endl
337 << "G4UnknownParticle is not defined in the physics list." << G4endl
338 << "Command ignored." << G4endl;
340 }
341}
342
344{
345 if(IsGoodForTrack(DP->GetDefinition())) return true;
346 auto* decayProducts
348 if(decayProducts != nullptr && decayProducts->entries()>0) return true;
349 G4cerr << G4endl
350 << "G4PrimaryTransformer: a shortlived primary particle is found"
351 << G4endl
352 << " without any valid decay table nor pre-assigned decay mode."
353 << G4endl;
354 G4Exception("G4PrimaryTransformer", "InvalidPrimary", JustWarning,
355 "This primary particle will be ignored.");
356 return false;
357}
358
361{
362 G4ParticleDefinition* partDef = pp->GetG4code();
363 if(partDef == nullptr)
364 {
365 partDef = particleTable->FindParticle(pp->GetPDGcode());
366 }
367 if(unknownParticleDefined && ((partDef == nullptr)||partDef->IsShortLived()))
368 {
369 partDef = unknown;
370 }
371 return partDef;
372}
373
375{
376 if(pd == nullptr)
377 { return false; }
378 if(!(pd->IsShortLived()))
379 { return true; }
380 //
381 // Following two lines should be removed if the user does not want to make
382 // shortlived primary particle with proper decay table to be converted into
383 // a track.
384 //
385 if(pd->GetDecayTable() != nullptr)
386 { return true; }
387
388 return false;
389}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
std::vector< G4Track * > G4TrackVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
double z() const
double x() const
double mag2() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
void SetPreAssignedDecayProducts(G4DecayProducts *aDecayProducts)
const G4DecayProducts * GetPreAssignedDecayProducts() const
G4ParticleDefinition * GetDefinition() const
G4PrimaryVertex * GetPrimaryVertex(G4int i=0) const
Definition G4Event.hh:142
G4int GetAtomicNumber() const
G4DecayTable * GetDecayTable() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
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
G4double GetPolZ() const
G4double GetPolX() const
G4PrimaryParticle * GetDaughter() const
G4bool CheckDynamicParticle(G4DynamicParticle *DP)
void SetUnknnownParticleDefined(G4bool vl)
virtual G4ParticleDefinition * GetDefinition(G4PrimaryParticle *pp)
G4ParticleTable * particleTable
G4ParticleDefinition * unknown
G4ParticleDefinition * opticalphoton
virtual G4bool IsGoodForTrack(G4ParticleDefinition *pd)
void SetDecayProducts(G4PrimaryParticle *mother, G4DynamicParticle *motherDP)
G4TrackVector * GimmePrimaries(G4Event *anEvent, G4int trackIDCounter=0)
void GenerateSingleTrack(G4PrimaryParticle *primaryParticle, G4double x0, G4double y0, G4double z0, G4double t0, G4double wv)
void GenerateTracks(G4PrimaryVertex *primaryVertex)
G4double GetT0() const
G4PrimaryVertex * GetNext() const
G4double GetWeight() const
G4double GetZ0() const
G4double GetX0() const
G4double GetY0() const
G4PrimaryParticle * GetPrimary(G4int i=0) const
#define DBL_MAX
Definition templates.hh:62