Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4WilsonAbrasionModel.hh
Go to the documentation of this file.
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 // * *
21 // * Parts of this code which have been developed by QinetiQ Ltd *
22 // * under contract to the European Space Agency (ESA) are the *
23 // * intellectual property of ESA. Rights to use, copy, modify and *
24 // * redistribute this software for general public use are granted *
25 // * in compliance with any licensing, distribution and development *
26 // * policy adopted by the Geant4 Collaboration. This code has been *
27 // * written by QinetiQ Ltd for the European Space Agency, under ESA *
28 // * contract 17191/03/NL/LvH (Aurora Programme). *
29 // * *
30 // * By using, copying, modifying or distributing the software (or *
31 // * any work based on the software) you agree to acknowledge its *
32 // * use in resulting scientific publications, and indicate your *
33 // * acceptance of all terms of the Geant4 Software license. *
34 // ********************************************************************
35 //
36 #ifndef G4WilsonAbrasionModel_h
37 #define G4WilsonAbrasionModel_h
38 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
39 //
40 // MODULE: G4WilsonAbrasionModel.hh
41 //
42 // Version: 1.0
43 // Date: 08/12/2009
44 // Author: P R Truscott
45 // Organisation: QinetiQ Ltd, UK
46 // Customer: ESA/ESTEC, NOORDWIJK
47 // Contract: 17191/03/NL/LvH
48 //
49 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
50 //
51 // CHANGE HISTORY
52 // --------------
53 //
54 // 6 October 2003, P R Truscott, QinetiQ Ltd, UK
55 // Created.
56 //
57 // 15 March 2004, P R Truscott, QinetiQ Ltd, UK
58 // Beta release
59 //
60 // 08 December 2009, P R Truscott, QinetiQ Ltd, Ltd
61 // ver 1.0
62 // Variable fradius defined. See .cc file for more details.
63 //
64 
65 #include "globals.hh"
66 
67 #include "G4HadronicInteraction.hh"
68 #include "G4WilsonAblationModel.hh"
69 #include "G4ExcitationHandler.hh"
70 #include "G4HadFinalState.hh"
71 #include "G4Track.hh"
72 #include "G4Nucleus.hh"
73 #include "G4Fragment.hh"
74 #include "G4HadProjectile.hh"
75 
76 
78 {
79  public:
80  G4WilsonAbrasionModel(G4bool useAblation1 = false);
83 
85 
87 
89  void SetVerboseLevel(G4int);
90  void SetUseAblation(G4bool);
96 
97  virtual void ModelDescription(std::ostream&) const;
98 
99  private:
100  void PrintWelcomeMessage();
101  G4Fragment* GetAbradedNucleons(G4int, G4double, G4double, G4double);
102  G4double GetNucleonInducedExcitation(G4double, G4double, G4double);
103  void SetConserveEnergy(G4bool);
104  G4bool GetConserveEnergy();
105 
106  private:
107  G4double r0sq;
108  G4double npK;
109  G4bool useAblation;
110  G4WilsonAblationModel* theAblation;
111  G4ExcitationHandler* theExcitationHandler;
112  G4ExcitationHandler* theExcitationHandlerx;
113  G4bool conserveEnergy;
114  G4bool conserveMomentum;
115  G4double B;
116  G4double third;
117  G4double fradius;
118 };
119 
120 inline void
122  {theExcitationHandler = aExcitationHandler;}
123 
125  {return theExcitationHandler;}
126 
128  {return useAblation;}
129 
130 inline void G4WilsonAbrasionModel::SetConserveEnergy(G4bool conserveEnergy1)
131  {conserveEnergy = conserveEnergy1;}
132 
133 inline G4bool G4WilsonAbrasionModel::GetConserveEnergy()
134  {return conserveEnergy;}
135 
137  {conserveMomentum = conserveMomentum1;}
138 
140  {return conserveMomentum;}
141 
143 {
144  verboseLevel = verboseLevel1;
145  if (useAblation) theAblation->SetVerboseLevel(verboseLevel);
146 }
147 
148 #endif
void SetExcitationHandler(G4ExcitationHandler *)
G4WilsonAbrasionModel(G4bool useAblation1=false)
int G4int
Definition: G4Types.hh:78
virtual void ModelDescription(std::ostream &) const
bool G4bool
Definition: G4Types.hh:79
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &, G4Nucleus &)
const G4WilsonAbrasionModel & operator=(G4WilsonAbrasionModel &right)
G4ExcitationHandler * GetExcitationHandler()
double G4double
Definition: G4Types.hh:76