source: trunk/source/processes/decay/include/G4Decay.hh @ 819

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

import all except CVS

File size: 7.3 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// $Id: G4Decay.hh,v 1.18 2007/07/23 23:13:04 kurasige Exp $
28// GEANT4 tag $Name: geant4-09-01-patch-02 $
29//
30//
31// ------------------------------------------------------------
32//      GEANT 4 class header file
33//
34//      History: first implementation, based on object model of
35//      7 July 1996 H.Kurashige
36// ------------------------------------------------------------
37//  New Physics scheme           18 Jan. 1997  H.Kurahige
38// ------------------------------------------------------------
39//   modified                     4  Feb. 1997  H.Kurahige
40//   modified                     8  Sep. 1997  H.Kurahige
41//   remove BuildPhysicsTable()   27 Nov. 1997   H.Kurashige
42//   modified for new ParticleChange 12 Mar. 1998  H.Kurashige
43//   added aPhysicsTable          2  Aug. 1998 H.Kurashige
44//   PreAssignedDecayTime         18 Jan. 2001 H.Kurashige
45//   Add External Decayer         23 Feb. 2001  H.Kurashige
46//   Remove PhysicsTable          12 Feb. 2002 H.Kurashige
47//   Fixed bug in PostStepGPIL
48//    in case of stopping during AlongStepDoIt 12 Mar. 2004 H.Kurashige
49//   Add GetRemainderLifeTime  10 Aug/2004 H.Kurashige
50//   Add DaughterPolarization     23 July 2008 H.Kurashige
51
52
53#ifndef G4Decay_h
54#define G4Decay_h 1
55
56#include "G4ios.hh"
57#include "globals.hh"
58#include "G4VRestDiscreteProcess.hh"
59#include "G4ParticleChangeForDecay.hh"
60class G4VExtDecayer;
61
62class G4Decay : public G4VRestDiscreteProcess
63{
64 // Class Description
65  //  This class is a decay process
66
67  public:
68    //  Constructors
69    G4Decay(const G4String& processName ="Decay");
70
71    //  Destructor
72    virtual ~G4Decay();
73
74  private:
75    //  copy constructor
76      G4Decay(const G4Decay &right);
77
78    //  Assignment Operation (generated)
79      G4Decay & operator=(const G4Decay &right);
80
81  public: //With Description
82     // G4Decay Process has both
83     // PostStepDoIt (for decay in flight)
84     //   and
85     // AtRestDoIt (for decay at rest)
86 
87     virtual G4VParticleChange *PostStepDoIt(
88                             const G4Track& aTrack,
89                             const G4Step& aStep
90                            );
91
92     virtual G4VParticleChange* AtRestDoIt(
93                             const G4Track& aTrack,
94                             const G4Step&  aStep
95                            );
96
97     virtual void BuildPhysicsTable(const G4ParticleDefinition&); 
98     // In G4Decay, thePhysicsTable stores values of
99    //    beta * std::sqrt( 1 - beta*beta)
100    //  as a function of normalized kinetic enregy (=Ekin/mass),
101    //  becasuse this table is universal for all particle types,
102
103
104    virtual G4bool IsApplicable(const G4ParticleDefinition&);
105    // returns "true" if the decay process can be applied to
106    // the particle type.
107 
108  protected: // With Description
109    virtual G4VParticleChange* DecayIt(
110                             const G4Track& aTrack,
111                             const G4Step&  aStep
112                            );
113    // The DecayIt() method returns by pointer a particle-change object,
114    // which has information of daughter particles.
115
116    // Set daughter polarization
117    //  NO OPERATION in the base class of G4Decay
118    virtual void DaughterPolarization(const G4Track& aTrack,
119                              G4DecayProducts* products);
120
121 public:
122    virtual G4double AtRestGetPhysicalInteractionLength(
123                             const G4Track& track,
124                             G4ForceCondition* condition
125                            );
126
127    virtual G4double PostStepGetPhysicalInteractionLength(
128                             const G4Track& track,
129                             G4double   previousStepSize,
130                             G4ForceCondition* condition
131                            );
132
133  protected: // With Description
134    // GetMeanFreePath returns ctau*beta*gamma for decay in flight
135    // GetMeanLifeTime returns ctau for decay at rest
136    virtual G4double GetMeanFreePath(const G4Track& aTrack,
137                              G4double   previousStepSize,
138                              G4ForceCondition* condition
139                             );
140
141    virtual G4double GetMeanLifeTime(const G4Track& aTrack,
142                              G4ForceCondition* condition
143                            );
144
145   public: //With Description
146     virtual void StartTracking(G4Track*);
147     virtual void EndTracking();
148      // inform Start/End of tracking for each track to the physics process
149
150   public: //With Description
151     void SetExtDecayer(G4VExtDecayer*);
152     const G4VExtDecayer* GetExtDecayer() const;
153     // Set/Get External Decayer
154   
155    G4double GetRemainderLifeTime() const; 
156    //Get Remainder of life time at rest decay
157
158  public:
159     void  SetVerboseLevel(G4int value);
160     G4int GetVerboseLevel() const;
161
162  protected:
163     G4int verboseLevel;
164     // controle flag for output message
165     //  0: Silent
166     //  1: Warning message
167     //  2: More
168
169  protected:
170    // HighestValue.
171    const G4double HighestValue;
172 
173    // Remainder of life time at rest
174    G4double                 fRemainderLifeTime;
175 
176    // ParticleChange for decay process
177    G4ParticleChangeForDecay fParticleChangeForDecay;
178   
179    // External Decayer
180    G4VExtDecayer*    pExtDecayer;
181};
182
183inline
184 void  G4Decay::SetVerboseLevel(G4int value){ verboseLevel = value; }
185
186inline
187 G4int G4Decay::GetVerboseLevel() const { return verboseLevel; }
188
189inline 
190  G4VParticleChange* G4Decay::AtRestDoIt(
191                             const G4Track& aTrack,
192                             const G4Step&  aStep
193                            )
194{
195  return DecayIt(aTrack, aStep);
196}
197
198inline 
199  G4VParticleChange* G4Decay::PostStepDoIt(
200                             const G4Track& aTrack,
201                             const G4Step&  aStep
202                            )
203{
204  return DecayIt(aTrack, aStep);
205}
206
207inline
208 void G4Decay::SetExtDecayer(G4VExtDecayer* val)
209{
210  pExtDecayer = val;
211}
212
213inline
214 const G4VExtDecayer* G4Decay::GetExtDecayer() const
215{
216  return pExtDecayer;
217}
218
219inline
220 G4double G4Decay::GetRemainderLifeTime() const 
221{
222  return fRemainderLifeTime;
223}
224
225#endif
226
227
228
229
230
231
232
233
234
235
Note: See TracBrowser for help on using the repository browser.