Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/lowenergy/src/G4Generator2BN.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

Diff markup

Differences between /processes/electromagnetic/lowenergy/src/G4Generator2BN.cc (Version 11.3.0) and /processes/electromagnetic/lowenergy/src/G4Generator2BN.cc (Version 10.3.p1)


  1 //                                                  1 //
  2 // *******************************************      2 // ********************************************************************
  3 // * License and Disclaimer                         3 // * License and Disclaimer                                           *
  4 // *                                                4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of th      5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided      6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License      7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/      8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.           9 // * include a list of copyright holders.                             *
 10 // *                                               10 // *                                                                  *
 11 // * Neither the authors of this software syst     11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing fin     12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warran     13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assum     14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file      15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitatio     16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                               17 // *                                                                  *
 18 // * This  code  implementation is the result      18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboratio     19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distri     20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  ag     21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publicati     22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Sof     23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // *******************************************     24 // ********************************************************************
 25 //                                                 25 //
                                                   >>  26 // $Id: G4Generator2BN.cc 95950 2016-03-03 10:42:48Z gcosmo $
 26 //                                                 27 //
 27 // -------------------------------------------     28 // -------------------------------------------------------------------
 28 //                                                 29 //
 29 // GEANT4 Class file                               30 // GEANT4 Class file
 30 //                                                 31 //
 31 //                                                 32 //
 32 // File name:     G4Generator2BN                   33 // File name:     G4Generator2BN
 33 //                                                 34 //
 34 // Author:        Andreia Trindade (andreia@li     35 // Author:        Andreia Trindade (andreia@lip.pt)
 35 //                Pedro Rodrigues  (psilva@lip     36 //                Pedro Rodrigues  (psilva@lip.pt)
 36 //                Luis Peralta     (luis@lip.p     37 //                Luis Peralta     (luis@lip.pt)
 37 //                                                 38 //
 38 // Creation date: 21 June 2003                     39 // Creation date: 21 June 2003
 39 //                                                 40 //
 40 // Modifications:                                  41 // Modifications:
 41 // 21 Jun 2003                                     42 // 21 Jun 2003                                 First implementation acording with new design
 42 // 05 Nov 2003  MGP                                43 // 05 Nov 2003  MGP                            Fixed compilation warning
 43 // 14 Mar 2004                                     44 // 14 Mar 2004                                 Code optimization
 44 //                                                 45 //
 45 // Class Description:                              46 // Class Description:
 46 //                                                 47 //
 47 // Concrete base class for Bremsstrahlung Angu     48 // Concrete base class for Bremsstrahlung Angular Distribution Generation
 48 // 2BN Distribution                                49 // 2BN Distribution
 49 //                                                 50 //
 50 // Class Description: End                          51 // Class Description: End
 51 //                                                 52 //
 52 // -------------------------------------------     53 // -------------------------------------------------------------------
 53 //                                                 54 //
 54                                                    55 
 55 #include "G4Generator2BN.hh"                       56 #include "G4Generator2BN.hh"
 56 #include "Randomize.hh"                            57 #include "Randomize.hh"
 57 #include "G4PhysicalConstants.hh"                  58 #include "G4PhysicalConstants.hh"
 58 #include "G4SystemOfUnits.hh"                      59 #include "G4SystemOfUnits.hh"
 59 #include "G4Exp.hh"                                60 #include "G4Exp.hh"
 60                                                    61 
 61 //                                                 62 //
 62                                                    63 
 63 G4double G4Generator2BN::Atab[320] =               64 G4double G4Generator2BN::Atab[320] =
 64          { 0.0264767, 0.0260006, 0.0257112, 0.     65          { 0.0264767, 0.0260006, 0.0257112, 0.0252475, 0.024792, 0.0243443, 0.0240726, 0.0236367,
 65            0.023209,  0.0227892, 0.0225362, 0.     66            0.023209,  0.0227892, 0.0225362, 0.0221284, 0.0217282,0.0214894, 0.0211005, 0.0207189,
 66            0.0204935, 0.0201227, 0.0197588, 0.     67            0.0204935, 0.0201227, 0.0197588, 0.019546,  0.0191923,0.0188454, 0.0186445, 0.0183072,
 67            0.0179763, 0.0177866, 0.0174649, 0.     68            0.0179763, 0.0177866, 0.0174649, 0.0172828, 0.0169702,0.0166634, 0.0164915, 0.0161933,
 68            0.0160283, 0.0157384, 0.0155801, 0.     69            0.0160283, 0.0157384, 0.0155801, 0.0152981, 0.0151463,0.0148721, 0.0147263, 0.0144598,
 69            0.01432,   0.0140607, 0.0139267, 0.     70            0.01432,   0.0140607, 0.0139267, 0.0136744, 0.0135459,0.0133005, 0.0131773, 0.0129385,
 70            0.0128205, 0.0125881, 0.012475,  0.     71            0.0128205, 0.0125881, 0.012475,  0.0122488, 0.0121406,0.0119204, 0.0118167, 0.0117158,
 71            0.0115032, 0.0114067, 0.0111995, 0.     72            0.0115032, 0.0114067, 0.0111995, 0.0111072, 0.0110175,0.0108173, 0.0107316, 0.0105365,
 72            0.0104547, 0.0102646, 0.0101865, 0.     73            0.0104547, 0.0102646, 0.0101865, 0.010111,  0.00992684,0.0098548,0.00967532,0.00960671,
 73            0.00943171,0.00936643,0.00930328,0.     74            0.00943171,0.00936643,0.00930328,0.0091337, 0.00907372,0.00890831,0.00885141,0.00869003,
 74            0.00863611,0.00858428,0.00842757,0.     75            0.00863611,0.00858428,0.00842757,0.00837854,0.0082256,0.00817931,0.00803,0.00798639,
 75            0.00784058,0.00779958,0.00776046,0.     76            0.00784058,0.00779958,0.00776046,0.00761866,0.00758201,0.00744346,0.00740928,0.00727384,
 76            0.00724201,0.00710969,0.00708004,0.     77            0.00724201,0.00710969,0.00708004,0.00695074,0.00692333,0.00679688,0.00677166,0.00664801,
 77            0.00662484,0.00650396,0.00648286,0.     78            0.00662484,0.00650396,0.00648286,0.00636458,0.00634545,0.00622977,0.00621258,0.00609936,
 78            0.00608412,0.00597331,0.00595991,0.     79            0.00608412,0.00597331,0.00595991,0.00585143,0.00583988,0.0057337,0.0057239,0.00561991,
 79            0.0056119, 0.00551005,0.00550377,0.     80            0.0056119, 0.00551005,0.00550377,0.00540399,0.00539938,0.00530162,0.00529872,0.00520292,
 80            0.0051091, 0.00510777,0.00501582,0.     81            0.0051091, 0.00510777,0.00501582,0.00501608,0.00492594,0.00492781,0.00483942,0.0048429,
 81            0.00475622,0.00476127,0.00467625,0.     82            0.00475622,0.00476127,0.00467625,0.00468287,0.00459947,0.00451785,0.00452581,0.00444573,
 82            0.00445522,0.00437664,0.00438768,0.     83            0.00445522,0.00437664,0.00438768,0.00431057,0.00432316,0.00424745,0.0042616,0.00418726,
 83            0.004203,  0.00413,   0.00405869,0.     84            0.004203,  0.00413,   0.00405869,0.00407563,0.00400561,0.00402414,0.00395536,0.00397553,
 84            0.00390795,0.00392975,0.00386339,0.     85            0.00390795,0.00392975,0.00386339,0.00379862,0.00382167,0.00375805,0.00378276,0.00372031,
 85            0.00374678,0.00368538,0.00371363,0.     86            0.00374678,0.00368538,0.00371363,0.00365335,0.00359463,0.0036242, 0.00356653,0.003598,
 86            0.00354139,0.00357481,0.00351921,0.     87            0.00354139,0.00357481,0.00351921,0.00355464,0.00350005,0.0034471,0.00348403,0.00343208,
 87            0.0034712, 0.00342026,0.00346165,0.     88            0.0034712, 0.00342026,0.00346165,0.00341172,0.00345548,0.00340657,0.00335944,0.00340491,
 88            0.00335885,0.00340692,0.00336191,0.     89            0.00335885,0.00340692,0.00336191,0.00341273,0.00336879,0.00342249,0.00337962,0.00333889,
 89            0.00339463,0.00335506,0.00341401,0.     90            0.00339463,0.00335506,0.00341401,0.00337558,0.00343797,0.00340067,0.00336584,0.00343059,
 90            0.0033969, 0.00346557,0.00343302,0.     91            0.0033969, 0.00346557,0.00343302,0.00350594,0.00347448,0.00344563,0.00352163,0.00349383,
 91            0.00357485,0.00354807,0.00352395,0.     92            0.00357485,0.00354807,0.00352395,0.00360885,0.00358571,0.00367661,0.00365446,0.00375194,
 92            0.00373078,0.00371234,0.00381532,0.     93            0.00373078,0.00371234,0.00381532,0.00379787,0.00390882,0.00389241,0.00387881,0.00399675,
 93            0.00398425,0.00411183,0.00410042,0.     94            0.00398425,0.00411183,0.00410042,0.00409197,0.00422843,0.00422123,0.00436974,0.0043637,
 94            0.00436082,0.00452075,0.00451934,0.     95            0.00436082,0.00452075,0.00451934,0.00452125,0.00469406,0.00469756,0.00488741,0.00489221,
 95            0.00490102,0.00510782,0.00511801,0.     96            0.00490102,0.00510782,0.00511801,0.00513271,0.0053589, 0.00537524,0.00562643,0.00564452,
 96            0.0056677, 0.00594482,0.00596999,0.     97            0.0056677, 0.00594482,0.00596999,0.0059999, 0.00630758,0.00634014,0.00637849,0.00672136,
 97            0.00676236,0.00680914,0.00719407,0.     98            0.00676236,0.00680914,0.00719407,0.0072439, 0.00730063,0.0077349, 0.00779494,0.00786293,
 98            0.00835577,0.0084276, 0.00850759,0.     99            0.00835577,0.0084276, 0.00850759,0.00907162,0.00915592,0.00924925,0.00935226,0.00999779,
 99            0.0101059, 0.0102249, 0.0109763, 0.    100            0.0101059, 0.0102249, 0.0109763, 0.0111003, 0.0112367, 0.0113862, 0.0122637, 0.0124196,
100            0.0125898, 0.0136311, 0.0138081, 0.    101            0.0125898, 0.0136311, 0.0138081, 0.0140011, 0.0142112, 0.0154536, 0.0156723, 0.0159099,
101            0.016168,  0.0176664, 0.0179339, 0.    102            0.016168,  0.0176664, 0.0179339, 0.0182246, 0.0185396, 0.020381,  0.0207026, 0.0210558,
102            0.0214374, 0.0237377, 0.0241275, 0.    103            0.0214374, 0.0237377, 0.0241275, 0.0245528, 0.0250106, 0.0255038, 0.0284158, 0.0289213,
103            0.0294621, 0.0300526, 0.0338619, 0.    104            0.0294621, 0.0300526, 0.0338619, 0.0344537, 0.0351108, 0.0358099, 0.036554,  0.0416399
104          };                                       105          };
105                                                   106 
106 G4double G4Generator2BN::ctab[320] =              107 G4double G4Generator2BN::ctab[320] =
107     { 0.482253, 0.482253, 0.489021, 0.489021,     108     { 0.482253, 0.482253, 0.489021, 0.489021, 0.489021, 0.489021, 0.495933,
108       0.495933, 0.495933, 0.495933, 0.502993,     109       0.495933, 0.495933, 0.495933, 0.502993, 0.502993, 0.502993, 0.510204,
109       0.510204, 0.510204, 0.517572, 0.517572,     110       0.510204, 0.510204, 0.517572, 0.517572, 0.517572, 0.5251,   0.5251,
110       0.5251,   0.532793, 0.532793, 0.532793,     111       0.5251,   0.532793, 0.532793, 0.532793, 0.540657, 0.540657, 0.548697,
111       0.548697, 0.548697, 0.556917, 0.556917,     112       0.548697, 0.548697, 0.556917, 0.556917, 0.565323, 0.565323, 0.573921,
112       0.573921, 0.582717, 0.582717, 0.591716,     113       0.573921, 0.582717, 0.582717, 0.591716, 0.591716, 0.600925, 0.600925,
113       0.610352, 0.610352, 0.620001, 0.620001,     114       0.610352, 0.610352, 0.620001, 0.620001, 0.629882, 0.629882, 0.64,
114       0.64,     0.650364, 0.650364, 0.660982,     115       0.64,     0.650364, 0.650364, 0.660982, 0.660982, 0.671862, 0.683013,
115       0.683013, 0.694444, 0.694444, 0.706165,     116       0.683013, 0.694444, 0.694444, 0.706165, 0.718184, 0.718184, 0.730514,
116       0.730514, 0.743163, 0.743163, 0.756144,     117       0.730514, 0.743163, 0.743163, 0.756144, 0.769468, 0.769468, 0.783147,
117       0.783147, 0.797194, 0.797194, 0.811622,     118       0.783147, 0.797194, 0.797194, 0.811622, 0.826446, 0.826446, 0.84168,
118       0.84168,  0.857339, 0.857339, 0.873439,     119       0.84168,  0.857339, 0.857339, 0.873439, 0.889996, 0.889996, 0.907029,
119       0.907029, 0.924556, 0.924556, 0.942596,     120       0.907029, 0.924556, 0.924556, 0.942596, 0.942596, 0.961169, 0.980296,
120       0.980296, 1,        1,        1.0203,       121       0.980296, 1,        1,        1.0203,   1.0203,   1.04123,  1.04123,
121       1.06281,  1.06281,  1.08507,  1.08507,      122       1.06281,  1.06281,  1.08507,  1.08507,  1.10803,  1.10803,  1.13173,
122       1.13173,  1.1562,   1.1562,   1.18147,      123       1.13173,  1.1562,   1.1562,   1.18147,  1.18147,  1.20758,  1.20758,
123       1.23457,  1.23457,  1.26247,  1.26247,      124       1.23457,  1.23457,  1.26247,  1.26247,  1.29132,  1.29132,  1.32118,
124       1.32118,  1.35208,  1.35208,  1.38408,      125       1.32118,  1.35208,  1.35208,  1.38408,  1.38408,  1.41723,  1.41723,
125       1.45159,  1.45159,  1.45159,  1.48721,      126       1.45159,  1.45159,  1.45159,  1.48721,  1.48721,  1.52416,  1.52416,
126       1.5625,   1.5625,   1.60231,  1.60231,      127       1.5625,   1.5625,   1.60231,  1.60231,  1.64366,  1.64366,  1.68663,
127       1.68663,  1.68663,  1.7313,   1.7313,       128       1.68663,  1.68663,  1.7313,   1.7313,   1.77778,  1.77778,  1.82615,
128       1.82615,  1.87652,  1.87652,  1.92901,      129       1.82615,  1.87652,  1.87652,  1.92901,  1.92901,  1.98373,  1.98373,
129       1.98373,  2.04082,  2.04082,  2.1004,       130       1.98373,  2.04082,  2.04082,  2.1004,   2.1004,   2.16263,  2.16263,
130       2.22767,  2.22767,  2.22767,  2.29568,      131       2.22767,  2.22767,  2.22767,  2.29568,  2.29568,  2.36686,  2.36686,
131       2.44141,  2.44141,  2.51953,  2.51953,      132       2.44141,  2.44141,  2.51953,  2.51953,  2.51953,  2.60146,  2.60146,
132       2.68745,  2.68745,  2.77778,  2.77778,      133       2.68745,  2.68745,  2.77778,  2.77778,  2.87274,  2.87274,  2.87274,
133       2.97265,  2.97265,  3.07787,  3.07787,      134       2.97265,  2.97265,  3.07787,  3.07787,  3.18878,  3.18878,  3.30579,
134       3.30579,  3.30579,  3.42936,  3.42936,      135       3.30579,  3.30579,  3.42936,  3.42936,  3.55999,  3.55999,  3.69822,
135       3.69822,  3.84468,  3.84468,  3.84468,      136       3.69822,  3.84468,  3.84468,  3.84468,  4,        4,        4.16493,
136       4.16493,  4.34028,  4.34028,  4.34028,      137       4.16493,  4.34028,  4.34028,  4.34028,  4.52694,  4.52694,  4.7259,
137       4.7259,   4.93827,  4.93827,  4.93827,      138       4.7259,   4.93827,  4.93827,  4.93827,  5.16529,  5.16529,  5.40833,
138       5.40833,  5.40833,  5.66893,  5.66893,      139       5.40833,  5.40833,  5.66893,  5.66893,  5.94884,  5.94884,  6.25,
139       6.25,     6.25,     6.57462,  6.57462,      140       6.25,     6.25,     6.57462,  6.57462,  6.92521,  6.92521,  6.92521,
140       7.3046,   7.3046,   7.71605,  7.71605,      141       7.3046,   7.3046,   7.71605,  7.71605,  7.71605,  8.16327,  8.16327,
141       8.65052,  8.65052,  8.65052,  9.18274,      142       8.65052,  8.65052,  8.65052,  9.18274,  9.18274,  9.18274,  9.76562,
142       9.76562,  10.4058,  10.4058,  10.4058,      143       9.76562,  10.4058,  10.4058,  10.4058,  11.1111,  11.1111,  11.1111,
143       11.8906,  11.8906,  12.7551,  12.7551,      144       11.8906,  11.8906,  12.7551,  12.7551,  12.7551,  13.7174,  13.7174,
144       13.7174,  14.7929,  14.7929,  14.7929,      145       13.7174,  14.7929,  14.7929,  14.7929,  16,       16,       16,
145       17.3611,  17.3611,  17.3611,  18.9036,      146       17.3611,  17.3611,  17.3611,  18.9036,  18.9036,  18.9036,  20.6612,
146       20.6612,  20.6612,  22.6757,  22.6757,      147       20.6612,  20.6612,  22.6757,  22.6757,  22.6757,  22.6757,  25,
147       25,       25,       27.7008,  27.7008,      148       25,       25,       27.7008,  27.7008,  27.7008,  27.7008,  30.8642,
148       30.8642,  30.8642,  34.6021,  34.6021,      149       30.8642,  30.8642,  34.6021,  34.6021,  34.6021,  34.6021,  39.0625,
149       39.0625,  39.0625,  39.0625,  44.4444,      150       39.0625,  39.0625,  39.0625,  44.4444,  44.4444,  44.4444,  44.4444,
150       51.0204,  51.0204,  51.0204,  51.0204,      151       51.0204,  51.0204,  51.0204,  51.0204,  59.1716,  59.1716,  59.1716,
151       59.1716,  59.1716,  69.4444,  69.4444,      152       59.1716,  59.1716,  69.4444,  69.4444,  69.4444,  69.4444,  82.6446,
152       82.6446,  82.6446,  82.6446,  82.6446,      153       82.6446,  82.6446,  82.6446,  82.6446,  100
153     };                                            154     };
154                                                   155 
155 //....oooOO0OOooo........oooOO0OOooo........oo << 
156                                                   156 
157 G4Generator2BN::G4Generator2BN(const G4String&    157 G4Generator2BN::G4Generator2BN(const G4String&)
158  : G4VEmAngularDistribution("AngularGen2BN")      158  : G4VEmAngularDistribution("AngularGen2BN")
159 {                                                 159 {
160   b = 1.2;                                        160   b = 1.2;
161   index_min = -300;                               161   index_min = -300;
162   index_max = 319;                                162   index_max = 319;
163                                                   163 
164   // Set parameters minimum limits Ekelectron     164   // Set parameters minimum limits Ekelectron = 250 eV and kphoton = 100 eV
165   kmin = 100*eV;                                  165   kmin = 100*eV;
166   Ekmin = 250*eV;                                 166   Ekmin = 250*eV;
167   kcut = 100*eV;                                  167   kcut = 100*eV;
168                                                   168 
169   // Increment Theta value for surface interpo    169   // Increment Theta value for surface interpolation
170   dtheta = 0.1*rad;                               170   dtheta = 0.1*rad;
171                                                   171 
172   nwarn = 0;                                      172   nwarn = 0;
173                                                   173 
174   // Construct Majorant Surface to 2BN cross-s    174   // Construct Majorant Surface to 2BN cross-section
175   // (we dont need this. Pre-calculated values    175   // (we dont need this. Pre-calculated values are used instead due to performance issues
176   // ConstructMajorantSurface();                  176   // ConstructMajorantSurface();
177 }                                                 177 }
178                                                   178 
179 //....oooOO0OOooo........oooOO0OOooo........oo << 179 //
                                                   >> 180 
                                                   >> 181 G4Generator2BN::~G4Generator2BN()
                                                   >> 182 {}
180                                                   183 
181 G4ThreeVector& G4Generator2BN::SampleDirection    184 G4ThreeVector& G4Generator2BN::SampleDirection(const G4DynamicParticle* dp,
182                  G4double out_energy,             185                  G4double out_energy,
183                  G4int,                           186                  G4int,
184                  const G4Material*)               187                  const G4Material*)
185 {                                                 188 {
186   G4double Ek  = dp->GetKineticEnergy();          189   G4double Ek  = dp->GetKineticEnergy();
187   G4double Eel = dp->GetTotalEnergy();            190   G4double Eel = dp->GetTotalEnergy();
188   if(Eel > 2*MeV) {                               191   if(Eel > 2*MeV) {
189     return fGenerator2BS.SampleDirection(dp, o    192     return fGenerator2BS.SampleDirection(dp, out_energy, 0);
190   }                                               193   }
191                                                   194 
192   G4double k   = Eel - out_energy;                195   G4double k   = Eel - out_energy;
193                                                   196 
194   G4double t;                                     197   G4double t;
195   G4double cte2;                                  198   G4double cte2;
196   G4double y, u;                                  199   G4double y, u;
197   G4double ds;                                    200   G4double ds;
198   G4double dmax;                                  201   G4double dmax;
199                                                   202 
                                                   >> 203   G4int trials = 0;
                                                   >> 204 
200   // find table index                             205   // find table index
201   G4int index = G4int(std::log10(Ek/MeV)*100)     206   G4int index = G4int(std::log10(Ek/MeV)*100) - index_min;
202   if(index > index_max) { index = index_max; }    207   if(index > index_max) { index = index_max; }
203   else if(index < 0) { index = 0; }               208   else if(index < 0) { index = 0; }
204                                                   209 
                                                   >> 210   //kmax = Ek;
                                                   >> 211   //kmin2 = kcut;
                                                   >> 212 
205   G4double c = ctab[index];                       213   G4double c = ctab[index];
206   G4double A = Atab[index];                       214   G4double A = Atab[index];
207   if(index < index_max) { A = std::max(A, Atab    215   if(index < index_max) { A = std::max(A, Atab[index+1]); }
208                                                   216 
209   do{                                             217   do{
210     // generate k accordimg to std::pow(k,-b)     218     // generate k accordimg to std::pow(k,-b)
                                                   >> 219     trials++;
211                                                   220 
212     // normalization constant                     221     // normalization constant
213     //  cte1 = (1-b)/(std::pow(kmax,(1-b))-std    222     //  cte1 = (1-b)/(std::pow(kmax,(1-b))-std::pow(kmin2,(1-b)));
214     //  y = G4UniformRand();                      223     //  y = G4UniformRand();
215     //  k = std::pow(((1-b)/cte1*y+std::pow(km    224     //  k = std::pow(((1-b)/cte1*y+std::pow(kmin2,(1-b))),(1/(1-b)));
216                                                   225 
217     // generate theta accordimg to theta/(1+c*    226     // generate theta accordimg to theta/(1+c*std::pow(theta,2)
218     // Normalization constant                     227     // Normalization constant
219     cte2 = 2*c/std::log(1+c*pi2);                 228     cte2 = 2*c/std::log(1+c*pi2);
220                                                   229 
221     y = G4UniformRand();                          230     y = G4UniformRand();
222     t = std::sqrt((G4Exp(2*c*y/cte2)-1)/c);       231     t = std::sqrt((G4Exp(2*c*y/cte2)-1)/c);
223     u = G4UniformRand();                          232     u = G4UniformRand();
224                                                   233 
225     // point acceptance algorithm                 234     // point acceptance algorithm
226     //fk = std::pow(k,-b);                        235     //fk = std::pow(k,-b);
227     //ft = t/(1+c*t*t);                           236     //ft = t/(1+c*t*t);
228     dmax = A*std::pow(k,-b)*t/(1+c*t*t);          237     dmax = A*std::pow(k,-b)*t/(1+c*t*t);
229     ds = Calculatedsdkdt(k,t,Eel);                238     ds = Calculatedsdkdt(k,t,Eel);
230                                                   239 
231     // violation point                            240     // violation point
232     if(ds > dmax && nwarn >= 20) {                241     if(ds > dmax && nwarn >= 20) {
233       ++nwarn;                                    242       ++nwarn;
234       G4cout << "### WARNING in G4Generator2BN    243       G4cout << "### WARNING in G4Generator2BN: Ekin(MeV)= " << Ek/MeV
235        << "  D(Ekin,k)/Dmax-1= " << (ds/dmax -    244        << "  D(Ekin,k)/Dmax-1= " << (ds/dmax - 1)
236        << "  results are not reliable!"           245        << "  results are not reliable!"
237        << G4endl;                                 246        << G4endl;
238       if(20 == nwarn) {                           247       if(20 == nwarn) {
239   G4cout << "   WARNING in G4Generator2BN is c    248   G4cout << "   WARNING in G4Generator2BN is closed" << G4endl;
240       }                                           249       }
241     }                                             250     }
242                                                   251 
243   } while(u*dmax > ds || t > CLHEP::pi);          252   } while(u*dmax > ds || t > CLHEP::pi);
244                                                   253 
245   G4double sint = std::sin(t);                    254   G4double sint = std::sin(t);
246   G4double phi  = twopi*G4UniformRand();          255   G4double phi  = twopi*G4UniformRand();
247                                                   256 
248   fLocalDirection.set(sint*std::cos(phi), sint    257   fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi),std::cos(t));
249   fLocalDirection.rotateUz(dp->GetMomentumDire    258   fLocalDirection.rotateUz(dp->GetMomentumDirection());
250                                                   259 
251   return fLocalDirection;                         260   return fLocalDirection;
252 }                                                 261 }
253                                                   262 
254 //....oooOO0OOooo........oooOO0OOooo........oo << 
255                                                << 
256 G4double G4Generator2BN::CalculateFkt(G4double    263 G4double G4Generator2BN::CalculateFkt(G4double k, G4double theta, G4double A, G4double c) const
257 {                                                 264 {
258   G4double Fkt_value = 0;                         265   G4double Fkt_value = 0;
259   Fkt_value = A*std::pow(k,-b)*theta/(1+c*thet    266   Fkt_value = A*std::pow(k,-b)*theta/(1+c*theta*theta);
260   return Fkt_value;                               267   return Fkt_value;
261 }                                                 268 }
262                                                   269 
263 //....oooOO0OOooo........oooOO0OOooo........oo << 
264                                                << 
265 G4double G4Generator2BN::Calculatedsdkdt(G4dou    270 G4double G4Generator2BN::Calculatedsdkdt(G4double kout, G4double theta, G4double Eel) const
266 {                                                 271 {
                                                   >> 272 
267   G4double dsdkdt_value = 0.;                     273   G4double dsdkdt_value = 0.;
268   G4double Z = 1;                                 274   G4double Z = 1;
269   // classic radius (in cm)                       275   // classic radius (in cm)
270   G4double r0 = 2.82E-13;                         276   G4double r0 = 2.82E-13;
271   //squared classic radius (in barn)              277   //squared classic radius (in barn)
272   G4double r02 = r0*r0*1.0E+24;                   278   G4double r02 = r0*r0*1.0E+24;
273                                                   279 
274   // Photon energy cannot be greater than elec    280   // Photon energy cannot be greater than electron kinetic energy
275   if(kout > (Eel-electron_mass_c2)){              281   if(kout > (Eel-electron_mass_c2)){
276     dsdkdt_value = 0;                             282     dsdkdt_value = 0;
277     return dsdkdt_value;                          283     return dsdkdt_value;
278   }                                               284   }
279                                                   285 
280   G4double E0 = Eel/electron_mass_c2;          << 286      G4double E0 = Eel/electron_mass_c2;
281   G4double k =  kout/electron_mass_c2;         << 287      G4double k =  kout/electron_mass_c2;
282   G4double E =  E0 - k;                        << 288      G4double E =  E0 - k;
283                                                << 289 
284   if(E <= 1*MeV ){                             << 290      if(E <= 1*MeV ){                                 // Kinematic limit at 1 MeV !!
285     dsdkdt_value = 0;                          << 291        dsdkdt_value = 0;
286     return dsdkdt_value;                       << 292        return dsdkdt_value;
287   }                                            << 293      }
288                                                << 294 
289   G4double p0 = std::sqrt(E0*E0-1);            << 295 
290   G4double p = std::sqrt(E*E-1);               << 296      G4double p0 = std::sqrt(E0*E0-1);
291   G4double LL = std::log((E*E0-1+p*p0)/(E*E0-1 << 297      G4double p = std::sqrt(E*E-1);
292   G4double delta0 = E0 - p0*std::cos(theta);   << 298      G4double LL = std::log((E*E0-1+p*p0)/(E*E0-1-p*p0));
293   G4double epsilon = std::log((E+p)/(E-p));    << 299      G4double delta0 = E0 - p0*std::cos(theta);
294   G4double Z2 = Z*Z;                           << 300      G4double epsilon = std::log((E+p)/(E-p));
295   G4double sintheta2 = std::sin(theta)*std::si << 301      G4double Z2 = Z*Z;
296   G4double E02 = E0*E0;                        << 302      G4double sintheta2 = std::sin(theta)*std::sin(theta);
297   G4double E2 = E*E;                           << 303      G4double E02 = E0*E0;
298   G4double p02 = E0*E0-1;                      << 304      G4double E2 = E*E;
299   G4double k2 = k*k;                           << 305      G4double p02 = E0*E0-1;
300   G4double delta02 = delta0*delta0;            << 306      G4double k2 = k*k;
301   G4double delta04 =  delta02* delta02;        << 307      G4double delta02 = delta0*delta0;
302   G4double Q = std::sqrt(p02+k2-2*k*p0*std::co << 308      G4double delta04 =  delta02* delta02;
303   G4double Q2 = Q*Q;                           << 309      G4double Q = std::sqrt(p02+k2-2*k*p0*std::cos(theta));
304   G4double epsilonQ = std::log((Q+p)/(Q-p));   << 310      G4double Q2 = Q*Q;
305                                                << 311      G4double epsilonQ = std::log((Q+p)/(Q-p));
306                                                << 312 
307   dsdkdt_value = Z2 * (r02/(8*pi*137)) * (1/k) << 313 
308     ( (8 * (sintheta2*(2*E02+1))/(p02*delta04) << 314      dsdkdt_value = Z2 * (r02/(8*pi*137)) * (1/k) * (p/p0) *
309       ((2*(5*E02+2*E*E0+3))/(p02 * delta02)) - << 315        ( (8 * (sintheta2*(2*E02+1))/(p02*delta04)) -
310       ((2*(p02-k2))/((Q2*delta02))) +          << 316          ((2*(5*E02+2*E*E0+3))/(p02 * delta02)) -
311       ((4*E)/(p02*delta0)) +                   << 317          ((2*(p02-k2))/((Q2*delta02))) +
312       (LL/(p*p0))*(                            << 318          ((4*E)/(p02*delta0)) +
313        ((4*E0*sintheta2*(3*k-p02*E))/(p02*delt << 319          (LL/(p*p0))*(
314        ((4*E02*(E02+E2))/(p02*delta02)) +      << 320                  ((4*E0*sintheta2*(3*k-p02*E))/(p02*delta04)) +
315        ((2-2*(7*E02-3*E*E0+E2))/(p02*delta02)) << 321                  ((4*E02*(E02+E2))/(p02*delta02)) +
316        (2*k*(E02+E*E0-1))/((p02*delta0))       << 322                  ((2-2*(7*E02-3*E*E0+E2))/(p02*delta02)) +
317        ) -                                     << 323                  (2*k*(E02+E*E0-1))/((p02*delta0))
318       ((4*epsilon)/(p*delta0)) +               << 324                  ) -
319       ((epsilonQ)/(p*Q))*                      << 325          ((4*epsilon)/(p*delta0)) +
320       (4/delta02-(6*k/delta0)-(2*k*(p02-k2))/( << 326          ((epsilonQ)/(p*Q))*
321       );                                       << 327          (4/delta02-(6*k/delta0)-(2*k*(p02-k2))/(Q2*delta0))
322                                                << 328          );
323   dsdkdt_value = dsdkdt_value*std::sin(theta); << 329 
324   return dsdkdt_value;                         << 
325 }                                              << 
326                                                   330 
327 //....oooOO0OOooo........oooOO0OOooo........oo << 331 
                                                   >> 332      dsdkdt_value = dsdkdt_value*std::sin(theta);
                                                   >> 333      return dsdkdt_value;
                                                   >> 334 }
328                                                   335 
329 void G4Generator2BN::ConstructMajorantSurface(    336 void G4Generator2BN::ConstructMajorantSurface()
330 {                                                 337 {
                                                   >> 338 
331   G4double Eel;                                   339   G4double Eel;
332   G4int vmax;                                     340   G4int vmax;
333   G4double Ek;                                    341   G4double Ek;
334   G4double k, theta, k0, theta0, thetamax;        342   G4double k, theta, k0, theta0, thetamax;
335   G4double ds, df, dsmax, dk, dt;                 343   G4double ds, df, dsmax, dk, dt;
336   G4double ratmin;                                344   G4double ratmin;
337   G4double ratio = 0;                             345   G4double ratio = 0;
338   G4double Vds, Vdf;                              346   G4double Vds, Vdf;
339   G4double A, c;                                  347   G4double A, c;
340                                                   348 
341   G4cout << "**** Constructing Majorant Surfac    349   G4cout << "**** Constructing Majorant Surface for 2BN Distribution ****" << G4endl;
342                                                   350 
343   if(kcut > kmin) kmin = kcut;                    351   if(kcut > kmin) kmin = kcut;
344                                                   352 
345   G4int i = 0;                                    353   G4int i = 0;
346   // Loop on energy index                         354   // Loop on energy index
347   for(G4int index = index_min; index < index_m    355   for(G4int index = index_min; index < index_max; index++){
348                                                   356 
349   G4double fraction = index/100.;                 357   G4double fraction = index/100.;
350   Ek = std::pow(10.,fraction);                    358   Ek = std::pow(10.,fraction);
351   Eel = Ek + electron_mass_c2;                    359   Eel = Ek + electron_mass_c2;
352                                                   360 
353   // find x-section maximum at k=kmin             361   // find x-section maximum at k=kmin
354   dsmax = 0.;                                     362   dsmax = 0.;
355   thetamax = 0.;                                  363   thetamax = 0.;
356                                                   364 
357   for(theta = 0.; theta < pi; theta = theta +     365   for(theta = 0.; theta < pi; theta = theta + dtheta){
358                                                   366 
359     ds = Calculatedsdkdt(kmin, theta, Eel);       367     ds = Calculatedsdkdt(kmin, theta, Eel);
360                                                   368 
361     if(ds > dsmax){                               369     if(ds > dsmax){
362       dsmax = ds;                                 370       dsmax = ds;
363       thetamax = theta;                           371       thetamax = theta;
364     }                                             372     }
365   }                                               373   }
366                                                   374 
367   // Compute surface paremeters at kmin           375   // Compute surface paremeters at kmin
368   if(Ek < kmin || thetamax == 0){                 376   if(Ek < kmin || thetamax == 0){
369     c = 0;                                        377     c = 0;
370     A = 0;                                        378     A = 0;
371   }else{                                          379   }else{
372     c = 1/(thetamax*thetamax);                    380     c = 1/(thetamax*thetamax);
373     A = 2*std::sqrt(c)*dsmax/(std::pow(kmin,-b    381     A = 2*std::sqrt(c)*dsmax/(std::pow(kmin,-b));
374   }                                               382   }
375                                                   383 
376   // look for correction factor to normalizati    384   // look for correction factor to normalization at kmin
377   ratmin = 1.;                                    385   ratmin = 1.;
378                                                   386 
379   // Volume under surfaces                        387   // Volume under surfaces
380   Vds = 0.;                                       388   Vds = 0.;
381   Vdf = 0.;                                       389   Vdf = 0.;
382   k0 = 0.;                                        390   k0 = 0.;
383   theta0 = 0.;                                    391   theta0 = 0.;
384                                                   392 
385   vmax = G4int(100.*std::log10(Ek/kmin));         393   vmax = G4int(100.*std::log10(Ek/kmin));
386                                                   394 
387   for(G4int v = 0; v < vmax; v++){                395   for(G4int v = 0; v < vmax; v++){
388     G4double fractionLocal = (v/100.);            396     G4double fractionLocal = (v/100.);
389     k = std::pow(10.,fractionLocal)*kmin;         397     k = std::pow(10.,fractionLocal)*kmin;
390                                                   398 
391     for(theta = 0.; theta < pi; theta = theta     399     for(theta = 0.; theta < pi; theta = theta + dtheta){
392       dk = k - k0;                                400       dk = k - k0;
393       dt = theta - theta0;                        401       dt = theta - theta0;
394       ds = Calculatedsdkdt(k,theta, Eel);         402       ds = Calculatedsdkdt(k,theta, Eel);
395       Vds = Vds + ds*dk*dt;                       403       Vds = Vds + ds*dk*dt;
396       df = CalculateFkt(k,theta, A, c);           404       df = CalculateFkt(k,theta, A, c);
397       Vdf = Vdf + df*dk*dt;                       405       Vdf = Vdf + df*dk*dt;
398                                                   406 
399       if(ds != 0.){                               407       if(ds != 0.){
400   if(df != 0.) ratio = df/ds;                     408   if(df != 0.) ratio = df/ds;
401       }                                           409       }
402                                                   410 
403       if(ratio < ratmin && ratio != 0.){          411       if(ratio < ratmin && ratio != 0.){
404   ratmin = ratio;                                 412   ratmin = ratio;
405       }                                           413       }
406     }                                             414     }
407   }                                               415   }
408                                                   416 
                                                   >> 417 
409   // sampling efficiency                          418   // sampling efficiency
410   Atab[i] = A/ratmin * 1.04;                      419   Atab[i] = A/ratmin * 1.04;
411   ctab[i] = c;                                    420   ctab[i] = c;
412                                                   421 
413 //  G4cout << Ek << " " << i << " " << index <    422 //  G4cout << Ek << " " << i << " " << index << " " << Atab[i] << " " << ctab[i] << " " << G4endl;
414   i++;                                            423   i++;
415   }                                               424   }
416 }                                              << 
417                                                   425 
418 //....oooOO0OOooo........oooOO0OOooo........oo << 426 }
419                                                   427 
420 void G4Generator2BN::PrintGeneratorInformation    428 void G4Generator2BN::PrintGeneratorInformation() const
421 {                                                 429 {
422   G4cout << "\n" << G4endl;                       430   G4cout << "\n" << G4endl;
423   G4cout << "Bremsstrahlung Angular Generator     431   G4cout << "Bremsstrahlung Angular Generator is 2BN Generator from 2BN Koch & Motz distribution (Rev Mod Phys 31(4), 920 (1959))" << G4endl;
424   G4cout << "\n" << G4endl;                       432   G4cout << "\n" << G4endl;
425 }                                                 433 }
426                                                   434