Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 //
27 // $Id$
28 //
29 //
30 
31 
32 #ifndef G4PrimaryParticle_h
33 #define G4PrimaryParticle_h 1
34 
35 #include "globals.hh"
36 #include "G4Allocator.hh"
37 #include "G4ThreeVector.hh"
38 
41 
42 // class description:
43 //
44 // This is a class which represents a primary particle.
45 // This is completely deferent class from G4Track or G4DynamicParticle.
46 // This class is designed with taking into account the possibility of
47 // making this class persistent, i.e. kept with G4Event class object
48 // to ODBMS. Thus this class is almost free from any other Geant4 classes.
49 // The only exception is a pointer to G4ParticleDefinition but it can be
50 // rebuilt by the PDGcode.
51 //
52 // Primary particles are stored in G4PrimaryVertex object with a form
53 // of linked list. Also, an object of this PrimaryParticle class can have
54 // one or more objects of this class as its daughters with a form of
55 // linked list.
56 // A parimary particle represented by this class object needs not to be
57 // a particle of type which Geant4 can simulate.
58 // case a) mother particle is not a particle Geant4 can simulate
59 // daughters associated to the mother will be examined.
60 // case b) mother particle is a perticle Geant4 can simulate
61 // daughters associated to the mother will be converted to G4Dynamic
62 // particle and be set to the mother G4Track. For this case, dauthers
63 // are used as the "pre-fixed" decay channel.
64 //
65 
67 {
68  public:
69  inline void *operator new(size_t);
70  inline void operator delete(void *aStackedTrack);
71 
72  public: // with description
74  G4PrimaryParticle(G4int Pcode);
76  G4double px,G4double py,G4double pz);
78  G4double px,G4double py,G4double pz,G4double E);
81  G4double px,G4double py,G4double pz);
83  G4double px,G4double py,G4double pz,G4double E);
84 
85  public:
86  virtual ~G4PrimaryParticle();
87 
88  // copy constructor and assignment operator
89  // nextParticle and daughterParticle is copied by object (i.e. deep copy)
90  // userInfo will nt be copied
93 
94  // equal operator returns 'true' only if the same object (i.e comarison by pointer value)
95  G4int operator==(const G4PrimaryParticle &right) const;
96  G4int operator!=(const G4PrimaryParticle &right) const;
97 
98  public: // with description
99  void Print() const;
100  // Print the properties of the particle.
101 
102 
103  public: // with description
104  // followings are get/set methods available.
105  // "trackID" will be set if this particle is sent to G4EventManager
106  // and converted to G4Track. Otherwise = -1.
107  // The mass and charge in G4ParticleDefinition will be used in default.
108  // "SetMass" and "SetCharge" methods are used to set dynamical mass and charge
109  // G4DynamicParticle."GetMass" and "GetCharge" methods will return
110  // those in G4ParticleDefinition if users do not set dynamical ones.
111 
112  G4int GetPDGcode() const;
113  void SetPDGcode(G4int Pcode);
115  void SetG4code(const G4ParticleDefinition * Gcode);
117  void SetParticleDefinition(const G4ParticleDefinition * pdef);
118  G4double GetMass() const;
119  void SetMass(G4double mas);
120  G4double GetCharge() const;
121  void SetCharge(G4double chg);
122  G4double GetKineticEnergy() const;
123  void SetKineticEnergy(G4double eKin);
124  const G4ThreeVector& GetMomentumDirection() const;
125  void SetMomentumDirection(const G4ThreeVector& p);
126  G4double GetTotalMomentum() const;
127  void Set4Momentum(G4double px, G4double py, G4double pz, G4double E);
128  G4double GetTotalEnergy() const;
129  void SetTotalEnergy(G4double eTot );
130  G4ThreeVector GetMomentum() const;
131  void SetMomentum(G4double px, G4double py, G4double pz);
132  G4double GetPx() const;
133  G4double GetPy() const;
134  G4double GetPz() const;
135  G4PrimaryParticle * GetNext() const;
136  void SetNext(G4PrimaryParticle * np);
137  G4PrimaryParticle * GetDaughter() const;
138  void SetDaughter(G4PrimaryParticle * np);
139  G4int GetTrackID() const;
140  void SetTrackID(G4int id);
142  void SetPolarization(const G4ThreeVector& pol);
143  void SetPolarization(G4double px,G4double py,G4double pz);
144  G4double GetPolX() const;
145  G4double GetPolY() const;
146  G4double GetPolZ() const;
147  G4double GetWeight() const;
148  void SetWeight(G4double w);
149  G4double GetProperTime() const;
150  void SetProperTime(G4double t);
153 
154  private:
155  G4int PDGcode;
156  const G4ParticleDefinition * G4code;
157 
158  G4ThreeVector direction;
159  G4double kinE;
160 
161  G4PrimaryParticle * nextParticle;
162  G4PrimaryParticle * daughterParticle;
163 
164  G4int trackID; // This will be set if this particle is
165  // sent to G4EventManager and converted to
166  // G4Track. Otherwise = -1.
167 
168  G4double mass;
169  G4double charge;
170  G4double polX;
171  G4double polY;
172  G4double polZ;
173  G4double Weight0;
174  G4double properTime;
176 
177 };
178 
179 #if defined G4PARTICLES_ALLOC_EXPORT
181 #else
183 #endif
184 
185 inline void * G4PrimaryParticle::operator new(size_t)
186 {
187  void * aPrimaryParticle;
188  aPrimaryParticle = (void *) aPrimaryParticleAllocator.MallocSingle();
189  return aPrimaryParticle;
190 }
191 
192 inline void G4PrimaryParticle::operator delete(void * aPrimaryParticle)
193 {
194  aPrimaryParticleAllocator.FreeSingle((G4PrimaryParticle *) aPrimaryParticle);
195 }
196 
198 { return mass; }
199 
201 { return charge; }
202 
204 { return PDGcode; }
205 
207 { return const_cast<G4ParticleDefinition*>(G4code); }
208 
210 { return G4code; }
211 
213 {
214  if (mass<0.) return kinE;
215  else return std::sqrt(kinE*(kinE+2.*mass));
216 }
217 
219 { return GetTotalMomentum()*direction;}
220 
222 { return direction;}
223 
225 { direction = p;}
226 
228 { return GetTotalMomentum()*direction.x();}
229 
231 { return GetTotalMomentum()*direction.y();}
232 
234 { return GetTotalMomentum()*direction.z();}
235 
237 {
238  if (mass<0.) return kinE;
239  else return kinE+mass;
240 }
241 
243 {
244  if (mass<0.) kinE = eTot;
245  else kinE = eTot - mass;
246 }
247 
249 { return kinE; }
250 
252 { kinE = eKin; }
253 
255 { return nextParticle; }
256 
258 { return daughterParticle; }
259 
261 { return trackID; }
262 
264 { return G4ThreeVector(polX,polY,polZ); }
265 
267 { return polX; }
268 
270 { return polY; }
271 
273 { return polZ; }
274 
276 { return Weight0; }
277 
279 { Weight0 = w; }
280 
282 { properTime = t; }
283 
285 { return properTime; }
286 
288 { userInfo = anInfo; }
289 
291 { return userInfo; }
292 
294 {
295  SetParticleDefinition(Gcode);
296 }
297 
299 {
300  if (nextParticle == 0) { nextParticle = np; }
301  else { nextParticle->SetNext(np); }
302 }
303 
305 {
306  if(daughterParticle == 0) { daughterParticle = np; }
307  else { daughterParticle->SetNext(np); }
308 }
309 
311 { trackID = id; }
312 
314 { mass = mas; }
315 
317 { charge = chg; }
318 
320 {
321  polX = px;
322  polY = py;
323  polZ = pz;
324 }
325 
327 {
328  polX = pol.x();
329  polY = pol.y();
330  polZ = pol.z();
331 }
332 
333 #endif
334