Geant4-11
G4KineticTrack.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//
28// -----------------------------------------------------------------------------
29// GEANT 4 class header file
30//
31// History: first implementation, A. Feliciello, 20th May 1998
32// -----------------------------------------------------------------------------
33
34#ifndef G4KineticTrack_h
35#define G4KineticTrack_h 1
36
38
39#include "globals.hh"
40#include "G4ios.hh"
41
42
43#include "Randomize.hh"
44#include "G4ThreeVector.hh"
45#include "G4LorentzVector.hh"
46#include "G4VKineticNucleon.hh"
47#include "G4Nucleon.hh"
49#include "G4VDecayChannel.hh"
50#include "G4Log.hh"
51
52// #include "G4Allocator.hh"
53
55
57{
58 public:
59
61
62 G4KineticTrack(const G4KineticTrack& right);
63
64 G4KineticTrack(const G4ParticleDefinition* aDefinition,
65 G4double aFormationTime,
66 const G4ThreeVector& aPosition,
67 const G4LorentzVector& a4Momentum);
69 const G4ThreeVector& aPosition,
70 const G4LorentzVector& a4Momentum);
71
73
75
76 G4bool operator==(const G4KineticTrack& right) const;
77
78 G4bool operator!=(const G4KineticTrack& right) const;
79/*
80 inline void *operator new(size_t);
81 inline void operator delete(void *aTrack);
82*/
84 void SetDefinition(const G4ParticleDefinition* aDefinition);
85
87 void SetFormationTime(G4double aFormationTime);
88
89 const G4ThreeVector& GetPosition() const;
90 void SetPosition(const G4ThreeVector aPosition);
91
92 const G4LorentzVector& Get4Momentum() const;
93 void Set4Momentum(const G4LorentzVector& a4Momentum);
94 void Update4Momentum(G4double aEnergy); // update E and p, not changing mass
95 void Update4Momentum(const G4ThreeVector & aMomentum); // idem
96 void SetTrackingMomentum(const G4LorentzVector& a4Momentum);
97 void UpdateTrackingMomentum(G4double aEnergy); // update E and p, not changing mass
98 void UpdateTrackingMomentum(const G4ThreeVector & aMomentum); // idem
99
101
103
104 void Hit();
105 void SetNucleon(G4Nucleon * aN) {theNucleon = aN;}
106
107 G4bool IsParticipant() const;
108
110
111 // LB move to public (before was private) LB
112 G4double* GetActualWidth() const;
113
114 G4double GetActualMass() const;
115 G4int GetnChannels() const;
116
117// position relativ to nucleus "state"
120
121 CascadeState SetState(const CascadeState new_state);
122 CascadeState GetState() const;
123 void SetProjectilePotential(const G4double aPotential);
125
126 void SetCreatorModelID(G4int id);
127 G4int GetCreatorModelID() const;
128
129 private:
130
131
132 void SetnChannels(const G4int aChannel);
133
134 void SetActualWidth(G4double* anActualWidth);
135
137
139 const G4double* m_ij) const;
140
141 G4double IntegrateCMMomentum(const G4double lowerLimit) const;
142
143 G4double IntegrateCMMomentum(const G4double lowerLimit ,const G4double polemass) const;
144
146
147 public:
148
149 G4double BrWig(const G4double Gamma,
150 const G4double rmass,
151 const G4double mass) const;
152
153private:
158public:
159 // friend G4double IntegrandFunction3 (G4double xmass);
160
161 // friend G4double IntegrandFunction4 (G4double xmass);
162
163 private:
164
166
168
170
174
176
178
180
182
183 // Temporary storage for daughter masses and widths
184 // (needed because Integrand Function cannot take > 1 argument)
187
189
191
193};
194
195// extern G4Allocator<G4KineticTrack> theKTAllocator;
196
197
198// Class G4KineticTrack
199/*
200inline void * G4KineticTrack::operator new(size_t)
201{
202 void * aT;
203 aT = (void *) theKTAllocator.MallocSingle();
204 return aT;
205}
206
207inline void G4KineticTrack::operator delete(void * aT)
208{
209 theKTAllocator.FreeSingle((G4KineticTrack *) aT);
210}
211*/
212
214{
215 return theDefinition;
216}
217
219{
220 theDefinition = aDefinition;
221}
222
224{
225 return theFormationTime;
226}
227
229{
230 theFormationTime = aFormationTime;
231}
232
234{
235 return thePosition;
236}
237
238inline void G4KineticTrack::SetPosition(const G4ThreeVector aPosition)
239{
240 thePosition = aPosition;
241}
242
244{
245 return theTotal4Momentum;
246}
247
249{
250 return the4Momentum;
251}
252
253inline void G4KineticTrack::Set4Momentum(const G4LorentzVector& a4Momentum)
254{
255// set the4Momentum and update theTotal4Momentum
256
257 theTotal4Momentum=a4Momentum;
260}
261
263{
264// update the4Momentum with aEnergy at constant mass (the4Momentum.mag()
265// updates theTotal4Momentum as well.
266 G4double newP(0);
268 if ( sqr(aEnergy) > mass2 )
269 {
270 newP = std::sqrt(sqr(aEnergy) - mass2 );
271 } else
272 {
273 aEnergy=std::sqrt(mass2);
274 }
276}
277
278inline void G4KineticTrack::Update4Momentum(const G4ThreeVector & aMomentum)
279{
280// update the4Momentum with aMomentum at constant mass (the4Momentum.mag()
281// updates theTotal4Momentum as well.
282 G4double newE=std::sqrt(theTotal4Momentum.mag2() + aMomentum.mag2());
283 Set4Momentum(G4LorentzVector(aMomentum, newE));
284}
285
287{
288// set the4Momentum and update theTotal4Momentum, keep the mass of aMomentum
289
290 the4Momentum = aMomentum;
292// keep mass of aMomentum for the total momentum
293 G4double mass2 = aMomentum.mag2();
295 theTotal4Momentum.setE(std::sqrt(mass2+p2));
296}
297
299{
300// update the4Momentum with aEnergy at constant mass (the4Momentum.mag()
301// updates theTotal4Momentum as well.
302 G4double newP(0);
304 if ( sqr(aEnergy) > mass2 )
305 {
306 newP = std::sqrt(sqr(aEnergy) - mass2 );
307 } else
308 {
309 aEnergy=std::sqrt(mass2);
310 }
312}
313
315{
316// update the4Momentum with aMomentum at constant mass (the4Momentum.mag()
317// updates theTotal4Momentum as well.
318 G4double newE=std::sqrt(theTotal4Momentum.mag2() + aMomentum.mag2());
319 SetTrackingMomentum(G4LorentzVector(aMomentum, newE));
320}
321
323{
324 return std::sqrt(std::abs(the4Momentum.mag2()));
325}
326
328{
329 return nChannels;
330}
331
332inline void G4KineticTrack::SetnChannels(const G4int numberOfChannels)
333{
334 nChannels = numberOfChannels;
335}
336
338{
339 return theActualWidth;
340}
341
342inline void G4KineticTrack::SetActualWidth(G4double* anActualWidth)
343{
344 theActualWidth = anActualWidth;
345}
346
347
348
350{
351 G4int index;
352 G4double theTotalActualWidth = 0.0;
353 for (index = nChannels - 1; index >= 0; index--)
354 {
355 theTotalActualWidth += theActualWidth[index];
356 }
357 return theTotalActualWidth;
358}
359
361{
362 G4double theTotalActualWidth = this->EvaluateTotalActualWidth();
363 G4double tau = CLHEP::hbar_Planck * (-1.0 / theTotalActualWidth);
364 G4double theResidualLifetime = tau * G4Log(G4UniformRand());
365 return theResidualLifetime*the4Momentum.gamma();
366}
367
369 const G4double* m_ij) const
370{
371 G4double theCMMomentum;
372 if((m_ij[0]+m_ij[1])<mass)
373 theCMMomentum = 1 / (2 * mass) *
374 std::sqrt (((mass * mass) - (m_ij[0] + m_ij[1]) * (m_ij[0] + m_ij[1])) *
375 ((mass * mass) - (m_ij[0] - m_ij[1]) * (m_ij[0] - m_ij[1])));
376 else
377 theCMMomentum=0.;
378
379 return theCMMomentum;
380}
381
382inline G4double G4KineticTrack::BrWig(const G4double Gamma, const G4double rmass, const G4double mass) const
383{
384 G4double Norm = CLHEP::twopi;
385 return (Gamma/((mass-rmass)*(mass-rmass)+Gamma*Gamma/4.))/Norm;
386}
387
388inline
390{
391 if(theNucleon)
392 {
393 theNucleon->Hit(1);
394 }
395}
396
397inline
399{
400 if(!theNucleon) return true;
401 return theNucleon->AreYouHit();
402}
403
404inline
406{
407 return theStateToNucleus;
408}
409
410inline
412{
414 theStateToNucleus=new_state;
415 return old_state;
416}
417
418inline
420{
421 theProjectilePotential = aPotential;
422}
423inline
425{
427}
428
429inline
431{
432 theCreatorModel = id;
433}
434inline
436{
437 return theCreatorModel;
438}
439
440#endif
441
442
443
G4double G4Log(G4double x)
Definition: G4Log.hh:226
CLHEP::HepLorentzVector G4LorentzVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
double mag2() const
Hep3Vector vect() const
G4double GetFormationTime() const
CascadeState SetState(const CascadeState new_state)
G4double * theDaughterMass
G4double GetProjectilePotential() const
G4double IntegrandFunction1(G4double xmass) const
void SetPosition(const G4ThreeVector aPosition)
G4double IntegrateCMMomentum2() const
G4KineticTrackVector * Decay()
G4KineticTrack & operator=(const G4KineticTrack &right)
void SetCreatorModelID(G4int id)
void SetnChannels(const G4int aChannel)
G4int GetCreatorModelID() const
G4LorentzVector the4Momentum
CascadeState GetState() const
void SetNucleon(G4Nucleon *aN)
G4bool operator!=(const G4KineticTrack &right) const
void SetProjectilePotential(const G4double aPotential)
G4ThreeVector thePosition
void Set4Momentum(const G4LorentzVector &a4Momentum)
G4LorentzVector theTotal4Momentum
G4double * theActualWidth
G4int GetnChannels() const
const G4ThreeVector & GetPosition() const
G4Nucleon * theNucleon
void SetTrackingMomentum(const G4LorentzVector &a4Momentum)
const G4ParticleDefinition * GetDefinition() const
G4bool operator==(const G4KineticTrack &right) const
G4bool IsParticipant() const
G4double theFormationTime
void SetDefinition(const G4ParticleDefinition *aDefinition)
const G4LorentzVector & GetTrackingMomentum() const
G4double * theDaughterWidth
G4double theProjectilePotential
G4double * GetActualWidth() const
G4double IntegrateCMMomentum(const G4double lowerLimit) const
G4double IntegrandFunction3(G4double xmass) const
void UpdateTrackingMomentum(G4double aEnergy)
G4double IntegrandFunction4(G4double xmass) const
G4double IntegrandFunction2(G4double xmass) const
void SetFormationTime(G4double aFormationTime)
void Update4Momentum(G4double aEnergy)
G4double SampleResidualLifetime()
G4double EvaluateCMMomentum(const G4double mass, const G4double *m_ij) const
G4double BrWig(const G4double Gamma, const G4double rmass, const G4double mass) const
G4double theActualMass
void SetActualWidth(G4double *anActualWidth)
CascadeState theStateToNucleus
G4LorentzVector theFermi3Momentum
const G4LorentzVector & Get4Momentum() const
const G4ParticleDefinition * theDefinition
G4double GetActualMass() const
G4double EvaluateTotalActualWidth()
G4bool AreYouHit() const
Definition: G4Nucleon.hh:98
void Hit(G4VSplitableHadron *aHit)
Definition: G4Nucleon.hh:91
static constexpr double hbar_Planck
static constexpr double twopi
Definition: SystemOfUnits.h:56
G4bool nucleon(G4int ityp)
T sqr(const T &x)
Definition: templates.hh:128