source: trunk/source/processes/hadronic/models/abla/include/G4Abla.hh @ 1355

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

update to last version 4.9.4

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