source: trunk/source/processes/hadronic/models/incl/include/G4InclDataDefs.hh@ 1036

Last change on this file since 1036 was 962, checked in by garnier, 17 years ago

update processes

File size: 15.4 KB
RevLine 
[819]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//
[962]26// $Id: G4InclDataDefs.hh,v 1.5 2008/06/25 17:20:04 kaitanie Exp $
[819]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// All data structures needed by INCL4 are defined here.
34
35#ifndef InclDataDefs_hh
36#define InclDataDefs_hh 1
37
38#define FSIZE 15
39/**
40 * Initial values of a hadronic cascade problem.
41 */
42class G4Calincl {
43public:
44 G4Calincl() {};
45 ~G4Calincl() {};
46
47 /**
48 * Here f is an array containing the following initial values:
49 * - f[0] : target mass number
50 * - f[1] : target charge number
51 * - f[2] : bullet energy
52 * - f[3] : minimum proton energy to leave the target (default: 0.0)
53 * - f[4] : nuclear potential (default: 45.0 MeV)
54 * - f[5] : time scale (default: 1.0)
55 * - f[6] : bullet type (1: proton, 2: neutron, 3: pi+, 4: pi0 5: pi-, 6:H2, 7: H3, 8: He3, 9: He4
56 * - f[7] : minimum neutron energy to leave the target (default: 0.0)
57 * - f[8] : target material identifier (G4Mat)
58 * - f[9] : not used
59 * - f[10] : not used
60 * - f[11] : not used
61 * - f[12] : not used
62 * - f[13] : not used
63 * - f[14] : not used
64 */
65 G4double f[FSIZE];
66
67 /**
68 * Number of events to be processed.
69 */
70 G4int icoup;
71};
72
73#define IGRAINESIZE 19
74/**
75 * Random seeds used by internal random number generators.
76 * @see G4Incl::standardRandom
77 * @see G4Incl::gaussianRandom
78 */
79class G4Hazard{
80public:
81 G4Hazard() {};
82 ~G4Hazard() {};
83
84 /**
85 * Random seed
86 */
87 G4long ial;
88
89 /**
90 * An array of random seeds.
91 */
92 G4long igraine[IGRAINESIZE];
93};
94
95#define MATSIZE 500
96#define MATGEOSIZE 6
97/**
98 * Target nuclei to be taken into account in the cascade problem.
99 */
100class G4Mat {
101public:
102 G4Mat() { };
103 ~G4Mat() { };
104
105 /**
106 * Charge numbers.
107 */
108 G4int zmat[MATSIZE];
109
110 /**
111 * Mass number
112 */
113 G4int amat[MATSIZE];
114
115 /**
116 *
117 */
118 G4double bmax_geo[MATGEOSIZE][MATSIZE];
119
120 /**
121 * Number of materials.
122 */
123 G4int nbmat;
124};
125
126#define LGNSIZE 9
127/**
128 * Properties of light nucleus used as a bullet.
129 */
130class G4LightGausNuc {
131public:
132 G4LightGausNuc() {};
133 ~G4LightGausNuc() {};
134
135 G4double rms1t[LGNSIZE];
136 G4double pf1t[LGNSIZE];
137 G4double pfln[LGNSIZE];
138 G4double tfln[LGNSIZE];
139 G4double vnuc[LGNSIZE];
140};
141
142#define LNSIZE 30
143/**
144 * Data of light nuclei.
145 */
146class G4LightNuc {
147public:
148 G4LightNuc() {};
149 ~G4LightNuc() {};
150
151 /**
152 * r
153 */
154 G4double r[LNSIZE];
155
156 /**
157 * a
158 */
159 G4double a[LNSIZE];
160};
161
162#define SAXWROWS 30
163#define SAXWCOLS 500
164/**
165 * Woods-Saxon density and its first derivative.
166 */
167class G4Saxw {
168public:
169 G4Saxw() {};
170 ~G4Saxw() {};
171
172 /**
173 * x
174 */
175 G4double x[SAXWROWS][SAXWCOLS];
176
177 /**
178 * y
179 */
180 G4double y[SAXWROWS][SAXWCOLS];
181
182 /**
183 * s
184 */
185 G4double s[SAXWROWS][SAXWCOLS];
186
187 /**
188 * imat
189 */
190 G4int imat;
191
192 /**
193 * n
194 */
195 G4int n;
196
197 /**
198 * k
199 */
200 G4int k;
201};
202
203/**
204 * Parameters for INCL4 model.
205 */
206class G4Ws {
207public:
208 G4Ws() {};
209 ~G4Ws() {};
210
211 /**
212 * r0
213 */
214 G4double r0;
215
216 /**
217 * adif
218 */
219 G4double adif;
220
221 /**
222 * Maximum radius of the nucleus
223 */
224 G4double rmaxws;
225
226 /**
227 * drws
228 */
229 G4double drws;
230
231 /**
232 * Shape of the surface of the nucleus:
233 * - -1: Woods-Saxon density with impact parameter dependence
234 * - 0: Woods-Saxon density without impact parameter dependence
235 * - 1: Sharp surface (hard sphere)
236 */
237 G4double nosurf;
238
239 /**
240 * Parameter related to the maximum radius of the nucleus.
241 *
242 * rmaxws = r0 + xfoisa*A
243 */
244 G4double xfoisa;
245
246 /**
247 * Pauli blocking used in the simulation:
248 * - 0: statistic Pauli blocking
249 * - 1: strict Pauli blocking
250 * - 2: no Pauli blocking
251 */
252 G4double npaulstr;
253
254 /**
255 * Maximum impact parameter
256 */
257 G4double bmax;
258};
259
260#define DTONSIZE 13
261/**
262 * Random seeds used by internal random number generators.
263 * @see G4Incl::standardRandom
264 * @see G4Incl::gaussianRandom
265 */
266class G4Dton {
267public:
268 G4Dton() {};
269 ~G4Dton() {};
270
271 G4double c[DTONSIZE];
272 G4double d[DTONSIZE];
273 G4double fn;
274};
275
276#define SPL2SIZE 100
277/**
278 * Random seeds used by internal random number generators.
279 * @see G4Incl::standardRandom
280 * @see G4Incl::gaussianRandom
281 */
282class G4Spl2 {
283public:
284 G4Spl2() {};
285 ~G4Spl2() {};
286
287 G4double x[SPL2SIZE];
288 G4double y[SPL2SIZE];
289 G4double a[SPL2SIZE];
290 G4double b[SPL2SIZE];
291 G4double c[SPL2SIZE];
292 G4int n;
293};
294
295// incl4.2.cc:
296
297//#define BL1SIZE 300
298#define BL1SIZE 3000
299/**
300 * Random seeds used by internal random number generators.
301 * @see G4Incl::standardRandom
302 * @see G4Incl::gaussianRandom
303 */
304class G4Bl1 {
305public:
306 G4Bl1() {};
307 ~G4Bl1() {};
308
309 G4double p1[BL1SIZE],p2[BL1SIZE],p3[BL1SIZE];
310 G4double eps[BL1SIZE];
311 G4int ind1[BL1SIZE],ind2[BL1SIZE];
312 G4double ta;
313};
314
315#define BL2CROISSIZE 19900
316#define BL2INDSIZE 19900
317/**
318 *
319 */
320class G4Bl2 {
321public:
322 G4Bl2() {};
323 ~G4Bl2() {};
324
325 /**
326 *
327 */
328 G4double crois[BL2CROISSIZE];
329
330 /**
331 *
332 */
333 G4int k;
334
335 /**
336 *
337 */
338 G4int ind[BL2INDSIZE];
339
340 /**
341 *
342 */
343 G4int jnd[BL2INDSIZE];
344};
345
346//#define BL3SIZE 300
347#define BL3SIZE 3000
348/**
349 *
350 */
351class G4Bl3 {
352public:
353 G4Bl3() {};
354 ~G4Bl3() {};
355
356 /**
357 * r1 and r2
358 */
359 G4double r1,r2;
360
361 /**
362 * Nucleon positions
363 */
364 G4double x1[BL3SIZE], x2[BL3SIZE],x3[BL3SIZE];
365
366 /**
367 * Mass numbers
368 */
369 G4int ia1,ia2;
370
371 /**
372 * rab2
373 */
374 G4double rab2;
375};
376
377/**
378 * G4Bl4
379 */
380class G4Bl4 {
381public:
382 G4Bl4() {};
383 ~G4Bl4() {};
384
385 /**
386 * tmax5
387 */
388 G4double tmax5;
389};
390
391//#define BL5SIZE 300
392#define BL5SIZE 3000
393/**
394 * G4Bl5
395 */
396class G4Bl5 {
397public:
398 G4Bl5() {};
399 ~G4Bl5() {};
400
401 /**
402 * tlg
403 */
404 G4double tlg[BL5SIZE];
405
406 /**
407 * nesc
408 */
409 G4int nesc[BL5SIZE];
410};
411
412/**
413 * G4Bl6
414 */
415class G4Bl6 {
416public:
417 G4Bl6() {};
418 ~G4Bl6() {};
419
420 /**
421 * xx10
422 */
423 G4double xx10;
424
425 /**
426 * isa
427 */
428 G4double isa;
429};
430
431/**
432 * G4Bl8
433 */
434class G4Bl8 {
435public:
436 G4Bl8() {};
437 ~G4Bl8() {};
438
439 /**
440 * rathr
441 */
442 G4double rathr;
443
444 /**
445 * ramass
446 */
447 G4double ramass;
448};
449
450//#define BL9SIZE 300
451#define BL9SIZE 3000
452/**
453 * G4Bl9
454 */
455class G4Bl9 {
456public:
457 G4Bl9() {
458 l1 = 0;
459 l2 = 0;
460 };
461 ~G4Bl9() {};
462
463 /**
464 * hel
465 */
466 G4double hel[BL9SIZE];
467
468 /**
469 * l1 and l2
470 */
471 G4int l1,l2;
472};
473
474/**
475 * G4Bl10
476 */
477class G4Bl10 {
478public:
479 G4Bl10() {};
480 ~G4Bl10() {};
481
482 /**
483 * ri4, rs4, r2i, r2s, pdummy, pf
484 */
485 G4double ri4,rs4,r2i,r2s,pdummy,pf;
486};
487
488/**
489 * G4Kind
490 */
491class G4Kind {
492public:
493 G4Kind() {};
494 ~G4Kind() {};
495
496 /**
497 * kindf7
498 */
499 G4int kindf7;
500};
501
502#define VARSIZE 3
503#define VAEPSSIZE 250
504#define VAAVM 1000
505/**
506 * Extra information on collisions between nucleons.
507 */
508class G4VarAvat {
509public:
510 G4VarAvat() {};
511 ~G4VarAvat() {};
512
513 /**
514 *
515 */
516 G4int kveux;
517
518 /**
519 *
520 */
521 G4double bavat;
522
523 /**
524 *
525 */
526 G4int nopartavat,ncolavat;
527
528 /**
529 *
530 */
531 G4double r1_in[VARSIZE],r1_first_avat[VARSIZE];
532
533 /**
534 *
535 */
536 G4double epsd[VAEPSSIZE],eps2[VAEPSSIZE],eps4[VAEPSSIZE],eps6[VAEPSSIZE],epsf[VAEPSSIZE];
537
538 /**
539 *
540 */
541 G4int nb_avat;
542
543 /**
544 *
545 */
546 G4double timeavat[VAAVM],l1avat[VAAVM],l2avat[VAAVM],jpartl1[VAAVM],jpartl2[VAAVM];
547
548 /**
549 *
550 */
551 G4double del1avat[VAAVM],del2avat[VAAVM],energyavat[VAAVM];
552
553 /**
554 *
555 */
556 G4double bloc_paul[VAAVM],bloc_cdpp[VAAVM],go_out[VAAVM];
557};
558
559#define VARNTPSIZE 255
560class G4VarNtp {
561public:
562 G4VarNtp() {};
563 ~G4VarNtp() {};
564
565 /**
[962]566 * Clear and initialize all variables and arrays.
567 */
568 void clear() {
569 particleIndex = 0;
570 projType = 0;
571 projEnergy = 0.0;
572 targetA = 0;
573 targetZ = 0;
574 massini = 0;
575 mzini = 0;
576 exini = 0;
577 pcorem = 0;
578 mcorem = 0;
579 pxrem = 0;
580 pyrem = 0;
581 pzrem = 0;
582 mulncasc = 0;
583 mulnevap = 0;
584 mulntot = 0;
585 bimpact = 0.0;
586 jremn = 0;
587 kfis = 0;
588 estfis = 0;
589 izfis = 0;
590 iafis = 0;
591 ntrack = 0;
592 for(G4int i = 0; i < VARNTPSIZE; i++) {
593 itypcasc[i] = 0;
594 avv[i] = 0;
595 zvv[i] = 0;
596 enerj[i] = 0.0;
597 plab[i] = 0.0;
598 tetlab[i] = 0.0;
599 philab[i] = 0.0;
600 full[i] = false;
601 }
602 }
603
604 void addParticle(G4double A, G4double Z, G4double E, G4double P, G4double theta, G4double phi) {
605 if(full[particleIndex]) {
606 G4cout <<"G4VarNtp: Error. Index i = " << particleIndex << " is already occupied by particle:" << G4endl;
607 G4cout <<"A = " << avv[particleIndex] << " Z = " << zvv[particleIndex] << G4endl;
608 G4cout <<"Tried to replace it with:" << G4endl;
609 G4cout <<"A = " << Z << " Z = " << Z << G4endl;
610 } else {
611 avv[particleIndex] = (int) A;
612 zvv[particleIndex] = (int) Z;
613 enerj[particleIndex] = E;
614 plab[particleIndex] = P;
615 tetlab[particleIndex] = theta;
616 philab[particleIndex] = phi;
617 full[particleIndex] = true;
618 ntrack = particleIndex + 1;
619 particleIndex++;
620 }
621 }
622
623 /**
624 * Baryon number conservation check.
625 */
626 G4int getTotalBaryonNumber() {
627 G4int baryonNumber = 0;
628 for(G4int i = 0; i < ntrack; i++) {
629 if(avv[i] > 0) {
630 baryonNumber += avv[i];
631 }
632 }
633 return baryonNumber;
634 }
635
636 /**
637 * Return total energy.
638 */
639 G4double getTotalEnergy() {
640 G4double energy = 0.0;
641 for(G4int i = 0; i < ntrack; i++) {
642 energy += std::sqrt(std::pow(plab[i], 2) + std::pow(getMass(i), 2)); // E^2 = p^2 + m^2
643 }
644
645 return energy;
646 }
647
648 /**
649 * Return total three momentum.
650 */
651 G4double getTotalThreeMomentum() {
652 G4double momentum = 0;
653 for(G4int i = 0; i < ntrack; i++) {
654 momentum += plab[i];
655 }
656 return momentum;
657 }
658
659 G4double getMomentumSum() {
660 G4double momentum = 0;
661 for(G4int i = 0; i < ntrack; i++) {
662 momentum += plab[i];
663 }
664 return momentum;
665 }
666
667 G4double getMass(G4int particle) {
668 const G4double protonMass = 938.272;
669 const G4double neutronMass = 939.565;
670 const G4double pionMass = 139.57;
671
672 G4double mass = 0.0;
673 if(avv[particle] == 1 && zvv[particle] == 1) mass = protonMass;
674 if(avv[particle] == 1 && zvv[particle] == 0) mass = neutronMass;
675 if(avv[particle] == -1) mass = pionMass;
676 if(avv[particle] > 1)
677 mass = avv[particle] * protonMass + zvv[particle] * neutronMass;
678 return mass;
679 }
680
681 /**
682 * Dump debugging output.
683 */
684 void dump() {
685 G4int nProton = 0, nNeutron = 0;
686 G4int nPiPlus = 0, nPiZero = 0, nPiMinus = 0;
687 G4int nH2 = 0, nHe3 = 0, nAlpha = 0;
688 G4int nFragments = 0;
689 G4int nParticles = 0;
690 G4cout <<"Particles produced in the event (" << ntrack << "):" << G4endl;
691 G4cout <<"A \t Z \t Ekin \t Ptot \t Theta \t Phi" << G4endl;
692 for(G4int i = 0; i < ntrack; i++) {
693 nParticles++;
694 if(avv[i] == 1 && zvv[i] == 1) nProton++; // Count multiplicities
695 if(avv[i] == 1 && zvv[i] == 0) nNeutron++;
696 if(avv[i] == -1 && zvv[i] == 1) nPiPlus++;
697 if(avv[i] == -1 && zvv[i] == 0) nPiZero++;
698 if(avv[i] == -1 && zvv[i] == -1) nPiMinus++;
699 if(avv[i] == 2 && zvv[i] == 1) nH2++;
700 if(avv[i] == 3 && zvv[i] == 2) nHe3++;
701 if(avv[i] == 4 && zvv[i] == 2) nAlpha++;
702 if( zvv[i] > 2) nFragments++;
703
704 G4cout << i << " \t " << avv[i] << " \t " << zvv[i] << " \t " << enerj[i] << " \t "
705 << plab[i] << " \t " << tetlab[i] << " \t " << philab[i] << G4endl;
706 }
707
708 G4cout <<"Summary of event: " << G4endl;
709 G4cout <<"Projectile type: " << projType <<" Energy: " << projEnergy << G4endl;
710 G4cout <<"Target A = " << targetA << " Z = " << targetZ << G4endl;
711 G4cout <<"Remnant from cascade: " << G4endl;
712 G4cout <<"A = " << massini << " Z = " << mzini << " excitation E = " << exini << G4endl;
713 G4cout <<"Particle multiplicities:" << G4endl;
714 G4cout <<"Protons: " << nProton << " Neutrons: " << nNeutron << G4endl;
715 G4cout <<"pi+: " << nPiPlus << " pi0: " << nPiZero << " pi-: " << nPiMinus << G4endl;
716 G4cout <<"H2: " << nH2 << " He3: " << nHe3 << " Alpha: " << nAlpha << G4endl;
717 G4cout <<"Nucleus fragments = " << nFragments << G4endl;
718 G4cout <<"Conservation laws:" << G4endl;
719 G4cout <<"Baryon number = " << getTotalBaryonNumber() << G4endl;
720 G4cout <<"Number of particles = " << nParticles << G4endl;
721 }
722
723 /**
724 * Projectile type.
725 */
726 G4int projType;
727
728 /**
729 * Projectile energy.
730 */
731 G4double projEnergy;
732
733 /**
734 * Target mass number.
735 */
736 G4int targetA;
737
738 /**
739 * Target charge number.
740 */
741 G4int targetZ;
742
743 /**
[819]744 * A of the remnant.
745 */
746 G4double massini;
747
748 /**
749 * Z of the remnant.
750 */
751 G4double mzini;
752
753 /**
754 * Excitation energy.
755 */
756 G4double exini;
757
[962]758 G4double pcorem, mcorem, pxrem, pyrem, pzrem;
759
[819]760 /**
761 * Cascade n multip.
762 */
763 G4int mulncasc;
764
765 /**
766 * Evaporation n multip.
767 */
768 G4int mulnevap;
769
770 /**
771 * Total n multip.
772 */
773 G4int mulntot;
774
775 /**
776 * Impact parameter.
777 */
778 G4double bimpact;
779
780 /**
781 * Remnant Intrinsic Spin.
782 */
783 G4int jremn;
784
785 /**
786 * Fission 1/0=Y/N.
787 */
788 G4int kfis;
789
790 /**
791 * Excit energy at fis.
792 */
793 G4double estfis;
794
795 /**
796 * Z of fiss nucleus.
797 */
798 G4int izfis;
799
800 /**
801 * A of fiss nucleus.
802 */
803 G4int iafis;
804
805 /**
806 * Number of particles.
807 */
808 G4int ntrack;
809
810 /**
[962]811 * The state of the index:
812 * true = reserved
813 * false = free
814 */
815 G4bool full[VARNTPSIZE];
816
817 /**
[819]818 * emitted in cascade (0) or evaporation (1).
819 */
820 G4int itypcasc[VARNTPSIZE];
821
822
823 /**
824 * A (-1 for pions).
825 */
826 G4int avv[VARNTPSIZE];
827
828 /**
829 * Z
830 */
831 G4int zvv[VARNTPSIZE];
832
833 /**
834 * Kinetic energy.
835 */
836 G4double enerj[VARNTPSIZE];
837
838 /**
839 * Momentum.
840 */
841 G4double plab[VARNTPSIZE];
842
843 /**
844 * Theta angle.
845 */
846 G4double tetlab[VARNTPSIZE];
847
848 /**
849 * Phi angle.
850 */
851 G4double philab[VARNTPSIZE];
[962]852
853private:
854 G4int particleIndex;
[819]855};
856
857/**
858 * Pauli blocking.
859 */
860class G4Paul {
861public:
862 G4Paul() {};
863 ~G4Paul() {};
864
865 /**
866 *
867 */
868 G4double ct0,ct1,ct2,ct3,ct4,ct5,ct6,pr,pr2,xrr,xrr2;
869
870 /**
871 *
872 */
873 G4double cp0,cp1,cp2,cp3,cp4,cp5,cp6;
874};
875
876
877#endif
Note: See TracBrowser for help on using the repository browser.