File indexing completed on 2025-01-18 09:58:58
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054 #ifndef G4PrimaryParticle_hh
0055 #define G4PrimaryParticle_hh 1
0056
0057 #include "G4Allocator.hh"
0058 #include "G4ThreeVector.hh"
0059 #include "globals.hh"
0060
0061 #include "pwdefs.hh"
0062
0063 class G4ParticleDefinition;
0064 class G4VUserPrimaryParticleInformation;
0065
0066 class G4PrimaryParticle
0067 {
0068 public:
0069
0070 G4PrimaryParticle();
0071 G4PrimaryParticle(G4int Pcode);
0072 G4PrimaryParticle(G4int Pcode, G4double px, G4double py, G4double pz);
0073 G4PrimaryParticle(G4int Pcode, G4double px, G4double py, G4double pz, G4double E);
0074 G4PrimaryParticle(const G4ParticleDefinition* Gcode);
0075 G4PrimaryParticle(const G4ParticleDefinition* Gcode, G4double px, G4double py, G4double pz);
0076 G4PrimaryParticle(const G4ParticleDefinition* Gcode, G4double px, G4double py, G4double pz,
0077 G4double E);
0078
0079
0080 virtual ~G4PrimaryParticle();
0081
0082
0083
0084
0085 G4PrimaryParticle(const G4PrimaryParticle& right);
0086 G4PrimaryParticle& operator=(const G4PrimaryParticle& right);
0087
0088
0089
0090 G4bool operator==(const G4PrimaryParticle& right) const;
0091 G4bool operator!=(const G4PrimaryParticle& right) const;
0092
0093
0094 inline void* operator new(std::size_t);
0095 inline void operator delete(void* aPrimaryParticle);
0096
0097
0098 void Print() const;
0099
0100
0101
0102
0103
0104
0105
0106
0107
0108
0109 inline G4int GetPDGcode() const;
0110 void SetPDGcode(G4int Pcode);
0111 inline G4ParticleDefinition* GetG4code() const;
0112 inline void SetG4code(const G4ParticleDefinition* Gcode);
0113 inline const G4ParticleDefinition* GetParticleDefinition() const;
0114 void SetParticleDefinition(const G4ParticleDefinition* pdef);
0115 inline G4double GetMass() const;
0116 inline void SetMass(G4double mas);
0117 inline G4double GetCharge() const;
0118 inline void SetCharge(G4double chg);
0119 inline G4double GetKineticEnergy() const;
0120 inline void SetKineticEnergy(G4double eKin);
0121 inline const G4ThreeVector& GetMomentumDirection() const;
0122 inline void SetMomentumDirection(const G4ThreeVector& p);
0123 inline G4double GetTotalMomentum() const;
0124 void Set4Momentum(G4double px, G4double py, G4double pz, G4double E);
0125 inline G4double GetTotalEnergy() const;
0126 inline void SetTotalEnergy(G4double eTot);
0127 inline G4ThreeVector GetMomentum() const;
0128 void SetMomentum(G4double px, G4double py, G4double pz);
0129 inline G4double GetPx() const;
0130 inline G4double GetPy() const;
0131 inline G4double GetPz() const;
0132 inline G4PrimaryParticle* GetNext() const;
0133 inline void SetNext(G4PrimaryParticle* np);
0134 inline void ClearNext();
0135 inline G4PrimaryParticle* GetDaughter() const;
0136 inline void SetDaughter(G4PrimaryParticle* np);
0137 inline G4int GetTrackID() const;
0138 inline void SetTrackID(G4int id);
0139 inline G4ThreeVector GetPolarization() const;
0140 inline void SetPolarization(const G4ThreeVector& pol);
0141 inline void SetPolarization(G4double px, G4double py, G4double pz);
0142 inline G4double GetPolX() const;
0143 inline G4double GetPolY() const;
0144 inline G4double GetPolZ() const;
0145 inline G4double GetWeight() const;
0146 inline void SetWeight(G4double w);
0147 inline G4double GetProperTime() const;
0148 inline void SetProperTime(G4double t);
0149 inline G4VUserPrimaryParticleInformation* GetUserInformation() const;
0150 inline void SetUserInformation(G4VUserPrimaryParticleInformation* anInfo);
0151
0152 private:
0153 const G4ParticleDefinition* G4code = nullptr;
0154
0155 G4ThreeVector direction;
0156 G4double kinE = 0.0;
0157
0158 G4PrimaryParticle* nextParticle = nullptr;
0159 G4PrimaryParticle* daughterParticle = nullptr;
0160
0161 G4double mass = -1.0;
0162 G4double charge = 0.0;
0163 G4double polX = 0.0;
0164 G4double polY = 0.0;
0165 G4double polZ = 0.0;
0166 G4double Weight0 = 1.0;
0167 G4double properTime = -1.0;
0168 G4VUserPrimaryParticleInformation* userInfo = nullptr;
0169
0170 G4int PDGcode = 0;
0171 G4int trackID = -1;
0172
0173
0174 };
0175
0176 extern G4PART_DLL G4Allocator<G4PrimaryParticle>*& aPrimaryParticleAllocator();
0177
0178
0179
0180
0181
0182 inline void* G4PrimaryParticle::operator new(std::size_t)
0183 {
0184 if (aPrimaryParticleAllocator() == nullptr) {
0185 aPrimaryParticleAllocator() = new G4Allocator<G4PrimaryParticle>;
0186 }
0187 return (void*)aPrimaryParticleAllocator()->MallocSingle();
0188 }
0189
0190 inline void G4PrimaryParticle::operator delete(void* aPrimaryParticle)
0191 {
0192 aPrimaryParticleAllocator()->FreeSingle((G4PrimaryParticle*)aPrimaryParticle);
0193 }
0194
0195 inline G4double G4PrimaryParticle::GetMass() const
0196 {
0197 return mass;
0198 }
0199
0200 inline G4double G4PrimaryParticle::GetCharge() const
0201 {
0202 return charge;
0203 }
0204
0205 inline G4int G4PrimaryParticle::GetPDGcode() const
0206 {
0207 return PDGcode;
0208 }
0209
0210 inline G4ParticleDefinition* G4PrimaryParticle::GetG4code() const
0211 {
0212 return const_cast<G4ParticleDefinition*>(G4code);
0213 }
0214
0215 inline const G4ParticleDefinition* G4PrimaryParticle::GetParticleDefinition() const
0216 {
0217 return G4code;
0218 }
0219
0220 inline G4double G4PrimaryParticle::GetTotalMomentum() const
0221 {
0222 if (mass < 0.) return kinE;
0223 return std::sqrt(kinE * (kinE + 2. * mass));
0224 }
0225
0226 inline G4ThreeVector G4PrimaryParticle::GetMomentum() const
0227 {
0228 return GetTotalMomentum() * direction;
0229 }
0230
0231 inline const G4ThreeVector& G4PrimaryParticle::GetMomentumDirection() const
0232 {
0233 return direction;
0234 }
0235
0236 inline void G4PrimaryParticle::SetMomentumDirection(const G4ThreeVector& p)
0237 {
0238 direction = p;
0239 }
0240
0241 inline G4double G4PrimaryParticle::GetPx() const
0242 {
0243 return GetTotalMomentum() * direction.x();
0244 }
0245
0246 inline G4double G4PrimaryParticle::GetPy() const
0247 {
0248 return GetTotalMomentum() * direction.y();
0249 }
0250
0251 inline G4double G4PrimaryParticle::GetPz() const
0252 {
0253 return GetTotalMomentum() * direction.z();
0254 }
0255
0256 inline G4double G4PrimaryParticle::GetTotalEnergy() const
0257 {
0258 if (mass < 0.) return kinE;
0259 return kinE + mass;
0260 }
0261
0262 inline void G4PrimaryParticle::SetTotalEnergy(G4double eTot)
0263 {
0264 if (mass < 0.)
0265 kinE = eTot;
0266 else
0267 kinE = eTot - mass;
0268 }
0269
0270 inline G4double G4PrimaryParticle::GetKineticEnergy() const
0271 {
0272 return kinE;
0273 }
0274
0275 inline void G4PrimaryParticle::SetKineticEnergy(G4double eKin)
0276 {
0277 kinE = eKin;
0278 }
0279
0280 inline G4PrimaryParticle* G4PrimaryParticle::GetNext() const
0281 {
0282 return nextParticle;
0283 }
0284
0285 inline G4PrimaryParticle* G4PrimaryParticle::GetDaughter() const
0286 {
0287 return daughterParticle;
0288 }
0289
0290 inline G4int G4PrimaryParticle::GetTrackID() const
0291 {
0292 return trackID;
0293 }
0294
0295 inline G4ThreeVector G4PrimaryParticle::GetPolarization() const
0296 {
0297 return G4ThreeVector(polX, polY, polZ);
0298 }
0299
0300 inline G4double G4PrimaryParticle::GetPolX() const
0301 {
0302 return polX;
0303 }
0304
0305 inline G4double G4PrimaryParticle::GetPolY() const
0306 {
0307 return polY;
0308 }
0309
0310 inline G4double G4PrimaryParticle::GetPolZ() const
0311 {
0312 return polZ;
0313 }
0314
0315 inline G4double G4PrimaryParticle::GetWeight() const
0316 {
0317 return Weight0;
0318 }
0319
0320 inline void G4PrimaryParticle::SetWeight(G4double w)
0321 {
0322 Weight0 = w;
0323 }
0324
0325 inline void G4PrimaryParticle::SetProperTime(G4double t)
0326 {
0327 properTime = t;
0328 }
0329
0330 inline G4double G4PrimaryParticle::GetProperTime() const
0331 {
0332 return properTime;
0333 }
0334
0335 inline void G4PrimaryParticle::SetUserInformation(G4VUserPrimaryParticleInformation* anInfo)
0336 {
0337 userInfo = anInfo;
0338 }
0339
0340 inline G4VUserPrimaryParticleInformation* G4PrimaryParticle::GetUserInformation() const
0341 {
0342 return userInfo;
0343 }
0344
0345 inline void G4PrimaryParticle::SetG4code(const G4ParticleDefinition* Gcode)
0346 {
0347 SetParticleDefinition(Gcode);
0348 }
0349
0350 inline void G4PrimaryParticle::SetNext(G4PrimaryParticle* np)
0351 {
0352 if (nextParticle == nullptr) {
0353 nextParticle = np;
0354 }
0355 else {
0356 nextParticle->SetNext(np);
0357 }
0358 }
0359
0360 inline void G4PrimaryParticle::ClearNext()
0361 {
0362 nextParticle = nullptr;
0363 }
0364
0365 inline void G4PrimaryParticle::SetDaughter(G4PrimaryParticle* np)
0366 {
0367 if (daughterParticle == nullptr) {
0368 daughterParticle = np;
0369 }
0370 else {
0371 daughterParticle->SetNext(np);
0372 }
0373 }
0374
0375 inline void G4PrimaryParticle::SetTrackID(G4int id)
0376 {
0377 trackID = id;
0378 }
0379
0380 inline void G4PrimaryParticle::SetMass(G4double mas)
0381 {
0382 mass = mas;
0383 }
0384
0385 inline void G4PrimaryParticle::SetCharge(G4double chg)
0386 {
0387 charge = chg;
0388 }
0389
0390 inline void G4PrimaryParticle::SetPolarization(G4double px, G4double py, G4double pz)
0391 {
0392 polX = px;
0393 polY = py;
0394 polZ = pz;
0395 }
0396
0397 inline void G4PrimaryParticle::SetPolarization(const G4ThreeVector& pol)
0398 {
0399 polX = pol.x();
0400 polY = pol.y();
0401 polZ = pol.z();
0402 }
0403
0404 #endif