Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4PrimaryParticle.hh
이 파일의 문서화 페이지로 가기
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: G4PrimaryParticle.hh 110257 2018-05-17 14:20:12Z gcosmo $
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 
39 #include "pwdefs.hh"
40 
43 
44 // class description:
45 //
46 // This is a class which represents a primary particle.
47 // This is completely deferent class from G4Track or G4DynamicParticle.
48 // This class is designed with taking into account the possibility of
49 // making this class persistent, i.e. kept with G4Event class object
50 // to ODBMS. Thus this class is almost free from any other Geant4 classes.
51 // The only exception is a pointer to G4ParticleDefinition but it can be
52 // rebuilt by the PDGcode.
53 //
54 // Primary particles are stored in G4PrimaryVertex object with a form
55 // of linked list. Also, an object of this PrimaryParticle class can have
56 // one or more objects of this class as its daughters with a form of
57 // linked list.
58 // A parimary particle represented by this class object needs not to be
59 // a particle of type which Geant4 can simulate.
60 // case a) mother particle is not a particle Geant4 can simulate
61 // daughters associated to the mother will be examined.
62 // case b) mother particle is a perticle Geant4 can simulate
63 // daughters associated to the mother will be converted to G4Dynamic
64 // particle and be set to the mother G4Track. For this case, dauthers
65 // are used as the "pre-fixed" decay channel.
66 //
67 
69 {
70  public:
71  inline void *operator new(size_t);
72  inline void operator delete(void *aStackedTrack);
73 
74  public: // with description
76  G4PrimaryParticle(G4int Pcode);
78  G4double px,G4double py,G4double pz);
83  G4double px,G4double py,G4double pz);
85  G4double px,G4double py,G4double pz,G4double E);
86 
87  public:
88  virtual ~G4PrimaryParticle();
89 
90  // copy constructor and assignment operator
91  // nextParticle and daughterParticle is copied by object (i.e. deep copy)
92  // userInfo will nt be copied
95 
96  // equal operator returns 'true' only if the same object (i.e comarison by pointer value)
97  G4int operator==(const G4PrimaryParticle &right) const;
98  G4int operator!=(const G4PrimaryParticle &right) const;
99 
100  public: // with description
101  void Print() const;
102  // Print the properties of the particle.
103 
104 
105  public: // with description
106  // followings are get/set methods available.
107  // "trackID" will be set if this particle is sent to G4EventManager
108  // and converted to G4Track. Otherwise = -1.
109  // The mass and charge in G4ParticleDefinition will be used in default.
110  // "SetMass" and "SetCharge" methods are used to set dynamical mass and charge
111  // G4DynamicParticle."GetMass" and "GetCharge" methods will return
112  // those in G4ParticleDefinition if users do not set dynamical ones.
113 
114  G4int GetPDGcode() const;
115  void SetPDGcode(G4int Pcode);
117  void SetG4code(const G4ParticleDefinition * Gcode);
119  void SetParticleDefinition(const G4ParticleDefinition * pdef);
120  G4double GetMass() const;
121  void SetMass(G4double mas);
122  G4double GetCharge() const;
123  void SetCharge(G4double chg);
124  G4double GetKineticEnergy() const;
125  void SetKineticEnergy(G4double eKin);
126  const G4ThreeVector& GetMomentumDirection() const;
127  void SetMomentumDirection(const G4ThreeVector& p);
128  G4double GetTotalMomentum() const;
129  void Set4Momentum(G4double px, G4double py, G4double pz, G4double E);
130  G4double GetTotalEnergy() const;
131  void SetTotalEnergy(G4double eTot );
132  G4ThreeVector GetMomentum() const;
133  void SetMomentum(G4double px, G4double py, G4double pz);
134  G4double GetPx() const;
135  G4double GetPy() const;
136  G4double GetPz() const;
137  G4PrimaryParticle * GetNext() const;
138  void SetNext(G4PrimaryParticle * np);
139  void ClearNext();
140  G4PrimaryParticle * GetDaughter() const;
141  void SetDaughter(G4PrimaryParticle * np);
142  G4int GetTrackID() const;
143  void SetTrackID(G4int id);
145  void SetPolarization(const G4ThreeVector& pol);
146  void SetPolarization(G4double px,G4double py,G4double pz);
147  G4double GetPolX() const;
148  G4double GetPolY() const;
149  G4double GetPolZ() const;
150  G4double GetWeight() const;
151  void SetWeight(G4double w);
152  G4double GetProperTime() const;
153  void SetProperTime(G4double t);
156 
157  private:
160 
163 
166 
167  G4int trackID; // This will be set if this particle is
168  // sent to G4EventManager and converted to
169  // G4Track. Otherwise = -1.
170 
179 
180 };
181 
183 
184 inline void * G4PrimaryParticle::operator new(size_t)
185 {
187  {
189  }
190  return (void *) aPrimaryParticleAllocator()->MallocSingle();
191 }
192 
193 inline void G4PrimaryParticle::operator delete(void * aPrimaryParticle)
194 {
195  aPrimaryParticleAllocator()->FreeSingle((G4PrimaryParticle *) aPrimaryParticle);
196 }
197 
199 { return mass; }
200 
202 { return charge; }
203 
205 { return PDGcode; }
206 
208 { return const_cast<G4ParticleDefinition*>(G4code); }
209 
211 { return G4code; }
212 
214 {
215  if (mass<0.) return kinE;
216  else return std::sqrt(kinE*(kinE+2.*mass));
217 }
218 
220 { return GetTotalMomentum()*direction;}
221 
223 { return direction;}
224 
226 { direction = p;}
227 
229 { return GetTotalMomentum()*direction.x();}
230 
232 { return GetTotalMomentum()*direction.y();}
233 
235 { return GetTotalMomentum()*direction.z();}
236 
238 {
239  if (mass<0.) return kinE;
240  else return kinE+mass;
241 }
242 
244 {
245  if (mass<0.) kinE = eTot;
246  else kinE = eTot - mass;
247 }
248 
250 { return kinE; }
251 
253 { kinE = eKin; }
254 
256 { return nextParticle; }
257 
259 { return daughterParticle; }
260 
262 { return trackID; }
263 
265 { return G4ThreeVector(polX,polY,polZ); }
266 
268 { return polX; }
269 
271 { return polY; }
272 
274 { return polZ; }
275 
277 { return Weight0; }
278 
280 { Weight0 = w; }
281 
283 { properTime = t; }
284 
286 { return properTime; }
287 
289 { userInfo = anInfo; }
290 
292 { return userInfo; }
293 
295 {
296  SetParticleDefinition(Gcode);
297 }
298 
300 {
301  if (nextParticle == 0) { nextParticle = np; }
302  else { nextParticle->SetNext(np); }
303 }
304 
306 {
307  nextParticle = nullptr;
308 }
309 
311 {
312  if(daughterParticle == 0) { daughterParticle = np; }
313  else { daughterParticle->SetNext(np); }
314 }
315 
317 { trackID = id; }
318 
320 { mass = mas; }
321 
323 { charge = chg; }
324 
326 {
327  polX = px;
328  polY = py;
329  polZ = pz;
330 }
331 
333 {
334  polX = pol.x();
335  polY = pol.y();
336  polZ = pol.z();
337 }
338 
339 #endif
G4ThreeVector GetMomentum() const
CLHEP::Hep3Vector G4ThreeVector
G4PART_DLL G4Allocator< G4PrimaryParticle > *& aPrimaryParticleAllocator()
G4int operator!=(const G4PrimaryParticle &right) const
G4PrimaryParticle & operator=(const G4PrimaryParticle &right)
const char * p
Definition: xmltok.h:285
G4double GetMass() const
G4double GetPolY() const
G4int GetTrackID() const
G4double GetPx() const
void SetUserInformation(G4VUserPrimaryParticleInformation *anInfo)
double z() const
void SetKineticEnergy(G4double eKin)
void SetMomentum(G4double px, G4double py, G4double pz)
G4ThreeVector direction
void SetWeight(G4double w)
G4PrimaryParticle * GetDaughter() const
G4double GetTotalMomentum() const
void SetDaughter(G4PrimaryParticle *np)
void SetTrackID(G4int id)
G4double GetWeight() const
G4double GetPy() const
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
G4PrimaryParticle * daughterParticle
G4int operator==(const G4PrimaryParticle &right) const
G4double GetPolX() const
G4double GetProperTime() const
G4ParticleDefinition * GetG4code() const
void SetParticleDefinition(const G4ParticleDefinition *pdef)
void SetPDGcode(G4int Pcode)
G4double GetPz() const
const G4ParticleDefinition * GetParticleDefinition() const
void SetPolarization(const G4ThreeVector &pol)
void SetTotalEnergy(G4double eTot)
void SetCharge(G4double chg)
G4ThreeVector GetPolarization() const
G4int GetPDGcode() const
void Set4Momentum(G4double px, G4double py, G4double pz, G4double E)
G4PrimaryParticle * nextParticle
G4PrimaryParticle * GetNext() const
G4double GetPolZ() const
int G4int
Definition: G4Types.hh:78
void SetMomentumDirection(const G4ThreeVector &p)
double x() const
const G4ThreeVector & GetMomentumDirection() const
void SetG4code(const G4ParticleDefinition *Gcode)
G4double GetCharge() const
G4double GetTotalEnergy() const
const G4ParticleDefinition * G4code
#define G4PART_DLL
Definition: pwdefs.hh:48
void SetProperTime(G4double t)
double y() const
void SetNext(G4PrimaryParticle *np)
void SetMass(G4double mas)
G4VUserPrimaryParticleInformation * userInfo
G4VUserPrimaryParticleInformation * GetUserInformation() const