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

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

File size: 21.8 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.12 2010/11/17 20:19:09 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#include "G4Nucleus.hh"
39#include "G4HadProjectile.hh"
40#include "G4ParticleTable.hh"
41#include "G4Track.hh"
42
43class G4InclFermi {
44public:
45  G4InclFermi() {
46    G4double hc = 197.328;
47    G4double fmp = 938.2796;
48    pf=1.37*hc;
49    pf2=pf*pf;
50    tf=std::sqrt(pf*pf+fmp*fmp)-fmp;
51  };
52  ~G4InclFermi() {};
53
54  G4double tf,pf,pf2;
55};
56
57#define max_a_proj 61
58
59/**
60 * (eps_c,p1_s,p2_s,p3_s,eps_c used to store the kinematics of
61 * nucleons for composit projectiles before entering the potential)
62 */
63class G4QuadvectProjo {
64public:
65  G4QuadvectProjo() {
66    for(G4int i = 0; i < max_a_proj; ++i) {
67      eps_c[i] = 0.0;
68      t_c[i] = 0.0;
69      p3_c[i] = 0.0;
70      p1_s[i] = 0.0;
71      p2_s[i] = 0.0;
72      p3_s[i] = 0.0;
73    }
74  };
75
76  ~G4QuadvectProjo() {};
77
78  G4double eps_c[max_a_proj],p3_c[max_a_proj],
79    p1_s[max_a_proj],p2_s[max_a_proj],p3_s[max_a_proj],
80    t_c[max_a_proj];
81};
82
83class G4VBe {
84public:
85  G4VBe()
86    :ia_be(0), iz_be(0),
87     rms_be(0.0), pms_be(0.0), bind_be(0.0)
88  { };
89
90  ~G4VBe() {};
91
92  G4int ia_be, iz_be;
93  G4double rms_be, pms_be, bind_be;
94};
95
96/**
97 * Projectile spectator
98 */
99class G4InclProjSpect {
100public:
101  G4InclProjSpect() {
102    //    G4cout <<"Projectile spectator data structure created!" << G4endl;
103    clear();
104  };
105  ~G4InclProjSpect() {};
106
107  void clear() {
108    for(G4int i = 0; i < 21; i++) tab[i] = 0.0;
109    for(G4int i = 0; i < 61; i++) n_projspec[i] = 0;
110    a_projspec = 0;
111    z_projspec = 0;
112    t_projspec = 0.0;
113    ex_projspec = 0.0;
114    p1_projspec = 0.0;
115    p2_projspec = 0.0;
116    p3_projspec = 0.0;
117    m_projspec = 0.0;
118  };
119
120  G4double tab[21];
121  G4int n_projspec[61];
122  G4int a_projspec,z_projspec;
123  G4double ex_projspec,t_projspec, p1_projspec, p2_projspec, p3_projspec, m_projspec;
124};
125
126#define IGRAINESIZE 19
127/**
128 * Random seeds used by internal random number generators.
129 * @see G4Incl::standardRandom
130 * @see G4Incl::gaussianRandom
131 */
132class G4Hazard{
133public:
134  G4Hazard() {
135    ial = 0;
136    for(G4int i = 0; i < IGRAINESIZE; ++i) {
137      igraine[i] = 0;
138    }
139  };
140
141  ~G4Hazard() {};
142
143  /**
144   * Random seed
145   */ 
146  G4long ial;
147
148  /**
149   * An array of random seeds.
150   */
151  G4long igraine[IGRAINESIZE];
152};
153
154#define MATSIZE 500
155#define MATGEOSIZE 6
156/**
157 * Target nuclei to be taken into account in the cascade problem.
158 */
159class G4Mat {
160public:
161  G4Mat() {
162    nbmat = 0;
163
164    for(G4int i = 0; i < MATSIZE; ++i) {
165      zmat[i] = 0;
166      amat[i] = 0;
167    }
168
169    for(G4int i = 0; i < MATGEOSIZE; ++i) {
170      for(G4int j = 0; j < MATSIZE; ++j) {
171        bmax_geo[i][j] = 0;
172      }
173    }
174};
175
176  ~G4Mat() { };
177
178  /**
179   * Charge numbers.
180   */
181  G4int zmat[MATSIZE];
182
183  /**
184   * Mass number
185   */ 
186  G4int amat[MATSIZE];
187
188  /**
189   *
190   */
191  G4double bmax_geo[MATGEOSIZE][MATSIZE];
192
193  /**
194   * Number of materials.
195   */
196  G4int nbmat;
197};
198
199#define LGNSIZE 9
200/**
201 * Properties of light nucleus used as a bullet.
202 */
203class G4LightGausNuc {
204public:
205  G4LightGausNuc() {
206    for(G4int i = 0; i < LGNSIZE; ++i) {
207      rms1t[i] = 0.0;
208      pf1t[i] = 0.0;
209      pfln[i] = 0.0;
210      tfln[i] = 0.0;
211      vnuc[i] = 0.0;
212    }
213  };
214
215  ~G4LightGausNuc() {};
216 
217  G4double rms1t[LGNSIZE];
218  G4double pf1t[LGNSIZE];
219  G4double pfln[LGNSIZE];
220  G4double tfln[LGNSIZE];
221  G4double vnuc[LGNSIZE];
222};
223
224#define LNSIZE 30
225/**
226 * Data of light nuclei.
227 */
228class G4LightNuc {
229public:
230  G4LightNuc() {
231    for(G4int i = 0; i < LNSIZE; ++i) {
232      r[i] = 0.0;
233      a[i] = 0.0;
234    }
235  };
236
237  ~G4LightNuc() {};
238
239  /**
240   * r
241   */
242  G4double r[LNSIZE];
243
244  /**
245   * a
246   */
247  G4double a[LNSIZE];
248};
249
250#define SAXWROWS 30
251#define SAXWCOLS 500
252/**
253 * Woods-Saxon density and its first derivative.
254 */
255class G4Saxw {
256public:
257  G4Saxw() {
258    for(G4int i = 0; i < SAXWROWS; ++i) {
259      for(G4int j = 0; j < SAXWCOLS; ++j) {
260        x[i][j] = 0.0;
261        y[i][j] = 0.0;
262        s[i][j] = 0.0;
263      }
264    }
265    imat = 0; n = 0; k = 0;
266  };
267
268  ~G4Saxw() {};
269 
270  /**
271   * x
272   */
273  G4double x[SAXWROWS][SAXWCOLS]; 
274
275  /**
276   * y
277   */
278  G4double y[SAXWROWS][SAXWCOLS]; 
279
280  /**
281   * s
282   */
283  G4double s[SAXWROWS][SAXWCOLS]; 
284
285  /**
286   * imat
287   */
288  G4int imat;
289
290  /**
291   * n
292   */
293  G4int n;
294
295  /**
296   * k
297   */
298  G4int k;
299};
300
301/**
302 * Parameters for INCL4 model.
303 */
304class G4Ws {
305public:
306  G4Ws() {
307    fneck = 0.0;
308    r0 = 0.0;
309    adif = 0.0;
310    rmaxws = 0.0;
311    drws = 0.0;
312    nosurf = 0.0;
313    xfoisa = 0.0;
314    bmax = 0.0;
315    npaulstr = 0.0;
316  };
317
318  ~G4Ws() {};
319 
320  /**
321   * r0
322   */
323  G4double r0;
324
325  /**
326   * adif
327   */
328  G4double adif;
329
330  /**
331   * Maximum radius of the nucleus
332   */
333  G4double rmaxws;
334
335  /**
336   * drws
337   */
338  G4double drws;
339
340  /**
341   * Shape of the surface of the nucleus:
342   * - -1: Woods-Saxon density with impact parameter dependence
343   * -  0: Woods-Saxon density without impact parameter dependence
344   * -  1: Sharp surface (hard sphere)
345   */
346  G4double nosurf;
347
348  /**
349   * Parameter related to the maximum radius of the nucleus.
350   *
351   * rmaxws = r0 + xfoisa*A
352   */
353  G4double xfoisa;
354
355  /**
356   * Pauli blocking used in the simulation:
357   * - 0: statistic Pauli blocking
358   * - 1: strict Pauli blocking
359   * - 2: no Pauli blocking
360   */
361  G4double npaulstr;
362
363  /**
364   * Maximum impact parameter
365   */
366  G4double bmax;
367
368  G4double fneck;
369};
370
371#define DTONSIZE 13
372/**
373 * Random seeds used by internal random number generators.
374 * @see G4Incl::standardRandom
375 * @see G4Incl::gaussianRandom
376 */
377class G4Dton {
378public:
379  G4Dton() {
380    fn = 0.0;
381    for(G4int i = 0; i < DTONSIZE; ++i) {
382      c[i] = 0.0;
383      d[i] = 0.0;
384    }
385  };
386
387  ~G4Dton() {};
388 
389  G4double c[DTONSIZE];
390  G4double d[DTONSIZE];
391  G4double fn;
392};
393
394#define SPL2SIZE 100
395/**
396 * Random seeds used by internal random number generators.
397 * @see G4Incl::standardRandom
398 * @see G4Incl::gaussianRandom
399 */
400class G4Spl2 {
401public:
402  G4Spl2() {
403    for(G4int i = 0; i < SPL2SIZE; ++i) {
404      x[i] = 0.0; y[i] = 0.0;
405      a[i] = 0.0; b[i] = 0.0; c[i] = 0.0;
406    }
407    n = 0;
408  };
409
410  ~G4Spl2() {};
411 
412  G4double x[SPL2SIZE];
413  G4double y[SPL2SIZE];
414  G4double a[SPL2SIZE];
415  G4double b[SPL2SIZE];
416  G4double c[SPL2SIZE];
417  G4int n;
418};
419
420// incl4.2.cc:
421
422//#define BL1SIZE 300
423#define BL1SIZE 3000
424/**
425 * Random seeds used by internal random number generators.
426 * @see G4Incl::standardRandom
427 * @see G4Incl::gaussianRandom
428 */
429class G4Bl1 {
430public:
431  G4Bl1() {
432    ta = 0.0;
433    for(G4int i = 0; i < BL1SIZE; ++i) {
434      p1[i] = 0.0; p2[i] = 0.0; p3[i] = 0.0; eps[i] = 0.0;
435      ind1[i] = 0; ind2[i] = 0;
436    }
437  };
438
439  ~G4Bl1() {};
440 
441  G4double p1[BL1SIZE],p2[BL1SIZE],p3[BL1SIZE];
442  G4double eps[BL1SIZE];
443  G4int ind1[BL1SIZE],ind2[BL1SIZE];
444  G4double ta;
445
446  void dump(G4int numberOfParticles) {
447    static G4int dumpNumber = 0;
448    G4cout <<"Dump number" << dumpNumber << " of particle 4-momenta (G4Bl1):" << G4endl;
449    G4cout <<"ta = " << ta << G4endl;
450    for(G4int i = 0; i < numberOfParticles; i++) {
451      G4cout <<"i = " << i << "   p1 = " << p1[i] << "   p2 = " << p2[i] << "   p3 = " << p3[i] << "   eps = " << eps[i] << G4endl;
452    }
453    dumpNumber++;
454  }
455};
456
457#define BL2SIZE 19900
458/**
459 *
460 */
461class G4Bl2 {
462public:
463  G4Bl2() {
464    k = 0;
465    for(G4int i = 0; i < BL2SIZE; ++i) {
466      crois[i] = 0.0;
467      ind[i] = 0;
468      jnd[i] = 0;
469    }   
470  };
471
472  ~G4Bl2() {};
473 
474  void dump() {
475    G4cout <<"Avatars: (number of avatars = " << k << ")" << G4endl;
476    for(G4int i = 0; i <= k; i++) {
477      G4cout <<"i = " << i << G4endl;
478      G4cout <<"crois[" << i << "] = " << crois[i] << G4endl;
479      G4cout <<"ind[" << i << "] = " << ind[i] << G4endl;
480      G4cout <<"jnd[" << i << "] = " << jnd[i] << G4endl;
481    }
482  }
483
484  /**
485   *
486   */
487  G4double crois[BL2SIZE];
488
489  /**
490   *
491   */
492  G4int k;
493
494  /**
495   *
496   */
497  G4int ind[BL2SIZE];
498
499  /**
500   *
501   */
502  G4int jnd[BL2SIZE];
503};
504
505//#define BL3SIZE 300
506#define BL3SIZE 3000
507/**
508 *
509 */
510class G4Bl3 {
511public:
512  G4Bl3() {
513    r1 = 0.0; r2 = 0.0;
514    ia1 = 0; ia2 = 0;
515    rab2 = 0.0;
516
517    for(G4int i = 0; i < BL3SIZE; ++i) {
518      x1[i] = 0.0; x2[i] = 0.0; x3[i] = 0.0;
519    }
520  };
521
522  ~G4Bl3() {};
523 
524  /**
525   * r1 and r2
526   */
527  G4double r1,r2;
528
529  /**
530   * Nucleon positions
531   */
532  G4double x1[BL3SIZE], x2[BL3SIZE],x3[BL3SIZE];
533
534  /**
535   * Mass numbers
536   */
537  G4int ia1,ia2;
538
539  /**
540   * rab2
541   */
542  G4double rab2;
543
544  void dump() {
545    static G4int dumpNumber = 0;
546    G4cout <<"Dump number" << dumpNumber << " of particle positions (G4Bl3):" << G4endl;
547    G4cout <<" ia1 = " << ia1 << G4endl;
548    G4cout <<" ia2 = " << ia2 << G4endl;
549    G4cout <<" rab2 = " << rab2 << G4endl;
550    for(G4int i = 0; i <= (ia1 + ia2); i++) {
551      G4cout <<"i = " << i << "   x1 = " << x1[i] << "   x2 = " << x2[i] << "   x3 = " << x3[i] << G4endl;
552    }
553    dumpNumber++;
554  }
555};
556
557/**
558 * G4Bl4
559 */
560class G4Bl4 {
561public:
562  G4Bl4()
563    :tmax5(0.0)
564  {};
565
566  ~G4Bl4() {};
567
568  /**
569   * tmax5
570   */
571  G4double tmax5;
572};
573
574//#define BL5SIZE 300
575#define BL5SIZE 3000
576/**
577 * G4Bl5
578 */
579class G4Bl5 {
580public:
581  G4Bl5() {
582    for(G4int i = 0; i < BL5SIZE; ++i) {
583      tlg[i] = 0.0;
584      nesc[i] = 0;
585    }
586  };
587
588  ~G4Bl5() {};
589 
590  /**
591   * tlg
592   */
593  G4double tlg[BL5SIZE];
594
595  /**
596   * nesc
597   */
598  G4int nesc[BL5SIZE];
599};
600
601/**
602 * G4Bl6
603 */
604class G4Bl6 {
605public:
606  G4Bl6()
607    :xx10(0.0), isa(0.0)
608  {};
609
610  ~G4Bl6() {};
611 
612  /**
613   * xx10
614   */
615  G4double xx10;
616
617  /**
618   * isa
619   */
620  G4double isa;
621};
622
623/**
624 * G4Bl8
625 */
626class G4Bl8 {
627public:
628  G4Bl8()
629    :rathr(0.0), ramass(0.0)
630  {};
631
632  ~G4Bl8() {};
633
634  /**
635   * rathr
636   */
637  G4double rathr;
638
639  /**
640   * ramass
641   */
642  G4double ramass;
643};
644
645//#define BL9SIZE 300
646#define BL9SIZE 3000
647/**
648 * G4Bl9
649 */
650class G4Bl9 {
651public:
652  G4Bl9() {
653    l1 = 0;
654    l2 = 0;
655    for(G4int i = 0; i < BL9SIZE; ++i) {
656      hel[i] = 0.0;
657    }
658  };
659  ~G4Bl9() {};
660
661  /**
662   * hel
663   */
664  G4double hel[BL9SIZE];
665
666  /**
667   * l1 and l2
668   */
669  G4int l1,l2;
670};
671
672/**
673 * G4Bl10
674 */
675class G4Bl10 {
676public:
677  G4Bl10()
678    :ri4(0.0), rs4(0.0), r2i(0.0), r2s(0.0), pdummy(0.0), pf(0.0)
679  {};
680  ~G4Bl10() {};
681
682  /**
683   * ri4, rs4, r2i, r2s, pdummy, pf
684   */
685  G4double ri4,rs4,r2i,r2s,pdummy,pf;
686};
687
688/**
689 * G4Kind
690 */
691class G4Kind {
692public:
693  G4Kind()
694    :kindf7(0)
695  {};
696  ~G4Kind() {};
697
698  /**
699   * kindf7
700   */
701  G4int kindf7;
702};
703
704/**
705 * Projectile parameters.
706 */
707class G4Bev {
708public:
709  /**
710   * Initialize all variables to zero.
711   */
712  G4Bev() {
713    ia_be = 0;
714    iz_be = 0;
715    rms_be = 0.0;
716    pms_be = 0.0;
717    bind_be = 0.0;
718  };
719  ~G4Bev() {};
720
721  /**
722   * Mass number.
723   */
724  G4int ia_be;
725
726  /**
727   * Charge number.
728   */
729  G4int iz_be;
730
731  /**
732   * rms
733   */
734  G4double rms_be;
735
736  /**
737   * pms
738   */
739  G4double pms_be;
740
741  /**
742   * bind
743   */
744  G4double bind_be;
745};
746
747#define VARSIZE 3
748#define VAEPSSIZE 250
749#define VAAVM 1000
750/**
751 * Extra information on collisions between nucleons.
752 */
753class G4VarAvat {
754public:
755  G4VarAvat() {
756    kveux = 0;
757    bavat = 0.0;
758    nopartavat = 0; ncolavat = 0;
759    nb_avat = 0;
760
761    for(G4int i = 0; i < VARSIZE; ++i) {
762      r1_in[i] = 0.0;
763      r1_first_avat[i] = 0.0;
764    }
765
766    for(G4int i = 0; i < VAEPSSIZE; ++i) {
767      epsd[i] = 0.0;
768      eps2[i] = 0.0;
769      eps4[i] = 0.0;
770      eps6[i] = 0.0;
771      epsf[i] = 0.0;
772    }
773
774    for(G4int i = 0; i < VAAVM; ++i) {
775      timeavat[i] = 0.0;
776      l1avat[i] = 0.0;
777      l2avat[i] = 0.0;
778      jpartl1[i] = 0.0;
779      jpartl2[i] = 0.0;
780      del1avat[i] = 0.0;
781      del2avat[i] = 0.0;
782      energyavat[i] = 0.0;
783      bloc_paul[i] = 0.0;
784      bloc_cdpp[i] = 0.0;
785      go_out[i] = 0.0;
786    }
787  };
788
789  ~G4VarAvat() {};
790
791  /**
792   *
793   */
794  G4int kveux;
795
796  /**
797   *
798   */
799  G4double bavat;
800
801  /**
802   *
803   */
804  G4int nopartavat,ncolavat;
805
806  /**
807   *
808   */
809  G4double r1_in[VARSIZE],r1_first_avat[VARSIZE];
810
811  /**
812   *
813   */
814  G4double epsd[VAEPSSIZE],eps2[VAEPSSIZE],eps4[VAEPSSIZE],eps6[VAEPSSIZE],epsf[VAEPSSIZE];
815
816  /**
817   *
818   */
819  G4int nb_avat;
820
821  /**
822   *
823   */
824  G4double timeavat[VAAVM],l1avat[VAAVM],l2avat[VAAVM],jpartl1[VAAVM],jpartl2[VAAVM];
825
826  /**
827   *
828   */
829  G4double del1avat[VAAVM],del2avat[VAAVM],energyavat[VAAVM];
830
831  /**
832   *
833   */
834  G4double bloc_paul[VAAVM],bloc_cdpp[VAAVM],go_out[VAAVM];
835};
836
837#define VARNTPSIZE 301
838class G4VarNtp {
839public:
840  G4VarNtp() {
841    clear();
842  };
843
844  ~G4VarNtp() {};
845
846  /**
847   * Clear and initialize all variables and arrays.
848   */
849  void clear() {
850    particleIndex = 0;
851    projType = 0;
852    projEnergy = 0.0;
853    targetA = 0;
854    targetZ = 0;
855    masp = 0.0; mzsp = 0.0; exsp = 0.0; mrem = 0.0;
856    // To be deleted?
857    spectatorA = 0;
858    spectatorZ = 0;
859    spectatorEx = 0.0;
860    spectatorM = 0.0;
861    spectatorT = 0.0;
862    spectatorP1 = 0.0;
863    spectatorP2 = 0.0;
864    spectatorP3 = 0.0;
865    massini = 0;
866    mzini = 0;
867    exini = 0;
868    pcorem = 0;
869    mcorem = 0;
870    pxrem = 0;
871    pyrem = 0;
872    pzrem = 0;
873    erecrem = 0;
874    mulncasc = 0;
875    mulnevap = 0;
876    mulntot = 0;
877    bimpact = 0.0;
878    jremn = 0;
879    kfis = 0;
880    estfis = 0;
881    izfis = 0;
882    iafis = 0;
883    ntrack = 0;
884    needsFermiBreakup = false;
885    for(G4int i = 0; i < VARNTPSIZE; i++) {
886      itypcasc[i] = 0;
887      avv[i] = 0;
888      zvv[i] = 0;
889      enerj[i] = 0.0;
890      plab[i] = 0.0;
891      tetlab[i] = 0.0;
892      philab[i] = 0.0;
893      full[i] = false;
894    }
895  }
896
897  /**
898   * Add a particle to the INCL/ABLA final output.
899   */
900  void addParticle(G4double A, G4double Z, G4double E, G4double P, G4double theta, G4double phi) {
901    if(full[particleIndex]) {
902      G4cout <<"G4VarNtp: Error. Index i = " << particleIndex << " is already occupied by particle:" << G4endl;
903      G4cout <<"A = " << avv[particleIndex] << " Z = " << zvv[particleIndex] << G4endl;
904      G4cout <<"Tried to replace it with:" << G4endl;
905      G4cout <<"A = " << Z << " Z = " << Z << G4endl;
906    } else {
907      avv[particleIndex] = (int) A;
908      zvv[particleIndex] = (int) Z;
909      enerj[particleIndex] = E;
910      plab[particleIndex] = P;
911      tetlab[particleIndex] = theta;
912      philab[particleIndex] = phi;
913      full[particleIndex] = true;
914      ntrack = particleIndex + 1;
915      particleIndex++;
916    }
917  }
918
919  /**
920   * Baryon number conservation check.
921   */
922  G4int getTotalBaryonNumber() {
923    G4int baryonNumber = 0;
924    for(G4int i = 0; i < ntrack; i++) {
925      if(avv[i] > 0) {
926        baryonNumber += avv[i];
927      }
928    }
929    return baryonNumber;
930  }
931
932  /**
933   * Return total energy.
934   */
935  G4double getTotalEnergy() {
936    G4double energy = 0.0;
937    for(G4int i = 0; i < ntrack; i++) {
938      energy += std::sqrt(std::pow(plab[i], 2) + std::pow(getMass(i), 2)); // E^2 = p^2 + m^2
939    }
940
941    return energy;
942  }
943
944  /**
945   * Return total three momentum.
946   */
947  G4double getTotalThreeMomentum() {
948    G4double momentum = 0;
949    for(G4int i = 0; i < ntrack; i++) {
950      momentum += plab[i];
951    }
952    return momentum;
953  }
954
955  G4double getMomentumSum() {
956    G4double momentum = 0;
957    for(G4int i = 0; i < ntrack; i++) {
958      momentum += plab[i];
959    }
960    return momentum;
961  }
962
963  G4double getMass(G4int particle) {
964    const G4double protonMass = 938.272;
965    const G4double neutronMass = 939.565;
966    const G4double pionMass = 139.57;
967
968    G4double mass = 0.0;
969    if(avv[particle] ==  1 && zvv[particle] ==  1) mass = protonMass;
970    if(avv[particle] ==  1 && zvv[particle] ==  0) mass = neutronMass;
971    if(avv[particle] == -1)                        mass = pionMass;
972    if(avv[particle] > 1)
973      mass = avv[particle] * protonMass + zvv[particle] * neutronMass;
974    return mass;
975  }
976
977  /**
978   * Dump debugging output.
979   */
980  void dump() {
981    G4int nProton = 0, nNeutron = 0;
982    G4int nPiPlus = 0, nPiZero = 0, nPiMinus = 0;
983    G4int nH2 = 0, nHe3 = 0, nAlpha = 0;
984    G4int nFragments = 0;
985    G4int nParticles = 0;
986    G4cout <<"Particles produced in the event (" << ntrack << "):" << G4endl;
987    G4cout <<"A \t Z \t Ekin \t Ptot \t Theta \t Phi" << G4endl;
988    for(G4int i = 0; i < ntrack; i++) {
989      nParticles++;
990      if(avv[i] ==  1 && zvv[i] ==  1) nProton++;  // Count multiplicities
991      if(avv[i] ==  1 && zvv[i] ==  0) nNeutron++;
992      if(avv[i] == -1 && zvv[i] ==  1) nPiPlus++;
993      if(avv[i] == -1 && zvv[i] ==  0) nPiZero++;
994      if(avv[i] == -1 && zvv[i] == -1) nPiMinus++;
995      if(avv[i] ==  2 && zvv[i] ==  1) nH2++;
996      if(avv[i] ==  3 && zvv[i] ==  2) nHe3++;
997      if(avv[i] ==  4 && zvv[i] ==  2) nAlpha++;
998      if(                zvv[i] >   2) nFragments++;
999
1000      G4cout << i << " \t " << avv[i] << " \t " << zvv[i] << " \t " << enerj[i] << " \t " 
1001             << plab[i] << " \t " << tetlab[i] << " \t " << philab[i] << G4endl;
1002    }
1003
1004    G4cout <<"Summary of event: " << G4endl;
1005    G4cout <<"Projectile type: " << projType <<" Energy: " << projEnergy << G4endl;
1006    G4cout <<"Target A = " << targetA << " Z = " << targetZ << G4endl;
1007    G4cout <<"Remnant from cascade: " << G4endl;
1008    G4cout <<"A = " << massini << " Z = " << mzini << " excitation E = " << exini << G4endl;
1009    G4cout <<"Particle multiplicities:" << G4endl;
1010    G4cout <<"Protons: " << nProton << " Neutrons:  " << nNeutron << G4endl;
1011    G4cout <<"pi+: " << nPiPlus << " pi0: " << nPiZero << " pi-: " << nPiMinus << G4endl;
1012    G4cout <<"H2: " << nH2 << " He3: " << nHe3 << " Alpha: " << nAlpha << G4endl;
1013    G4cout <<"Nucleus fragments = " << nFragments << G4endl;
1014    G4cout <<"Conservation laws:" << G4endl;
1015    G4cout <<"Baryon number = " <<  getTotalBaryonNumber() << G4endl;
1016    G4cout <<"Number of particles = " << nParticles << G4endl;
1017  }
1018
1019  /**
1020   * Projectile type.
1021   */
1022  G4int projType;
1023
1024  /**
1025   * Projectile energy.
1026   */
1027  G4double projEnergy;
1028
1029  /**
1030   * Target mass number.
1031   */
1032  G4int targetA;
1033
1034  /**
1035   * Target charge number.
1036   */
1037  G4int targetZ;
1038
1039  /**
1040   * Projectile spectator A, Z, Eex;
1041   */
1042  G4double masp, mzsp, exsp, mrem;
1043
1044  /**
1045   * Spectator nucleus mass number for light ion projectile support.
1046   */
1047  G4int spectatorA;
1048
1049  /**
1050   * Spectator nucleus charge number for light ion projectile support.
1051   */
1052  G4int spectatorZ;
1053
1054  /**
1055   * Spectator nucleus excitation energy for light ion projectile support.
1056   */
1057  G4double spectatorEx;
1058
1059  /**
1060   * Spectator nucleus mass.
1061   */
1062  G4double spectatorM;
1063
1064  /**
1065   * Spectator nucleus kinetic energy.
1066   */
1067  G4double spectatorT;
1068
1069  /**
1070   * Spectator nucleus momentum x-component.
1071   */
1072  G4double spectatorP1;
1073
1074  /**
1075   * Spectator nucleus momentum y-component.
1076   */
1077  G4double spectatorP2;
1078
1079  /**
1080   * Spectator nucleus momentum z-component.
1081   */
1082  G4double spectatorP3;
1083
1084  /**
1085   * A of the remnant.
1086   */
1087  G4double massini;
1088
1089  /**
1090   * Z of the remnant.
1091   */
1092  G4double mzini;
1093
1094  /**
1095   * Excitation energy.
1096   */
1097  G4double exini;
1098
1099  G4double pcorem, mcorem, pxrem, pyrem, pzrem, erecrem;
1100
1101  /**
1102   * Cascade n multip.
1103   */
1104  G4int mulncasc;
1105
1106  /**
1107   * Evaporation n multip.
1108   */
1109  G4int mulnevap;
1110
1111  /**
1112   * Total n multip.
1113   */
1114  G4int mulntot;
1115
1116  /**
1117   * Impact parameter.
1118   */
1119  G4double bimpact;
1120
1121  /**
1122   * Remnant Intrinsic Spin.
1123   */
1124  G4int jremn;
1125
1126  /**
1127   * Fission 1/0=Y/N.
1128   */
1129  G4int kfis;
1130
1131  /**
1132   * Excit energy at fis.
1133   */
1134  G4double estfis;
1135
1136  /**
1137   * Z of fiss nucleus.
1138   */
1139  G4int izfis;
1140
1141  /**
1142   * A of fiss nucleus.
1143   */
1144  G4int iafis;
1145
1146  /**
1147   * Number of particles.
1148   */
1149  G4int ntrack;
1150
1151  /**
1152   * The state of the index:
1153   * true = reserved
1154   * false = free
1155   */
1156  G4bool full[VARNTPSIZE];
1157
1158  /**
1159   * Does this nucleus require Fermi break-up treatment? Only
1160   * applicable when used together with Geant4.
1161   * true = do fermi break-up (and skip ABLA part)
1162   * false = use ABLA
1163   */
1164  G4bool needsFermiBreakup;
1165
1166  /**
1167   * emitted in cascade (0) or evaporation (1).
1168   */
1169  G4int itypcasc[VARNTPSIZE];
1170
1171 
1172  /**
1173   * A (-1 for pions).
1174   */
1175  G4int avv[VARNTPSIZE];
1176
1177  /**
1178   * Z
1179   */
1180  G4int zvv[VARNTPSIZE];
1181
1182  /**
1183   * Kinetic energy.
1184   */
1185  G4double enerj[VARNTPSIZE];
1186
1187  /**
1188   * Momentum.
1189   */
1190  G4double plab[VARNTPSIZE];
1191
1192  /**
1193   * Theta angle.
1194   */
1195  G4double tetlab[VARNTPSIZE];
1196
1197  /**
1198   * Phi angle.
1199   */
1200  G4double philab[VARNTPSIZE];
1201
1202private:
1203  G4int particleIndex;
1204};
1205
1206/**
1207 * Pauli blocking.
1208 */
1209class G4Paul {
1210public:
1211  G4Paul()
1212    :ct0(0.0), ct1(0.0), ct2(0.0), ct3(0.0), ct4(0.0), ct5(0.0), ct6(0.0),
1213     pr(0.0), pr2(0.0), xrr(0.0), xrr2(0.0),
1214     cp0(0.0), cp1(0.0), cp2(0.0), cp3(0.0), cp4(0.0), cp5(0.0), cp6(0.0)
1215  {};
1216
1217  ~G4Paul() {};
1218 
1219  /**
1220   *
1221   */
1222  G4double ct0,ct1,ct2,ct3,ct4,ct5,ct6,pr,pr2,xrr,xrr2;
1223
1224  /**
1225   *
1226   */
1227  G4double cp0,cp1,cp2,cp3,cp4,cp5,cp6;
1228};
1229
1230
1231#endif
Note: See TracBrowser for help on using the repository browser.