source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/include/G4QDiffraction.hh @ 836

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

import all except CVS

File size: 5.5 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// $Id: G4QDiffraction.hh,v 1.1 2007/08/23 15:58:43 mkossov Exp $
27// GEANT4 tag $Name:  $
28//
29//      ---------------- G4QDiffraction header ----------------
30//                 by Mikhail Kossov, Aug 2007.
31//  Header of G4QDiffraction class (hadron+A) of the CHIPS Simulation Branch in GEANT4
32// -------------------------------------------------------------------------------
33// This is a unique CHIPS class for the Hadron-Nuclear Diffractive Interaction Prosesses
34// -------------------------------------------------------------------------------
35// At present (Aug-07) it is based on the G4QDiffractionRatio class and is not tested.
36// The normalization is based on the temporary G4QProtonNuclearCrossSection class
37// -------------------------------------------------------------------------------
38
39#ifndef G4QDiffraction_hh
40#define G4QDiffraction_hh
41
42// GEANT4 Headers
43#include "globals.hh"
44#include "G4ios.hh"
45#include "Randomize.hh"
46#include "G4VDiscreteProcess.hh"
47#include "G4Track.hh"
48#include "G4Step.hh"
49#include "G4ParticleTypes.hh"
50#include "G4VParticleChange.hh"
51#include "G4ParticleDefinition.hh"
52#include "G4DynamicParticle.hh"
53#include "G4NucleiPropertiesTable.hh"
54#include "G4ThreeVector.hh"
55#include "G4LorentzVector.hh"
56
57// CHIPS Headers
58#include "G4QDiffractionRatio.hh"
59#include "G4QProtonNuclearCrossSection.hh"
60#include "G4QIsotope.hh"
61#include "G4QPDGToG4Particle.hh"
62#include "G4QCHIPSWorld.hh"
63#include "G4QHadronVector.hh"
64#include <vector>
65
66class G4QDiffraction : public G4VDiscreteProcess
67{
68public:
69
70  // Constructor
71  G4QDiffraction(const G4String& processName ="CHIPS_DiffractiveInteraction");
72
73  // Destructor
74  ~G4QDiffraction();
75
76  G4bool IsApplicable(const G4ParticleDefinition& particle);
77
78  G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize,
79                           G4ForceCondition* condition);
80  // It returns the MeanFreePath of the process for the current track :
81  // (energy, material)
82  // The previousStepSize and G4ForceCondition* are not used.
83  // This function overloads a virtual function of the base class.                   
84  // It is invoked by the ProcessManager of the Particle.
85 
86
87  G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep); 
88  // It computes the final state of the process (at end of step),
89  // returned as a ParticleChange object.                           
90  // This function overloads a virtual function of the base class.
91  // It is invoked by the ProcessManager of the Particle.
92
93
94  G4LorentzVector GetEnegryMomentumConservation();
95
96  G4int GetNumberOfNeutronsInTarget();
97
98private:
99
100  // Hide assignment operator as private
101  G4QDiffraction& operator=(const G4QDiffraction &right);
102
103  // Copy constructor
104  G4QDiffraction(const G4QDiffraction&);
105
106  // Calculate Cross-Section of the Diffraction Reaction (p is in GeV @@ units)
107  G4double CalculateXS(G4double p, G4int Z, G4int N, G4int pPDG);
108
109                // BODY
110  // Static Parameters --------------------------------------------------------------------
111  static G4int    nPartCWorld; // The#of particles for hadronization (limit of A of fragm.)
112  //--------------------------------- End of static parameters ---------------------------
113  // Working parameters
114  G4VQCrossSection* theCS;
115  G4LorentzVector EnMomConservation;                  // Residual of Energy/Momentum Cons.
116  G4int nOfNeutrons;                                  // #of neutrons in the target nucleus
117
118  // Modifires for the reaction
119  G4double Time;                                      // Time shift of the capture reaction
120  G4double EnergyDeposition;                          // Energy deposited in the reaction
121  static std::vector <G4int> ElementZ;                // Z of the element(i) in theLastCalc
122  static std::vector <G4double> ElProbInMat;          // SumProbabilityElements in Material
123  static std::vector <std::vector<G4int>*> ElIsoN;    // N of isotope(j) of Element(i)
124  static std::vector <std::vector<G4double>*> IsoProbInEl;// SumProbabIsotopes in Element i
125};
126#endif
127
Note: See TracBrowser for help on using the repository browser.