source: trunk/source/processes/hadronic/models/incl/include/G4Abla.hh @ 1201

Last change on this file since 1201 was 962, checked in by garnier, 15 years ago

update processes

File size: 10.6 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: G4Abla.hh,v 1.11 2008/06/25 17:20:03 kaitanie Exp $
27// Translation of INCL4.2/ABLA V3
28// Pekka Kaitaniemi, HIP (translation)
29// Christelle Schmidt, IPNL (fission code)
30// Alain Boudard, CEA (contact person INCL/ABLA)
31// Aatos Heikkinen, HIP (project coordination)
32
33#include "globals.hh"
34
35#include "G4InclRandomNumbers.hh"
36#include "G4AblaDataDefs.hh"
37#include "G4InclDataDefs.hh"
38#include "G4AblaFissionBase.hh"
39
40#ifndef G4Abla_hh
41#define G4Abla_hh 1
42
43/**
44 *  Class containing ABLA de-excitation code.
45 */
46
47class G4Abla {
48
49public:
50  /**
51   * Basic constructor.
52   */
53  G4Abla();
54
55  /**
56   * This constructor is used by standalone test driver and the Geant4 interface.
57   *
58   * @param aHazard random seeds
59   * @param aVolant data structure for ABLA output
60   * @param aVarNtp data structure for transfering ABLA output to Geant4 interface
61   */
62  G4Abla(G4Hazard *aHazard, G4Volant *aVolant, G4VarNtp *aVarntp);
63
64  /**
65   * Constructor that is to be used only for testing purposes.
66   * @param aHazard random seeds
67   * @param aVolant data structure for ABLA output   
68   */
69  G4Abla(G4Hazard *hazard, G4Volant *volant);
70
71  /**
72   * Basic destructor.
73   */
74  ~G4Abla();
75
76  /**
77   * Set verbosity level.
78   */
79  void setVerboseLevel(G4int level) {
80    verboseLevel = level;
81  }
82
83  /**
84   * Get the internal output data structure pointer.
85   */
86  G4Volant* getVolant() {
87    return volant;
88  }
89
90  /**
91   * Main interface to the de-excitation code.
92   *
93   * @param nucleusA mass number of the nucleus
94   * @param nucleusZ charge number of the nucleus
95   * @param nucleusMass mass of the nucleus
96   * @param excitationEnergy excitation energy of the nucleus
97   * @param angularMomentum angular momentum of the nucleus (produced as output by INCL4)
98   * @param recoilEnergy recoil energy of the nucleus
99   * @param momX momentum x-component
100   * @param momY momentum y-component
101   * @param momZ momentum z-component
102   * @param eventnumber number of the event
103   */
104  void breakItUp(G4double nucleusA, G4double nucleusZ, G4double nucleusMass, G4double excitationEnergy,
105                 G4double angularMomentum, G4double recoilEnergy, G4double momX, G4double momY, G4double momZ,
106                 G4int eventnumber);
107
108  // Evaporation
109public:
110  /**
111   * Initialize ABLA evaporation code.
112   *
113   */
114  void initEvapora();
115
116  /**
117   * Coefficient of collective enhancement including damping                         
118   * Input: z,a,bet,sig,u                                                 
119   * Output: qr - collective enhancement factor                           
120   * See  junghans et al., nucl. phys. a 629 (1998) 635                   
121   * @param z charge number
122   * @param a mass number
123   * @param bet beta deformation
124   * @param sig perpendicular spin cut-off factor
125   * @param u Energy
126   * @return Coefficient of collective enhancement   
127   */
128  void qrot(G4double z, G4double a, G4double bet, G4double sig, G4double u, G4double *qr);
129
130  /**
131   * Model de la goutte liquide de c. f. weizsacker.
132   * usually an obsolete option
133   */
134  void mglw(G4double a, G4double z, G4double *el);
135
136  /**
137   * Mglms
138   */
139  void mglms(G4double a, G4double z, G4int refopt4, G4double *el);
140
141  /**
142   *
143   */
144  G4double spdef(G4int a, G4int z, G4int optxfis);
145
146  /**
147   * Calculation of fissility parameter
148   */
149  G4double fissility(int a,int z, int optxfis);
150
151  /**
152   * Main evaporation routine.
153   */
154  void evapora(G4double zprf, G4double aprf, G4double ee, G4double jprf, 
155               G4double *zf_par, G4double *af_par, G4double *mtota_par,
156               G4double *pleva_par, G4double *pxeva_par, G4double *pyeva_par,
157               G4int *ff_par, G4int *inttype_par, G4int *inum_par);
158
159  /**
160   * Calculation of particle emission probabilities.
161   */
162  void direct(G4double zprf,G4double a, G4double ee, G4double jprf, 
163              G4double *probp_par, G4double *probn_par, G4double *proba_par, 
164              G4double *probf_par, G4double *ptotl_par, G4double *sn_par, G4double *sbp_par, G4double *sba_par, G4double *ecn_par, 
165              G4double *ecp_par,G4double *eca_par, G4double *bp_par, G4double *ba_par, G4int inttype, G4int inum, G4int itest);
166
167  /**
168   * Level density parameters.
169   */
170  void densniv(G4double a, G4double z, G4double ee, G4double esous, G4double *dens, G4double bshell, G4double bs, G4double bk, 
171               G4double *temp, G4int optshp, G4int optcol, G4double defbet);
172
173  /**
174   * This subroutine calculates the fission barriers                                                                 
175   * of the liquid-drop model of Myers and Swiatecki (1967).                                                                 
176   * Analytic parameterization of Dahlinger 1982
177   * replaces tables. Barrier heights from Myers and Swiatecki                                                               
178   */
179  G4double bfms67(G4double zms, G4double ams);
180
181  /**
182   * This subroutine calculates the ordinary legendre polynomials of   
183   * order 0 to n-1 of argument x and stores them in the vector pl.   
184   * They are calculated by recursion relation from the first two     
185   * polynomials.                                                     
186   * Written by A.J.Sierk  LANL  t-9  February, 1984                   
187   */
188  void lpoly(G4double x, G4int n, G4double pl[]);
189
190  /**
191   * This function will calculate the liquid-drop nuclear mass for spheri
192   * configuration according to the preprint NUCLEAR GROUND-STATE       
193   * MASSES and DEFORMATIONS by P. Mo"ller et al. from August 16, 1993 p.
194   * All constants are taken from this publication for consistency.     
195   */
196  G4double eflmac(G4int ia, G4int iz, G4int flag, G4int optshp);
197
198  /**
199   * Procedure for calculating the pairing correction to the binding   
200   * energy of a specific nucleus.
201   */
202  void appariem(G4double a, G4double z, G4double *del);
203
204  /**
205   * PROCEDURE FOR CALCULATING THE PARITY OF THE NUMBER N.             
206   * RETURNS -1 IF N IS ODD AND +1 IF N IS EVEN                       
207   */
208  void parite(G4double n, G4double *par);
209
210  /**
211   * RISE TIME IN WHICH THE FISSION WIDTH HAS REACHED     
212   * 90 PERCENT OF ITS FINAL VALUE
213   */
214  G4double tau(G4double bet, G4double homega, G4double ef, G4double t);
215
216  /**
217   * KRAMERS FAKTOR  - REDUCTION OF THE FISSION PROBABILITY       
218   * INDEPENDENT OF EXCITATION ENERGY
219   */
220  G4double cram(G4double bet, G4double homega);
221
222  /**
223   * CALCULATION OF THE SURFACE BS OR CURVATURE BK OF A NUCLEUS       
224   * RELATIVE TO THE SPHERICAL CONFIGURATION                           
225   * BASED ON  MYERS, DROPLET MODEL FOR ARBITRARY SHAPES               
226   */
227  G4double bipol(int iflag, G4double y);
228
229  /**
230   * THIS SUBROUTINE RETURNS THE BARRIER HEIGHT BFIS, THE             
231   * GROUND-STATE ENERGY SEGS, IN MEV, AND THE ANGULAR MOMENTUM       
232   * AT WHICH THE FISSION BARRIER DISAPPEARS, LMAX, IN UNITS OF       
233   * H-BAR, WHEN CALLED WITH INTEGER AGUMENTS IZ, THE ATOMIC           
234   * NUMBER, IA, THE ATOMIC MASS NUMBER, AND IL, THE ANGULAR           
235   * MOMENTUM IN UNITS OF H-BAR. (PLANCK'S CONSTANT DIVIDED BY         
236   * 2*PI).                                                           
237   */
238  void barfit(G4int iz, G4int ia, G4int il, G4double *sbfis, G4double *segs, G4double *selmax);
239
240  /**
241   * Random numbers.
242   */
243  G4double haz(G4int k);
244  void standardRandom(G4double *rndm, G4long *seed);
245
246  /**
247   * TIRAGE ALEATOIRE DANS UNE EXPONENTIELLLE : Y=EXP(-X/T)
248   */ 
249  G4double expohaz(G4int k, G4double T);
250
251  /**
252   * DISTRIBUTION DE MAXWELL
253   */
254  G4double fd(G4double E);
255
256  /**
257   *FONCTION INTEGRALE DE FD(E)
258   */
259  G4double f(G4double E);
260
261  /**
262   * tirage aleatoire dans une maxwellienne
263   */
264  G4double fmaxhaz(G4double T);
265
266  /**
267   *
268   */
269  G4double pace2(G4double a, G4double z);
270
271  /**
272   *
273   */
274  void guet(G4double *x_par, G4double *z_par, G4double *find_par);
275   
276public:
277  // Coordinate system transformations:
278  void lorab(G4double gam, G4double eta, G4double ein, G4double pin[],
279             G4double *eout, G4double pout[]);
280
281  void translab(G4double gamrem, G4double etrem, G4double csrem[4], G4int nopart, G4int ndec);
282  void translabpf(G4double masse1, G4double t1, G4double p1, G4double ctet1,
283                  G4double phi1, G4double gamrem, G4double etrem, G4double R[][4],
284                  G4double *plab1, G4double *gam1, G4double *eta1, G4double csdir[]);
285
286  void rotab(G4double R[4][4], G4double pin[4], G4double pout[4]);
287
288  // Utils
289  G4int min(G4int a, G4int b);
290  G4double min(G4double a, G4double b);
291  G4int max(G4int a, G4int b);
292  G4double max(G4double a, G4double b);
293
294  G4int nint(G4double number);
295  G4int secnds(G4int x);
296  G4int mod(G4int a, G4int b);
297  G4double dmod(G4double a, G4double b);
298  G4double dint(G4double a);
299  G4int idint(G4double a);
300  G4int idnint(G4double value);
301  G4double utilabs(G4double a);
302  G4double dmin1(G4double a, G4double b, G4double c);
303  G4Ec2sub* getFrldmTable() {
304    return ec2sub;
305  }
306
307private:
308  G4int verboseLevel;
309  G4int ilast;
310
311  G4AblaFissionBase *fissionModel;
312  G4InclRandomInterface *randomGenerator;
313  G4Pace *pace;
314  G4Hazard *hazard;
315  G4Ald *ald;
316  G4Eenuc *eenuc;
317  G4Ec2sub *ec2sub;
318  G4Ecld *ecld; 
319  G4Fb *fb;
320  G4Fiss *fiss;
321  G4Opt *opt;
322  G4Volant *volant;
323  G4VarNtp *varntp; 
324};
325
326#endif
Note: See TracBrowser for help on using the repository browser.