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

Last change on this file since 1199 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 15.7 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: G4InclDataDefs.hh,v 1.7 2009/11/18 10:43:14 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// 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 void dump() {
326 G4cout <<"Avatars: (number of avatars = " << k << ")" << G4endl;
327 for(G4int i = 0; i <= k; i++) {
328 G4cout <<"i = " << i << G4endl;
329 G4cout <<"crois[" << i << "] = " << crois[i] << G4endl;
330 G4cout <<"ind[" << i << "] = " << ind[i] << G4endl;
331 G4cout <<"jnd[" << i << "] = " << jnd[i] << G4endl;
332 }
333 }
334
335 /**
336 *
337 */
338 G4double crois[BL2CROISSIZE];
339
340 /**
341 *
342 */
343 G4int k;
344
345 /**
346 *
347 */
348 G4int ind[BL2INDSIZE];
349
350 /**
351 *
352 */
353 G4int jnd[BL2INDSIZE];
354};
355
356//#define BL3SIZE 300
357#define BL3SIZE 3000
358/**
359 *
360 */
361class G4Bl3 {
362public:
363 G4Bl3() {};
364 ~G4Bl3() {};
365
366 /**
367 * r1 and r2
368 */
369 G4double r1,r2;
370
371 /**
372 * Nucleon positions
373 */
374 G4double x1[BL3SIZE], x2[BL3SIZE],x3[BL3SIZE];
375
376 /**
377 * Mass numbers
378 */
379 G4int ia1,ia2;
380
381 /**
382 * rab2
383 */
384 G4double rab2;
385};
386
387/**
388 * G4Bl4
389 */
390class G4Bl4 {
391public:
392 G4Bl4() {};
393 ~G4Bl4() {};
394
395 /**
396 * tmax5
397 */
398 G4double tmax5;
399};
400
401//#define BL5SIZE 300
402#define BL5SIZE 3000
403/**
404 * G4Bl5
405 */
406class G4Bl5 {
407public:
408 G4Bl5() {};
409 ~G4Bl5() {};
410
411 /**
412 * tlg
413 */
414 G4double tlg[BL5SIZE];
415
416 /**
417 * nesc
418 */
419 G4int nesc[BL5SIZE];
420};
421
422/**
423 * G4Bl6
424 */
425class G4Bl6 {
426public:
427 G4Bl6() {};
428 ~G4Bl6() {};
429
430 /**
431 * xx10
432 */
433 G4double xx10;
434
435 /**
436 * isa
437 */
438 G4double isa;
439};
440
441/**
442 * G4Bl8
443 */
444class G4Bl8 {
445public:
446 G4Bl8() {};
447 ~G4Bl8() {};
448
449 /**
450 * rathr
451 */
452 G4double rathr;
453
454 /**
455 * ramass
456 */
457 G4double ramass;
458};
459
460//#define BL9SIZE 300
461#define BL9SIZE 3000
462/**
463 * G4Bl9
464 */
465class G4Bl9 {
466public:
467 G4Bl9() {
468 l1 = 0;
469 l2 = 0;
470 };
471 ~G4Bl9() {};
472
473 /**
474 * hel
475 */
476 G4double hel[BL9SIZE];
477
478 /**
479 * l1 and l2
480 */
481 G4int l1,l2;
482};
483
484/**
485 * G4Bl10
486 */
487class G4Bl10 {
488public:
489 G4Bl10() {};
490 ~G4Bl10() {};
491
492 /**
493 * ri4, rs4, r2i, r2s, pdummy, pf
494 */
495 G4double ri4,rs4,r2i,r2s,pdummy,pf;
496};
497
498/**
499 * G4Kind
500 */
501class G4Kind {
502public:
503 G4Kind() {};
504 ~G4Kind() {};
505
506 /**
507 * kindf7
508 */
509 G4int kindf7;
510};
511
512#define VARSIZE 3
513#define VAEPSSIZE 250
514#define VAAVM 1000
515/**
516 * Extra information on collisions between nucleons.
517 */
518class G4VarAvat {
519public:
520 G4VarAvat() {};
521 ~G4VarAvat() {};
522
523 /**
524 *
525 */
526 G4int kveux;
527
528 /**
529 *
530 */
531 G4double bavat;
532
533 /**
534 *
535 */
536 G4int nopartavat,ncolavat;
537
538 /**
539 *
540 */
541 G4double r1_in[VARSIZE],r1_first_avat[VARSIZE];
542
543 /**
544 *
545 */
546 G4double epsd[VAEPSSIZE],eps2[VAEPSSIZE],eps4[VAEPSSIZE],eps6[VAEPSSIZE],epsf[VAEPSSIZE];
547
548 /**
549 *
550 */
551 G4int nb_avat;
552
553 /**
554 *
555 */
556 G4double timeavat[VAAVM],l1avat[VAAVM],l2avat[VAAVM],jpartl1[VAAVM],jpartl2[VAAVM];
557
558 /**
559 *
560 */
561 G4double del1avat[VAAVM],del2avat[VAAVM],energyavat[VAAVM];
562
563 /**
564 *
565 */
566 G4double bloc_paul[VAAVM],bloc_cdpp[VAAVM],go_out[VAAVM];
567};
568
569#define VARNTPSIZE 255
570class G4VarNtp {
571public:
572 G4VarNtp() {};
573 ~G4VarNtp() {};
574
575 /**
576 * Clear and initialize all variables and arrays.
577 */
578 void clear() {
579 particleIndex = 0;
580 projType = 0;
581 projEnergy = 0.0;
582 targetA = 0;
583 targetZ = 0;
584 massini = 0;
585 mzini = 0;
586 exini = 0;
587 pcorem = 0;
588 mcorem = 0;
589 pxrem = 0;
590 pyrem = 0;
591 pzrem = 0;
592 mulncasc = 0;
593 mulnevap = 0;
594 mulntot = 0;
595 bimpact = 0.0;
596 jremn = 0;
597 kfis = 0;
598 estfis = 0;
599 izfis = 0;
600 iafis = 0;
601 ntrack = 0;
602 for(G4int i = 0; i < VARNTPSIZE; i++) {
603 itypcasc[i] = 0;
604 avv[i] = 0;
605 zvv[i] = 0;
606 enerj[i] = 0.0;
607 plab[i] = 0.0;
608 tetlab[i] = 0.0;
609 philab[i] = 0.0;
610 full[i] = false;
611 }
612 }
613
614 void addParticle(G4double A, G4double Z, G4double E, G4double P, G4double theta, G4double phi) {
615 if(full[particleIndex]) {
616 G4cout <<"G4VarNtp: Error. Index i = " << particleIndex << " is already occupied by particle:" << G4endl;
617 G4cout <<"A = " << avv[particleIndex] << " Z = " << zvv[particleIndex] << G4endl;
618 G4cout <<"Tried to replace it with:" << G4endl;
619 G4cout <<"A = " << Z << " Z = " << Z << G4endl;
620 } else {
621 avv[particleIndex] = (int) A;
622 zvv[particleIndex] = (int) Z;
623 enerj[particleIndex] = E;
624 plab[particleIndex] = P;
625 tetlab[particleIndex] = theta;
626 philab[particleIndex] = phi;
627 full[particleIndex] = true;
628 ntrack = particleIndex + 1;
629 particleIndex++;
630 }
631 }
632
633 /**
634 * Baryon number conservation check.
635 */
636 G4int getTotalBaryonNumber() {
637 G4int baryonNumber = 0;
638 for(G4int i = 0; i < ntrack; i++) {
639 if(avv[i] > 0) {
640 baryonNumber += avv[i];
641 }
642 }
643 return baryonNumber;
644 }
645
646 /**
647 * Return total energy.
648 */
649 G4double getTotalEnergy() {
650 G4double energy = 0.0;
651 for(G4int i = 0; i < ntrack; i++) {
652 energy += std::sqrt(std::pow(plab[i], 2) + std::pow(getMass(i), 2)); // E^2 = p^2 + m^2
653 }
654
655 return energy;
656 }
657
658 /**
659 * Return total three momentum.
660 */
661 G4double getTotalThreeMomentum() {
662 G4double momentum = 0;
663 for(G4int i = 0; i < ntrack; i++) {
664 momentum += plab[i];
665 }
666 return momentum;
667 }
668
669 G4double getMomentumSum() {
670 G4double momentum = 0;
671 for(G4int i = 0; i < ntrack; i++) {
672 momentum += plab[i];
673 }
674 return momentum;
675 }
676
677 G4double getMass(G4int particle) {
678 const G4double protonMass = 938.272;
679 const G4double neutronMass = 939.565;
680 const G4double pionMass = 139.57;
681
682 G4double mass = 0.0;
683 if(avv[particle] == 1 && zvv[particle] == 1) mass = protonMass;
684 if(avv[particle] == 1 && zvv[particle] == 0) mass = neutronMass;
685 if(avv[particle] == -1) mass = pionMass;
686 if(avv[particle] > 1)
687 mass = avv[particle] * protonMass + zvv[particle] * neutronMass;
688 return mass;
689 }
690
691 /**
692 * Dump debugging output.
693 */
694 void dump() {
695 G4int nProton = 0, nNeutron = 0;
696 G4int nPiPlus = 0, nPiZero = 0, nPiMinus = 0;
697 G4int nH2 = 0, nHe3 = 0, nAlpha = 0;
698 G4int nFragments = 0;
699 G4int nParticles = 0;
700 G4cout <<"Particles produced in the event (" << ntrack << "):" << G4endl;
701 G4cout <<"A \t Z \t Ekin \t Ptot \t Theta \t Phi" << G4endl;
702 for(G4int i = 0; i < ntrack; i++) {
703 nParticles++;
704 if(avv[i] == 1 && zvv[i] == 1) nProton++; // Count multiplicities
705 if(avv[i] == 1 && zvv[i] == 0) nNeutron++;
706 if(avv[i] == -1 && zvv[i] == 1) nPiPlus++;
707 if(avv[i] == -1 && zvv[i] == 0) nPiZero++;
708 if(avv[i] == -1 && zvv[i] == -1) nPiMinus++;
709 if(avv[i] == 2 && zvv[i] == 1) nH2++;
710 if(avv[i] == 3 && zvv[i] == 2) nHe3++;
711 if(avv[i] == 4 && zvv[i] == 2) nAlpha++;
712 if( zvv[i] > 2) nFragments++;
713
714 G4cout << i << " \t " << avv[i] << " \t " << zvv[i] << " \t " << enerj[i] << " \t "
715 << plab[i] << " \t " << tetlab[i] << " \t " << philab[i] << G4endl;
716 }
717
718 G4cout <<"Summary of event: " << G4endl;
719 G4cout <<"Projectile type: " << projType <<" Energy: " << projEnergy << G4endl;
720 G4cout <<"Target A = " << targetA << " Z = " << targetZ << G4endl;
721 G4cout <<"Remnant from cascade: " << G4endl;
722 G4cout <<"A = " << massini << " Z = " << mzini << " excitation E = " << exini << G4endl;
723 G4cout <<"Particle multiplicities:" << G4endl;
724 G4cout <<"Protons: " << nProton << " Neutrons: " << nNeutron << G4endl;
725 G4cout <<"pi+: " << nPiPlus << " pi0: " << nPiZero << " pi-: " << nPiMinus << G4endl;
726 G4cout <<"H2: " << nH2 << " He3: " << nHe3 << " Alpha: " << nAlpha << G4endl;
727 G4cout <<"Nucleus fragments = " << nFragments << G4endl;
728 G4cout <<"Conservation laws:" << G4endl;
729 G4cout <<"Baryon number = " << getTotalBaryonNumber() << G4endl;
730 G4cout <<"Number of particles = " << nParticles << G4endl;
731 }
732
733 /**
734 * Projectile type.
735 */
736 G4int projType;
737
738 /**
739 * Projectile energy.
740 */
741 G4double projEnergy;
742
743 /**
744 * Target mass number.
745 */
746 G4int targetA;
747
748 /**
749 * Target charge number.
750 */
751 G4int targetZ;
752
753 /**
754 * A of the remnant.
755 */
756 G4double massini;
757
758 /**
759 * Z of the remnant.
760 */
761 G4double mzini;
762
763 /**
764 * Excitation energy.
765 */
766 G4double exini;
767
768 G4double pcorem, mcorem, pxrem, pyrem, pzrem;
769
770 /**
771 * Cascade n multip.
772 */
773 G4int mulncasc;
774
775 /**
776 * Evaporation n multip.
777 */
778 G4int mulnevap;
779
780 /**
781 * Total n multip.
782 */
783 G4int mulntot;
784
785 /**
786 * Impact parameter.
787 */
788 G4double bimpact;
789
790 /**
791 * Remnant Intrinsic Spin.
792 */
793 G4int jremn;
794
795 /**
796 * Fission 1/0=Y/N.
797 */
798 G4int kfis;
799
800 /**
801 * Excit energy at fis.
802 */
803 G4double estfis;
804
805 /**
806 * Z of fiss nucleus.
807 */
808 G4int izfis;
809
810 /**
811 * A of fiss nucleus.
812 */
813 G4int iafis;
814
815 /**
816 * Number of particles.
817 */
818 G4int ntrack;
819
820 /**
821 * The state of the index:
822 * true = reserved
823 * false = free
824 */
825 G4bool full[VARNTPSIZE];
826
827 /**
828 * emitted in cascade (0) or evaporation (1).
829 */
830 G4int itypcasc[VARNTPSIZE];
831
832
833 /**
834 * A (-1 for pions).
835 */
836 G4int avv[VARNTPSIZE];
837
838 /**
839 * Z
840 */
841 G4int zvv[VARNTPSIZE];
842
843 /**
844 * Kinetic energy.
845 */
846 G4double enerj[VARNTPSIZE];
847
848 /**
849 * Momentum.
850 */
851 G4double plab[VARNTPSIZE];
852
853 /**
854 * Theta angle.
855 */
856 G4double tetlab[VARNTPSIZE];
857
858 /**
859 * Phi angle.
860 */
861 G4double philab[VARNTPSIZE];
862
863private:
864 G4int particleIndex;
865};
866
867/**
868 * Pauli blocking.
869 */
870class G4Paul {
871public:
872 G4Paul() {};
873 ~G4Paul() {};
874
875 /**
876 *
877 */
878 G4double ct0,ct1,ct2,ct3,ct4,ct5,ct6,pr,pr2,xrr,xrr2;
879
880 /**
881 *
882 */
883 G4double cp0,cp1,cp2,cp3,cp4,cp5,cp6;
884};
885
886
887#endif
Note: See TracBrowser for help on using the repository browser.