source: trunk/source/processes/hadronic/models/photolepton_hadron/muon_nuclear/include/G4MuNuclearInteraction.hh @ 1340

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

update ti head

File size: 5.0 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: G4MuNuclearInteraction.hh,v 1.7 2009/03/04 19:09:20 vnivanch Exp $
28// GEANT4 tag $Name: geant4-09-03-ref-09 $
29//
30// $Id:
31// ------------------------------------------------------------
32//      GEANT 4 class header file
33//
34//      History: first implementation, based on object model of
35//      2nd December 1995, G.Cosmo
36//      -------- G4MuNuclearInteraction physics process ---------
37//                by Laszlo Urban, May 1998
38// ************************************************************
39
40#ifndef G4MuNuclearInteraction_h
41#define G4MuNuclearInteraction_h 1
42
43#include "G4ios.hh"
44#include "globals.hh"
45#include "Randomize.hh"
46#include "G4VDiscreteProcess.hh"
47#include "G4Track.hh"
48#include "G4Step.hh"
49#include "G4MuonMinus.hh"
50#include "G4MuonPlus.hh"
51#include "G4PionZero.hh"
52#include "G4OrderedTable.hh"
53#include "G4PhysicsTable.hh"
54#include "G4ElementTable.hh"
55#include "G4PhysicsLogVector.hh"
56#include "G4ParametrizedHadronicVertex.hh"
57#include "G4HadronicProcessType.hh"
58 
59class G4MuNuclearInteraction : public G4VDiscreteProcess
60 
61{ 
62public:
63 
64  G4MuNuclearInteraction(const G4String& processName = "muNuclear");
65 
66  virtual ~G4MuNuclearInteraction();
67
68  void SetPhysicsTableBining(G4double lowE, G4double highE, G4int nBins);
69
70  virtual G4bool IsApplicable(const G4ParticleDefinition&);
71
72  virtual void PreparePhysicsTable(const G4ParticleDefinition& ParticleType);
73
74  virtual void BuildPhysicsTable(const G4ParticleDefinition& ParticleType);
75
76  virtual void PrintInfoDefinition() ;
77
78  virtual G4VParticleChange *PostStepDoIt(const G4Track& track,
79                                          const G4Step& Step  ) ;
80
81protected:
82
83  virtual G4double GetMeanFreePath(const G4Track& track,
84                                   G4double previousStepSize,
85                                   G4ForceCondition* condition ) ;
86 
87  virtual G4double ComputeMeanFreePath( const G4ParticleDefinition* ParticleType,
88                                        G4double KineticEnergy,
89                                        const G4Material* aMaterial);
90
91  void ComputePartialSumSigma(  const G4ParticleDefinition* ParticleType,
92                                G4double KineticEnergy,
93                                const G4Material* aMaterial);
94
95  virtual G4double ComputeMicroscopicCrossSection(
96                                const G4ParticleDefinition* ParticleType,
97                                G4double KineticEnergy,
98                                G4double AtomicNumber,
99                                G4double AtomicMass);
100
101  virtual G4double ComputeDMicroscopicCrossSection(
102                                const G4ParticleDefinition* ParticleType,
103                                G4double KineticEnergy,
104                                G4double AtomicNumber,
105                                G4double AtomicMass,
106                                G4double epsilon);
107
108private:
109
110  G4MuNuclearInteraction & operator=(const G4MuNuclearInteraction &right);
111  G4MuNuclearInteraction(const G4MuNuclearInteraction&);
112
113  G4Element* SelectRandomAtom(G4Material* aMaterial) const;
114
115  void MakeSamplingTables( const G4ParticleDefinition* ParticleType );
116
117private:
118
119  G4PhysicsTable* theMeanFreePathTable;
120  G4PhysicsTable* theCrossSectionTable ;       
121
122  G4OrderedTable PartialSumSigma;     
123
124  G4double LowestKineticEnergy; 
125  G4double HighestKineticEnergy;   
126  G4int TotBin;     
127
128  //cut from R.P. Kokoulin
129  const G4double CutFixed ;
130  // for the atomic weight conversion
131  G4double GramPerMole ;
132
133  const G4MuonMinus* theMuonMinus;
134  const G4MuonPlus* theMuonPlus;
135  const G4PionZero* thePionZero;
136
137  // tables for sampling ..............
138  static G4int nzdat,ntdat,NBIN ;
139  static G4double zdat[5],adat[5],tdat[8] ;
140  static G4double ya[1001],proba[5][8][1001] ;
141     
142  G4ParametrizedHadronicVertex theHadronicVertex;
143};
144 
145#endif
Note: See TracBrowser for help on using the repository browser.