Geant4-11
G4PrimaryParticle.hh
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// G4PrimaryParticle
27//
28// Class description:
29//
30// This class represents a primary particle.
31// G4PrimaryParticle is a completely different object from G4Track or
32// G4DynamicParticle. G4PrimaryParticle is designed to take into account
33// the possibility of making the object persistent, i.e. kept with G4Event
34// class object within an ODBMS. This class is therefore almost independent
35// from any other Geant4 object. The only exception is a pointer to
36// G4ParticleDefinition, which can be rebuilt by its PDGcode.
37//
38// Primary particles are stored in G4PrimaryVertex through a form of linked
39// list. A G4PrimaryParticle object can have one or more objects of this
40// class as its daughters as a form of linked list.
41// A primary particle represented by this class object needs not to be
42// a particle type which Geant4 can simulate:
43// Case a) mother particle is not a particle Geant4 can simulate
44// daughters associated to the mother will be examined.
45// Case b) mother particle is a perticle Geant4 can simulate
46// daughters associated to the mother will be converted to
47// G4DynamicParticle and be set to the mother G4Track.
48// For this case, daugthers are used as the "pre-fixed"
49// decay channel.
50
51// Authors: G.Cosmo, 2 December 1995 - Design, based on object model
52// M.Asai, 29 January 1996 - First implementation
53// --------------------------------------------------------------------
54#ifndef G4PrimaryParticle_hh
55#define G4PrimaryParticle_hh 1
56
57#include "globals.hh"
58#include "G4Allocator.hh"
59#include "G4ThreeVector.hh"
60
61#include "pwdefs.hh"
62
65
67{
68 public:
69
73 G4double px, G4double py, G4double pz);
75 G4double px,G4double py, G4double pz, G4double E);
78 G4double px,G4double py,G4double pz);
80 G4double px,G4double py, G4double pz, G4double E);
81 // Constructors
82
83 virtual ~G4PrimaryParticle();
84 // Destructor
85
88 // Copy constructor and assignment operator
89 // NOTE: nextParticle and daughterParticle are copied by object
90 // (i.e. deep copy); userInfo will not be copied
91
92 G4bool operator==(const G4PrimaryParticle& right) const;
93 G4bool operator!=(const G4PrimaryParticle& right) const;
94 // Equality operator returns 'true' only if same object
95 // (i.e. comparison by pointer value)
96
97 inline void* operator new(std::size_t);
98 inline void operator delete(void* aPrimaryParticle);
99 // Overloaded new/delete operators
100
101 void Print() const;
102 // Print the properties of the particle
103
104
105 // Followings are the available accessors/modifiers.
106 // "trackID" will be set if this particle is sent to G4EventManager
107 // and converted to G4Track. Otherwise = -1.
108 // The mass and charge in G4ParticleDefinition will be used by default.
109 // "SetMass" and "SetCharge" methods are used to set dynamical mass and
110 // charge of G4DynamicParticle.
111 // "GetMass" and "GetCharge" methods will return those in
112 // G4ParticleDefinition if users do not set dynamical ones
113
114 inline G4int GetPDGcode() const;
115 void SetPDGcode(G4int Pcode);
116 inline G4ParticleDefinition* GetG4code() const;
117 inline void SetG4code(const G4ParticleDefinition* Gcode);
118 inline const G4ParticleDefinition* GetParticleDefinition() const;
120 inline G4double GetMass() const;
121 inline void SetMass(G4double mas);
122 inline G4double GetCharge() const;
123 inline void SetCharge(G4double chg);
124 inline G4double GetKineticEnergy() const;
125 inline void SetKineticEnergy(G4double eKin);
126 inline const G4ThreeVector& GetMomentumDirection() const;
127 inline void SetMomentumDirection(const G4ThreeVector& p);
128 inline G4double GetTotalMomentum() const;
129 void Set4Momentum(G4double px, G4double py, G4double pz, G4double E);
130 inline G4double GetTotalEnergy() const;
131 inline void SetTotalEnergy(G4double eTot);
132 inline G4ThreeVector GetMomentum() const;
133 void SetMomentum(G4double px, G4double py, G4double pz);
134 inline G4double GetPx() const;
135 inline G4double GetPy() const;
136 inline G4double GetPz() const;
137 inline G4PrimaryParticle* GetNext() const;
138 inline void SetNext(G4PrimaryParticle* np);
139 inline void ClearNext();
140 inline G4PrimaryParticle* GetDaughter() const;
141 inline void SetDaughter(G4PrimaryParticle* np);
142 inline G4int GetTrackID() const;
143 inline void SetTrackID(G4int id);
144 inline G4ThreeVector GetPolarization() const;
145 inline void SetPolarization(const G4ThreeVector& pol);
146 inline void SetPolarization(G4double px, G4double py, G4double pz);
147 inline G4double GetPolX() const;
148 inline G4double GetPolY() const;
149 inline G4double GetPolZ() const;
150 inline G4double GetWeight() const;
151 inline void SetWeight(G4double w);
152 inline G4double GetProperTime() const;
153 inline void SetProperTime(G4double t);
156
157 private:
158
159 const G4ParticleDefinition* G4code = nullptr;
160
163
166
167 G4double mass = -1.0;
175
177 G4int trackID = -1; // This will be set if this particle is
178 // sent to G4EventManager and converted to
179 // G4Track. Otherwise = -1
180};
181
183
184// ------------------------
185// Inline methods
186// ------------------------
187
188inline
189void* G4PrimaryParticle::operator new(std::size_t)
190{
191 if (aPrimaryParticleAllocator() == nullptr)
192 {
194 }
195 return (void*) aPrimaryParticleAllocator()->MallocSingle();
196}
197
198inline
199void G4PrimaryParticle::operator delete(void* aPrimaryParticle)
200{
202 ->FreeSingle((G4PrimaryParticle*) aPrimaryParticle);
203}
204
205inline
207{
208 return mass;
209}
210
211inline
213{
214 return charge;
215}
216
217inline
219{
220 return PDGcode;
221}
222
223inline
225{
226 return const_cast<G4ParticleDefinition*>(G4code);
227}
228
229inline
231{
232 return G4code;
233}
234
235inline
237{
238 if (mass<0.) return kinE;
239 else return std::sqrt(kinE*(kinE+2.*mass));
240}
241
242inline
244{
246}
247
248inline
250{
251 return direction;
252}
253
254inline
256{
257 direction = p;
258}
259
260inline
262{
263 return GetTotalMomentum()*direction.x();
264}
265
266inline
268{
269 return GetTotalMomentum()*direction.y();
270}
271
272inline
274{
275 return GetTotalMomentum()*direction.z();
276}
277
279{
280 if (mass<0.) return kinE;
281 else return kinE+mass;
282}
283
284inline
286{
287 if (mass<0.) kinE = eTot;
288 else kinE = eTot - mass;
289}
290
291inline
293{
294 return kinE;
295}
296
297inline
299{
300 kinE = eKin;
301}
302
303inline
305{
306 return nextParticle;
307}
308
309inline
311{
312 return daughterParticle;
313}
314
315inline
317{
318 return trackID;
319}
320
321inline
323{
324 return G4ThreeVector(polX,polY,polZ);
325}
326
327inline
329{
330 return polX;
331}
332
333inline
335{
336 return polY;
337}
338
339inline
341{
342 return polZ;
343}
344
345inline
347{
348 return Weight0;
349}
350
351inline
353{
354 Weight0 = w;
355}
356
357inline
359{
360 properTime = t;
361}
362
363inline
365{
366 return properTime;
367}
368
369inline
372{
373 userInfo = anInfo;
374}
375
376inline
379{
380 return userInfo;
381}
382
383inline
385{
387}
388
389inline
391{
392 if (nextParticle == nullptr) { nextParticle = np; }
393 else { nextParticle->SetNext(np); }
394}
395
396inline
398{
399 nextParticle = nullptr;
400}
401
402inline
404{
405 if(daughterParticle == nullptr) { daughterParticle = np; }
406 else { daughterParticle->SetNext(np); }
407}
408
409inline
411{
412 trackID = id;
413}
414
415inline
417{
418 mass = mas;
419}
420
421inline
423{
424 charge = chg;
425}
426
427inline
429{
430 polX = px;
431 polY = py;
432 polZ = pz;
433}
434
435inline
437{
438 polX = pol.x();
439 polY = pol.y();
440 polZ = pol.z();
441}
442
443#endif
G4PART_DLL G4Allocator< G4PrimaryParticle > *& aPrimaryParticleAllocator()
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
double z() const
double x() const
double y() const
void SetPDGcode(G4int Pcode)
G4double GetTotalEnergy() const
G4double GetWeight() const
void SetTotalEnergy(G4double eTot)
G4double GetCharge() const
void SetG4code(const G4ParticleDefinition *Gcode)
G4double GetKineticEnergy() const
G4bool operator==(const G4PrimaryParticle &right) const
void SetCharge(G4double chg)
G4VUserPrimaryParticleInformation * userInfo
G4VUserPrimaryParticleInformation * GetUserInformation() const
G4double GetProperTime() const
const G4ThreeVector & GetMomentumDirection() const
G4PrimaryParticle * daughterParticle
G4double GetPolY() const
const G4ParticleDefinition * G4code
void SetPolarization(const G4ThreeVector &pol)
void SetTrackID(G4int id)
void SetNext(G4PrimaryParticle *np)
G4ThreeVector GetPolarization() const
G4PrimaryParticle * GetNext() const
G4PrimaryParticle & operator=(const G4PrimaryParticle &right)
void SetKineticEnergy(G4double eKin)
G4bool operator!=(const G4PrimaryParticle &right) const
void SetWeight(G4double w)
void Set4Momentum(G4double px, G4double py, G4double pz, G4double E)
G4double GetPy() const
void SetMomentum(G4double px, G4double py, G4double pz)
G4int GetTrackID() const
G4double GetMass() const
G4ThreeVector direction
G4double GetTotalMomentum() const
void SetMomentumDirection(const G4ThreeVector &p)
G4ThreeVector GetMomentum() const
G4int GetPDGcode() const
void SetMass(G4double mas)
G4double GetPz() const
void SetProperTime(G4double t)
void SetUserInformation(G4VUserPrimaryParticleInformation *anInfo)
G4double GetPx() const
G4double GetPolZ() const
void SetParticleDefinition(const G4ParticleDefinition *pdef)
void SetDaughter(G4PrimaryParticle *np)
G4double GetPolX() const
G4PrimaryParticle * GetDaughter() const
const G4ParticleDefinition * GetParticleDefinition() const
G4PrimaryParticle * nextParticle
G4ParticleDefinition * GetG4code() const
#define G4PART_DLL
Definition: pwdefs.hh:45