source: trunk/source/processes/hadronic/models/util/include/G4KineticTrack.hh @ 1340

Last change on this file since 1340 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 12.2 KB
Line 
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// $Id: G4KineticTrack.hh,v 1.0 1998/05/20
29// -----------------------------------------------------------------------------
30//      GEANT 4 class header file
31//
32//      History: first implementation, A. Feliciello, 20th May 1998
33// -----------------------------------------------------------------------------
34
35#ifndef G4KineticTrack_h
36#define G4KineticTrack_h 1
37
38#include "globals.hh"
39#include "G4ios.hh"
40
41
42#include "Randomize.hh"
43#include "G4ThreeVector.hh"
44#include "G4LorentzVector.hh"
45#include "G4VKineticNucleon.hh"
46#include "G4Nucleon.hh"
47#include "G4ParticleDefinition.hh"
48#include "G4VDecayChannel.hh"
49
50// #include "G4Allocator.hh"
51
52class G4KineticTrackVector;
53
54
55
56
57
58class G4KineticTrack : public G4VKineticNucleon
59{
60  public:
61     
62      G4KineticTrack();
63
64      G4KineticTrack(const G4KineticTrack& right);
65
66      G4KineticTrack(G4ParticleDefinition* aDefinition,
67                     G4double aFormationTime,
68                     G4ThreeVector aPosition, 
69                     G4LorentzVector& a4Momentum);
70      G4KineticTrack(G4Nucleon * nucleon,
71                     G4ThreeVector aPosition,
72                     G4LorentzVector& a4Momentum);
73
74      ~G4KineticTrack();
75
76      const G4KineticTrack& operator=(const G4KineticTrack& right);
77
78      G4int operator==(const G4KineticTrack& right) const;
79
80      G4int operator!=(const G4KineticTrack& right) const;
81/*
82      inline void *operator new(size_t);
83      inline void operator delete(void *aTrack);
84*/
85      G4ParticleDefinition* GetDefinition() const;
86      void SetDefinition(G4ParticleDefinition* aDefinition);
87
88      G4double GetFormationTime() const;
89      void SetFormationTime(G4double aFormationTime);
90
91      const G4ThreeVector& GetPosition() const;
92      void SetPosition(const G4ThreeVector aPosition);
93     
94      const G4LorentzVector& Get4Momentum() const;
95      void Set4Momentum(const G4LorentzVector& a4Momentum);
96      void Update4Momentum(G4double aEnergy);                   // update E and p, not changing mass
97      void Update4Momentum(const G4ThreeVector & aMomentum);    // idem
98      void SetTrackingMomentum(const G4LorentzVector& a4Momentum);
99      void UpdateTrackingMomentum(G4double aEnergy);                    // update E and p, not changing mass
100      void UpdateTrackingMomentum(const G4ThreeVector & aMomentum);     // idem
101
102      const G4LorentzVector& GetTrackingMomentum() const;
103     
104      G4double SampleResidualLifetime();
105     
106      void Hit(); 
107      void SetNucleon(G4Nucleon * aN) {theNucleon = aN;}
108           
109      G4bool IsParticipant() const; 
110
111      G4KineticTrackVector* Decay();
112     
113  // LB move to public (before was private) LB
114      G4double* GetActualWidth() const;
115
116      G4double GetActualMass() const;
117      G4int GetnChannels() const;
118     
119//   position relativ to nucleus "state"
120      enum CascadeState {undefined, outside, going_in, inside, 
121                         going_out, gone_out, captured, miss_nucleus };
122     
123      CascadeState SetState(const CascadeState new_state);
124      CascadeState GetState() const;
125      void SetProjectilePotential(const G4double aPotential);
126      G4double GetProjectilePotential() const;
127
128     
129  private:
130
131
132      void SetnChannels(const G4int aChannel);
133
134      void SetActualWidth(G4double* anActualWidth); 
135     
136      G4double EvaluateTotalActualWidth();
137
138      G4double EvaluateCMMomentum (const G4double mass,
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
145      G4double IntegrateCMMomentum2() const;
146     
147  public:
148     
149      G4double BrWig(const G4double Gamma, 
150                     const G4double rmass, 
151                     const G4double mass) const;
152
153private:
154      G4double IntegrandFunction1 (G4double xmass) const;
155      G4double IntegrandFunction2 (G4double xmass) const;
156      G4double IntegrandFunction3 (G4double xmass) const;
157      G4double IntegrandFunction4 (G4double xmass) const;
158public:
159  //   friend G4double IntegrandFunction3 (G4double xmass);
160
161  //   friend G4double IntegrandFunction4 (G4double xmass);
162     
163  // LB new variable created LB
164      G4int chosench;
165
166
167  private:
168 
169      G4ParticleDefinition* theDefinition;
170
171      G4double theFormationTime;
172
173      G4ThreeVector thePosition;
174
175      G4LorentzVector the4Momentum;
176      G4LorentzVector theFermi3Momentum;
177      G4LorentzVector theTotal4Momentum;
178     
179      G4Nucleon * theNucleon;
180     
181      G4int nChannels;
182     
183      G4double theActualMass;
184           
185      G4double* theActualWidth;
186
187     // Temporary storage for daughter masses and widths
188      // (needed because Integrand Function cannot take > 1 argument)
189      G4double* theDaughterMass;
190      G4double* theDaughterWidth;
191
192      CascadeState theStateToNucleus;
193
194      G4double theProjectilePotential;
195};
196
197// extern G4Allocator<G4KineticTrack> theKTAllocator;
198
199
200// Class G4KineticTrack
201/*
202inline void * G4KineticTrack::operator new(size_t)
203{
204  void * aT;
205  aT = (void *) theKTAllocator.MallocSingle();
206  return aT;
207}
208
209inline void G4KineticTrack::operator delete(void * aT)
210{
211  theKTAllocator.FreeSingle((G4KineticTrack *) aT);
212}
213*/
214
215inline G4ParticleDefinition* G4KineticTrack::GetDefinition() const
216{
217  return theDefinition;
218}
219
220inline void G4KineticTrack::SetDefinition(G4ParticleDefinition* aDefinition)
221{
222  theDefinition = aDefinition;
223}
224
225
226
227inline G4double G4KineticTrack::GetFormationTime() const
228{
229  return theFormationTime;
230}
231
232inline void G4KineticTrack::SetFormationTime(G4double aFormationTime)
233{
234  theFormationTime = aFormationTime;
235}
236
237
238
239inline const G4ThreeVector& G4KineticTrack::GetPosition() const
240{
241  return thePosition;
242}
243
244inline void G4KineticTrack::SetPosition(const G4ThreeVector aPosition)
245{
246  thePosition = aPosition;
247}
248
249
250inline const G4LorentzVector& G4KineticTrack::Get4Momentum() const
251{
252  return theTotal4Momentum;
253}
254
255inline const G4LorentzVector& G4KineticTrack::GetTrackingMomentum() const
256{
257   return the4Momentum;
258}
259
260inline void G4KineticTrack::Set4Momentum(const G4LorentzVector& a4Momentum)
261{
262//  set the4Momentum and update theTotal4Momentum
263
264  theTotal4Momentum=a4Momentum;
265  the4Momentum = theTotal4Momentum;
266  theFermi3Momentum=G4LorentzVector(0);
267}
268
269inline void G4KineticTrack::Update4Momentum(G4double aEnergy)
270{
271// update the4Momentum with aEnergy at constant mass (the4Momentum.mag() 
272//   updates theTotal4Momentum as well.
273  G4double newP(0);
274  G4double mass2=theTotal4Momentum.mag2();
275  if ( sqr(aEnergy) > mass2 )
276  {
277      newP = std::sqrt(sqr(aEnergy) - mass2 );
278  } else
279  {
280      aEnergy=std::sqrt(mass2);
281  }
282  Set4Momentum(G4LorentzVector(newP*the4Momentum.vect().unit(), aEnergy));
283}
284
285inline void G4KineticTrack::Update4Momentum(const G4ThreeVector & aMomentum)
286{
287// update the4Momentum with aMomentum at constant mass (the4Momentum.mag() 
288//   updates theTotal4Momentum as well.
289  G4double newE=std::sqrt(theTotal4Momentum.mag2() + aMomentum.mag2());
290  Set4Momentum(G4LorentzVector(aMomentum, newE));
291}
292
293inline void G4KineticTrack::SetTrackingMomentum(const G4LorentzVector& aMomentum)
294{
295//  set the4Momentum and update theTotal4Momentum, keep the mass of aMomentum
296
297  the4Momentum = aMomentum;
298  theTotal4Momentum=the4Momentum+theFermi3Momentum;
299//     keep mass of aMomentum for the total momentum
300  G4double m2 = aMomentum.mag2();
301  G4double p2=theTotal4Momentum.vect().mag2();
302  theTotal4Momentum.setE(std::sqrt(m2+p2));
303}
304
305inline void G4KineticTrack::UpdateTrackingMomentum(G4double aEnergy)
306{
307// update the4Momentum with aEnergy at constant mass (the4Momentum.mag() 
308//   updates theTotal4Momentum as well.
309  G4double newP(0);
310  G4double mass2=theTotal4Momentum.mag2();
311  if ( sqr(aEnergy) > mass2 )
312  {
313      newP = std::sqrt(sqr(aEnergy) - mass2 );
314  } else
315  {
316      aEnergy=std::sqrt(mass2);
317  }
318  SetTrackingMomentum(G4LorentzVector(newP*the4Momentum.vect().unit(), aEnergy));
319}
320
321inline void G4KineticTrack::UpdateTrackingMomentum(const G4ThreeVector & aMomentum)
322{
323// update the4Momentum with aMomentum at constant mass (the4Momentum.mag() 
324//   updates theTotal4Momentum as well.
325  G4double newE=std::sqrt(theTotal4Momentum.mag2() + aMomentum.mag2());
326  SetTrackingMomentum(G4LorentzVector(aMomentum, newE));
327}
328
329
330
331
332inline G4double G4KineticTrack::GetActualMass() const
333{
334  return std::sqrt(std::abs(the4Momentum.mag2()));
335}
336
337
338
339inline G4int G4KineticTrack::GetnChannels() const
340{
341  return nChannels;
342}
343
344inline void G4KineticTrack::SetnChannels(const G4int numberOfChannels)
345{
346  nChannels = numberOfChannels;
347}
348
349
350
351inline G4double* G4KineticTrack::GetActualWidth() const
352{
353  return theActualWidth;
354}
355
356inline void G4KineticTrack::SetActualWidth(G4double* anActualWidth)
357{
358  theActualWidth = anActualWidth;
359}
360
361
362
363inline G4double G4KineticTrack::EvaluateTotalActualWidth()
364{
365 G4int index;
366 G4double theTotalActualWidth = 0.0;
367 for (index = nChannels - 1; index >= 0; index--)
368    {
369     theTotalActualWidth += theActualWidth[index];
370    }
371 return theTotalActualWidth;
372}
373
374
375
376inline G4double G4KineticTrack::SampleResidualLifetime()
377{
378 G4double theTotalActualWidth = this->EvaluateTotalActualWidth();
379 G4double tau = hbar_Planck * (-1.0 / theTotalActualWidth);
380 G4double theResidualLifetime = tau * std::log(G4UniformRand());
381 return theResidualLifetime*the4Momentum.gamma();
382}
383
384
385
386inline G4double G4KineticTrack::EvaluateCMMomentum(const G4double m, 
387                                                 const G4double* m_ij) const
388{
389  G4double theCMMomentum;
390  if((m_ij[0]+m_ij[1])<m)
391   theCMMomentum = 1 / (2 * m) * 
392          std::sqrt (((m * m) - (m_ij[0] + m_ij[1]) * (m_ij[0] + m_ij[1])) *
393                ((m * m) - (m_ij[0] - m_ij[1]) * (m_ij[0] - m_ij[1])));
394  else
395   theCMMomentum=0.;
396
397 return theCMMomentum;
398}     
399
400inline G4double G4KineticTrack::BrWig(const G4double Gamma, const G4double rmass, const G4double mass) const 
401{               
402  G4double Norm = twopi;
403  return (Gamma/((mass-rmass)*(mass-rmass)+Gamma*Gamma/4.))/Norm;
404}
405     
406inline     
407void G4KineticTrack::Hit() 
408{
409  if(theNucleon) 
410  {
411    theNucleon->Hit(1);
412  }
413}
414
415inline
416G4bool G4KineticTrack::IsParticipant() const 
417{ 
418  if(!theNucleon) return true;
419  return theNucleon->AreYouHit(); 
420}
421
422inline 
423G4KineticTrack::CascadeState G4KineticTrack::GetState() const
424{
425        return theStateToNucleus;
426}
427
428inline
429G4KineticTrack::CascadeState G4KineticTrack::SetState(const CascadeState new_state)
430{
431        CascadeState old_state=theStateToNucleus;
432        theStateToNucleus=new_state;
433        return old_state;
434}
435
436inline
437void G4KineticTrack::SetProjectilePotential(G4double aPotential)
438{
439        theProjectilePotential = aPotential;
440}
441inline
442G4double G4KineticTrack::GetProjectilePotential() const
443{
444        return theProjectilePotential;
445}
446
447#endif
448
449
450
Note: See TracBrowser for help on using the repository browser.