Geant4 10.7.0
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
45{
48}
49
51{
52}
53
55{
57 if(unknown != nullptr)
58 { unknownParticleDefined = true; }
59 else
60 { unknownParticleDefined = false; }
61 opticalphoton = particleTable->FindParticle("opticalphoton");
62 if(opticalphoton != nullptr)
63 { opticalphotonDefined = true; }
64 else
65 { opticalphotonDefined = false; }
66}
67
70{
71 trackID = trackIDCounter;
72
73 for(auto tr : TV) delete tr;
74 TV.clear();
75
76 // Loop over vertices
77 //
78 G4PrimaryVertex* nextVertex = anEvent->GetPrimaryVertex();
79 while(nextVertex != nullptr) // Loop checking 12.28.2015 M.Asai
80 {
81 GenerateTracks(nextVertex);
82 nextVertex = nextVertex->GetNext();
83 }
84 return &TV;
85}
86
88{
89 G4double X0 = primaryVertex->GetX0();
90 G4double Y0 = primaryVertex->GetY0();
91 G4double Z0 = primaryVertex->GetZ0();
92 G4double T0 = primaryVertex->GetT0();
93 G4double WV = primaryVertex->GetWeight();
94
95#ifdef G4VERBOSE
96 if(verboseLevel>2)
97 {
98 primaryVertex->Print();
99 }
100 else if (verboseLevel==1)
101 {
102 G4cout << "G4PrimaryTransformer::PrimaryVertex ("
103 << X0 / mm << "(mm),"
104 << Y0 / mm << "(mm),"
105 << Z0 / mm << "(mm),"
106 << T0 / nanosecond << "(nsec))" << G4endl;
107 }
108#endif
109
110 G4PrimaryParticle* primaryParticle = primaryVertex->GetPrimary();
111 while( primaryParticle != nullptr ) // Loop checking 12.28.2015 M.Asai
112 {
113 GenerateSingleTrack( primaryParticle, X0, Y0, Z0, T0, WV );
114 primaryParticle = primaryParticle->GetNext();
115 }
116}
117
119GenerateSingleTrack( G4PrimaryParticle* primaryParticle,
120 G4double x0, G4double y0, G4double z0,
121 G4double t0, G4double wv)
122{
123 G4ParticleDefinition* partDef = GetDefinition(primaryParticle);
124 if(!IsGoodForTrack(partDef))
125 { // The particle cannot be converted to G4Track, check daughters
126#ifdef G4VERBOSE
127 if(verboseLevel>2)
128 {
129 G4cout << "Primary particle (PDGcode " << primaryParticle->GetPDGcode()
130 << ") --- Ignored" << G4endl;
131 }
132#endif
133 G4PrimaryParticle* daughter = primaryParticle->GetDaughter();
134 while(daughter != nullptr) // Loop checking 12.28.2015 M.Asai
135 {
136 GenerateSingleTrack(daughter,x0,y0,z0,t0,wv);
137 daughter = daughter->GetNext();
138 }
139 }
140 else // The particle is defined in GEANT4
141 {
142 // Create G4DynamicParticle object
143#ifdef G4VERBOSE
144 if(verboseLevel>1)
145 {
146 G4cout << "Primary particle (" << partDef->GetParticleName()
147 << ") --- Transferred with momentum "
148 << primaryParticle->GetMomentum()
149 << G4endl;
150 }
151#endif
152 G4DynamicParticle* DP =
153 new G4DynamicParticle(partDef,
154 primaryParticle->GetMomentumDirection(),
155 primaryParticle->GetKineticEnergy());
157 && primaryParticle->GetPolarization().mag2()==0.)
158 {
159 if(nWarn<10)
160 {
161 G4Exception("G4PrimaryTransformer::GenerateSingleTrack",
162 "ZeroPolarization", JustWarning,
163 "Polarization of the optical photon is null.\
164 Random polarization is assumed.");
165 G4cerr << "This warning message is issued up to 10 times." << G4endl;
166 ++nWarn;
167 }
168
169 G4double angle = G4UniformRand() * 360.0*deg;
170 G4ThreeVector normal (1., 0., 0.);
171 G4ThreeVector kphoton = DP->GetMomentumDirection();
172 G4ThreeVector product = normal.cross(kphoton);
173 G4double modul2 = product*product;
174
175 G4ThreeVector e_perpend (0., 0., 1.);
176 if (modul2 > 0.) e_perpend = (1./std::sqrt(modul2))*product;
177 G4ThreeVector e_paralle = e_perpend.cross(kphoton);
178
179 G4ThreeVector polar = std::cos(angle)*e_paralle
180 + std::sin(angle)*e_perpend;
181 DP->SetPolarization(polar.x(),polar.y(),polar.z());
182 }
183 else
184 {
185 DP->SetPolarization(primaryParticle->GetPolX(),
186 primaryParticle->GetPolY(),
187 primaryParticle->GetPolZ());
188 }
189 if(primaryParticle->GetProperTime()>=0.0)
190 {
191 DP->SetPreAssignedDecayProperTime(primaryParticle->GetProperTime());
192 }
193
194 // Set Mass if it is specified
195 //
196 G4double pmas = primaryParticle->GetMass();
197 if(pmas>=0.) { DP->SetMass(pmas); }
198
199 // Set Charge if it is specified
200 //
201 if (primaryParticle->GetCharge()<DBL_MAX)
202 {
203 if (partDef->GetAtomicNumber() <0)
204 {
205 DP->SetCharge(primaryParticle->GetCharge());
206 }
207 else // ions
208 {
209 G4int iz = partDef->GetAtomicNumber();
210 G4int iq = static_cast<G4int>(primaryParticle->GetCharge()/eplus);
211 G4int n_e = iz - iq;
212 if (n_e>0) DP->AddElectron(0,n_e);
213 }
214 }
215 // Set decay products to the DynamicParticle
216 //
217 SetDecayProducts( primaryParticle, DP );
218
219 // Set primary particle
220 //
221 DP->SetPrimaryParticle(primaryParticle);
222
223 // Set PDG code if it is different from G4ParticleDefinition
224 //
225 if(partDef->GetPDGEncoding()==0 && primaryParticle->GetPDGcode()!=0)
226 {
227 DP->SetPDGcode(primaryParticle->GetPDGcode());
228 }
229
230 // Check the particle is properly constructed
231 //
232 if(!CheckDynamicParticle(DP))
233 {
234 delete DP;
235 return;
236 }
237
238 // Create G4Track object
239 //
240 G4Track* track = new G4Track(DP,t0,G4ThreeVector(x0,y0,z0));
241
242 // Set trackID and let primary particle know it
243 //
244 ++trackID;
245 track->SetTrackID(trackID);
246 primaryParticle->SetTrackID(trackID);
247
248 // Set parentID to 0 as a primary particle
249 //
250 track->SetParentID(0);
251
252 // Set weight ( vertex weight * particle weight )
253 //
254 track->SetWeight(wv*(primaryParticle->GetWeight()));
255
256 // Store it to G4TrackVector
257 //
258 TV.push_back( track );
259 }
260}
261
264{
265 G4PrimaryParticle* daughter = mother->GetDaughter();
266 if(daughter == nullptr) return;
267 G4DecayProducts* decayProducts
269 if(decayProducts == nullptr)
270 {
271 decayProducts = new G4DecayProducts(*motherDP);
272 motherDP->SetPreAssignedDecayProducts(decayProducts);
273 }
274 while(daughter != nullptr)
275 {
276 G4ParticleDefinition* partDef = GetDefinition(daughter);
277 if(!IsGoodForTrack(partDef))
278 {
279#ifdef G4VERBOSE
280 if(verboseLevel>2)
281 {
282 G4cout << " >> Decay product (PDGcode " << daughter->GetPDGcode()
283 << ") --- Ignored" << G4endl;
284 }
285#endif
286 SetDecayProducts(daughter,motherDP);
287 }
288 else
289 {
290#ifdef G4VERBOSE
291 if(verboseLevel>1)
292 {
293 G4cout << " >> Decay product (" << partDef->GetParticleName()
294 << ") --- Attached with momentum " << daughter->GetMomentum()
295 << G4endl;
296 }
297#endif
299 = new G4DynamicParticle(partDef,daughter->GetMomentum());
300 DP->SetPrimaryParticle(daughter);
301
302 // Decay proper time for daughter
303 //
304 if(daughter->GetProperTime()>=0.0)
305 {
307 }
308
309 // Set Charge and Mass is specified
310 //
311 if (daughter->GetCharge()<DBL_MAX)
312 {
313 DP->SetCharge(daughter->GetCharge());
314 }
315 G4double pmas = daughter->GetMass();
316 if(pmas>=0.)
317 {
318 DP->SetMass(pmas);
319 }
320
321 // Set Polarization
322 //
323 DP->SetPolarization(daughter->GetPolX(),
324 daughter->GetPolY(),
325 daughter->GetPolZ());
326 decayProducts->PushProducts(DP);
327 SetDecayProducts(daughter,DP);
328
329 // Check the particle is properly constructed
330 //
331 if(!CheckDynamicParticle(DP))
332 {
333 delete DP;
334 return;
335 }
336 }
337 daughter = daughter->GetNext();
338 }
339}
340
342{
345 {
346 G4cerr << "unknownParticleDefined cannot be set true because" << G4endl
347 << "G4UnknownParticle is not defined in the physics list." << G4endl
348 << "Command ignored." << G4endl;
350 }
351}
352
354{
355 if(IsGoodForTrack(DP->GetDefinition())) return true;
356 G4DecayProducts* decayProducts
358 if(decayProducts != nullptr && decayProducts->entries()>0) return true;
359 G4cerr << G4endl
360 << "G4PrimaryTransformer: a shortlived primary particle is found"
361 << G4endl
362 << " without any valid decay table nor pre-assigned decay mode."
363 << G4endl;
364 G4Exception("G4PrimaryTransformer", "InvalidPrimary", JustWarning,
365 "This primary particle will be ignored.");
366 return false;
367}
368
371{
372 G4ParticleDefinition* partDef = pp->GetG4code();
373 if(partDef == nullptr)
374 {
375 partDef = particleTable->FindParticle(pp->GetPDGcode());
376 }
377 if(unknownParticleDefined && ((!partDef)||partDef->IsShortLived()))
378 {
379 partDef = unknown;
380 }
381 return partDef;
382}
383
385{
386 if(pd == nullptr)
387 { return false; }
388 else if(!(pd->IsShortLived()))
389 { return true; }
390 //
391 // Following two lines should be removed if the user does not want to make
392 // shortlived primary particle with proper decay table to be converted into
393 // a track.
394 //
395 else if(pd->GetDecayTable() != nullptr)
396 { return true; }
397
398 return false;
399}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
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:57
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
G4int entries() const
G4int PushProducts(G4DynamicParticle *aParticle)
void SetPreAssignedDecayProducts(G4DecayProducts *aDecayProducts)
void SetPolarization(const G4ThreeVector &)
void SetCharge(G4double charge)
void SetPDGcode(G4int c)
const G4ThreeVector & GetMomentumDirection() const
void SetMass(G4double mass)
const G4DecayProducts * GetPreAssignedDecayProducts() const
G4ParticleDefinition * GetDefinition() const
void SetPrimaryParticle(G4PrimaryParticle *p)
void SetPreAssignedDecayProperTime(G4double)
void AddElectron(G4int orbit, G4int number=1)
G4PrimaryVertex * GetPrimaryVertex(G4int i=0) const
Definition: G4Event.hh:137
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
G4int GetPDGcode() 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
void Print() const
G4PrimaryVertex * GetNext() const
G4double GetWeight() const
G4double GetZ0() const
G4double GetX0() const
G4double GetY0() const
G4PrimaryParticle * GetPrimary(G4int i=0) const
void SetWeight(G4double aValue)
void SetTrackID(const G4int aValue)
void SetParentID(const G4int aValue)
#define DBL_MAX
Definition: templates.hh:62