source: trunk/source/processes/management/include/G4VDiscreteProcess.hh @ 1346

Last change on this file since 1346 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 6.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//
28// $Id: G4VDiscreteProcess.hh,v 1.9 2007/11/15 04:10:18 kurasige Exp $
29// GEANT4 tag $Name: geant4-09-04-beta-01 $
30//
31//
32// ------------------------------------------------------------
33//      GEANT 4 class header file
34//
35//      History: first implementation, based on object model of
36//      2nd December 1995, G.Cosmo
37//      add G4VDiscreteProcess(const G4String&) 24 Jul 1996, Hisaya kurashige
38//
39// Class Description 
40//  Abstract class which defines the public behavior of
41//  discrete physics interactions.
42//
43// ------------------------------------------------------------
44//   New Physics scheme           18 Dec. 1996  H.Kurahige
45// ------------------------------------------------------------
46//   modified for new ParticleChange 12 Mar. 1998  H.Kurashige
47//   Fixed a bug in PostStepGetPhysicalInteractionLength 
48//                                15 Apr. 2002 H.Kurashige
49//
50
51#ifndef G4VDiscreteProcess_h
52#define G4VDiscreteProcess_h 1
53
54#include "globals.hh"
55#include "G4ios.hh"
56
57#include "G4VProcess.hh"
58
59class G4VDiscreteProcess : public G4VProcess
60{
61  //  Abstract class which defines the public behavior of
62  //  discrete physics interactions.
63  public:     
64
65      G4VDiscreteProcess(const G4String& ,
66                         G4ProcessType   aType = fNotDefined );
67      G4VDiscreteProcess(G4VDiscreteProcess &);
68
69      virtual ~G4VDiscreteProcess();
70
71  public :// with description
72      virtual G4double PostStepGetPhysicalInteractionLength(
73                             const G4Track& track,
74                             G4double   previousStepSize,
75                             G4ForceCondition* condition
76                            );
77
78      virtual G4VParticleChange* PostStepDoIt(
79                             const G4Track& ,
80                             const G4Step& 
81                            );
82
83     //  no operation in  AtRestDoIt and  AlongStepDoIt
84      virtual G4double AlongStepGetPhysicalInteractionLength(
85                             const G4Track&,
86                             G4double  ,
87                             G4double  ,
88                             G4double& ,
89                             G4GPILSelection*
90                            ){ return -1.0; };
91
92      virtual G4double AtRestGetPhysicalInteractionLength(
93                             const G4Track& ,
94                             G4ForceCondition* 
95                            ) { return -1.0; };
96
97     //  no operation in  AtRestDoIt and  AlongStepDoIt
98      virtual G4VParticleChange* AtRestDoIt(
99                             const G4Track& ,
100                             const G4Step&
101                            ) {return 0;};
102
103      virtual G4VParticleChange* AlongStepDoIt(
104                             const G4Track& ,
105                             const G4Step& 
106                            ) {return 0;};
107 
108  protected:// with description
109     virtual G4double GetMeanFreePath(const G4Track& aTrack,
110                             G4double   previousStepSize,
111                             G4ForceCondition* condition
112                                                               )=0;
113      //  Calculates from the macroscopic cross section a mean
114      //  free path, the value is returned in units of distance.
115 
116  private:
117  // hide default constructor and assignment operator as private
118      G4VDiscreteProcess();
119      G4VDiscreteProcess & operator=(const G4VDiscreteProcess &right);
120
121};
122
123// -----------------------------------------
124//  inlined function members implementation
125// -----------------------------------------
126#include "G4Step.hh"
127#include "G4Track.hh"
128#include "G4MaterialTable.hh"
129#include "G4VParticleChange.hh"
130
131inline G4double G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(
132                             const G4Track& track,
133                             G4double   previousStepSize,
134                             G4ForceCondition* condition
135                            )
136{
137  if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0)) {
138    // beggining of tracking (or just after DoIt of this process)
139    ResetNumberOfInteractionLengthLeft();
140  } else if ( previousStepSize > 0.0) {
141    // subtract NumberOfInteractionLengthLeft
142    SubtractNumberOfInteractionLengthLeft(previousStepSize);
143  } else {
144    // zero step
145    //  DO NOTHING
146  }
147
148  // condition is set to "Not Forced"
149  *condition = NotForced;
150
151  // get mean free path
152  currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
153
154  G4double value;
155  if (currentInteractionLength <DBL_MAX) {
156    value = theNumberOfInteractionLengthLeft * currentInteractionLength;
157  } else {
158    value = DBL_MAX;
159  }
160#ifdef G4VERBOSE
161  if (verboseLevel>1){
162    G4cout << "G4VDiscreteProcess::PostStepGetPhysicalInteractionLength ";
163    G4cout << "[ " << GetProcessName() << "]" <<G4endl;
164    track.GetDynamicParticle()->DumpInfo();
165    G4cout << " in Material  " <<  track.GetMaterial()->GetName() <<G4endl;
166    G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
167  }
168#endif
169  return value;
170}
171
172inline G4VParticleChange* G4VDiscreteProcess::PostStepDoIt(
173                             const G4Track& ,
174                             const G4Step& 
175                            )
176{ 
177//  clear NumberOfInteractionLengthLeft
178    ClearNumberOfInteractionLengthLeft();
179
180    return pParticleChange;
181}
182
183
184#endif
185
186
187
Note: See TracBrowser for help on using the repository browser.