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

Last change on this file since 908 was 819, checked in by garnier, 17 years ago

import all except CVS

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