Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/particle_hp/src/G4NRESP71M03.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 ]

  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 #include "G4NRESP71M03.hh"
 27 
 28 #include "G4Alpha.hh"
 29 #include "G4Geantino.hh"
 30 #include "G4IonTable.hh"
 31 #include "G4Neutron.hh"
 32 #include "G4PhysicalConstants.hh"
 33 #include "G4RotationMatrix.hh"
 34 #include "G4SystemOfUnits.hh"
 35 #include "Randomize.hh"
 36 
 37 // Energies for which angular distributions are given.
 38 const G4double G4NRESP71M03::BEN2[ndist] = {
 39   5700.,  8000.,  8640.,  8990.,  9220.,  9410.,  9830.,  10400., 10800., 11250., 11460.,
 40   11870., 12140., 12320., 12570., 12940., 13420., 13760., 14020., 14200., 14440., 14620.,
 41   14820., 15050., 15660., 15980., 16470., 16940., 17970., 18000., 19000., 20000.};
 42 // Angular distribution of alpha particles for each energy value in BEN2.
 43 const G4double G4NRESP71M03::B2[ndist][nrhos] = {
 44   {0.000,  2839.,  4027.,  4951.,  5736.,  6437.,  7074.,  7671.,  8230.,  8765.,  9273.,
 45    9767.,  10241., 10703., 11152., 11595., 12026., 12453., 12871., 13285., 13697., 14102.,
 46    14507., 14909., 15308., 15711., 16110., 16509., 16911., 17316., 17721., 18133., 18547.,
 47    18965., 19393., 19823., 20266., 20715., 21177., 21651., 22145., 22654., 23188., 23747.,
 48    24344., 24981., 25682., 26467., 27391., 28579., 31416.}, /* En=5700 keV */
 49   {0.000,  1771.,  2528.,  3129.,  3653.,  4128.,  4574.,  4998.,  5403.,  5802.,  6192.,
 50    6578.,  6961.,  7345.,  7725.,  8108.,  8494.,  8887.,  9280.,  9676.,  10078., 10483.,
 51    10891., 11306., 11721., 12142., 12563., 12987., 13411., 13841., 14272., 14705., 15142.,
 52    15585., 16034., 16490., 16958., 17438., 17935., 18453., 18994., 19572., 20187., 20857.,
 53    21595., 22424., 23373., 24482., 25820., 27539., 31416.}, /* En=8000 keV */
 54   {0.000,  3518.,  4863.,  5805.,  6540.,  7143.,  7665.,  8121.,  8535.,  8912.,  9261.,
 55    9588.,  9896.,  10191., 10470., 10744., 11004., 11259., 11507., 11749., 11988., 12223.,
 56    12456., 12685., 12915., 13144., 13373., 13603., 13832., 14064., 14300., 14539., 14784.,
 57    15032., 15286., 15547., 15821., 16100., 16399., 16710., 17043., 17401., 17790., 18221.,
 58    18711., 19283., 19980., 20897., 22251., 24554., 31416.}, /* En=8640 keV */
 59   {0.000,  2287.,  3364.,  4319.,  5277.,  6327.,  7536.,  8868.,  10043., 10964., 11696.,
 60    12305., 12839., 13317., 13757., 14168., 14558., 14931., 15293., 15645., 15993., 16339.,
 61    16681., 17024., 17372., 17724., 18086., 18453., 18833., 19226., 19638., 20065., 20511.,
 62    20976., 21463., 21965., 22481., 23005., 23527., 24048., 24560., 25069., 25569., 26068.,
 63    26568., 27080., 27614., 28183., 28820., 29612., 31416.}, /* En=8990 keV */
 64   {0.000,  1690.,  2456.,  3097.,  3697.,  4294.,  4919.,  5614.,  6433.,  7517.,  9076.,
 65    10722., 11906., 12795., 13515., 14134., 14686., 15189., 15654., 16094., 16515., 16920.,
 66    17313., 17699., 18082., 18463., 18849., 19239., 19638., 20049., 20480., 20932., 21406.,
 67    21912., 22446., 23002., 23571., 24133., 24677., 25195., 25682., 26147., 26590., 27023.,
 68    27451., 27878., 28321., 28789., 29314., 29955., 31416.}, /* En=9220 keV */
 69   {0.000,  1674.,  2434.,  3078.,  3685.,  4297.,  4951.,  5695.,  6619.,  7913.,  9597.,
 70    11014., 12069., 12918., 13640., 14278., 14853., 15378., 15865., 16323., 16757., 17171.,
 71    17574., 17966., 18353., 18733., 19116., 19502., 19892., 20294., 20712., 21146., 21601.,
 72    22085., 22591., 23122., 23665., 24212., 24743., 25248., 25735., 26194., 26637., 27067.,
 73    27492., 27919., 28355., 28820., 29336., 29973., 31416.}, /* En=9410 keV */
 74   {0.000,  2302.,  3245.,  3967.,  4577.,  5117.,  5617.,  6085.,  6534.,  6971.,  7401.,
 75    7838.,  8278.,  8739.,  9229.,  9767.,  10382., 11140., 12186., 13637., 14818., 15604.,
 76    16198., 16688., 17121., 17514., 17881., 18230., 18569., 18902., 19232., 19568., 19911.,
 77    20269., 20646., 21051., 21497., 22009., 22619., 23364., 24253., 25148., 25921., 26571.,
 78    27146., 27677., 28189., 28710., 29267., 29936., 31416.}, /* En=9830 keV */
 79   {0.000,  2195.,  2934.,  3458.,  3879.,  4244.,  4567.,  4866.,  5142.,  5406.,  5658.,
 80    5899.,  6138.,  6371.,  6600.,  6832.,  7062.,  7291.,  7527.,  7766.,  8011.,  8265.,
 81    8529.,  8809.,  9104.,  9424.,  9773.,  10159., 10602., 11124., 11762., 12572., 13643.,
 82    15060., 16801., 18312., 19352., 20140., 20794., 21378., 21925., 22449., 22971., 23505.,
 83    24067., 24674., 25355., 26153., 27137., 28443., 31416.}, /* En=10400 keV */
 84   {0.000,  1712.,  2406.,  2937.,  3380.,  3773.,  4134.,  4470.,  4787.,  5095.,  5394.,
 85    5689.,  5978.,  6270.,  6565.,  6864.,  7175.,  7495.,  7831.,  8190.,  8579.,  9003.,
 86    9484.,  10034., 10675., 11422., 12258., 13147., 14055., 14994., 15996., 17077., 18164.,
 87    19138., 19971., 20696., 21353., 21965., 22553., 23128., 23703., 24281., 24865., 25462.,
 88    26068., 26684., 27316., 27972., 28689., 29543., 31416.}, /* En=10800 keV */
 89   {0.000,  1853.,  2651.,  3289.,  3851.,  4366.,  4853.,  5324.,  5786.,  6245.,  6707.,
 90    7172.,  7646.,  8139.,  8645.,  9179.,  9735.,  10326., 10945., 11592., 12261., 12943.,
 91    13621., 14287., 14928., 15544., 16128., 16691., 17228., 17746., 18249., 18736., 19213.,
 92    19682., 20143., 20602., 21061., 21523., 21984., 22452., 22933., 23423., 23929., 24457.,
 93    25013., 25603., 26247., 26964., 27803., 28874., 31416.}, /* En=11250 keV */
 94   {0.000,  1959.,  2770.,  3397.,  3934.,  4417.,  4866.,  5294.,  5709.,  6118.,  6529.,
 95    6948.,  7382.,  7840.,  8333.,  8876.,  9484.,  10167., 10902., 11629., 12292., 12877.,
 96    13398., 13871., 14309., 14725., 15125., 15517., 15906., 16298., 16697., 17109., 17540.,
 97    17992., 18473., 18982., 19519., 20077., 20646., 21217., 21783., 22344., 22906., 23477.,
 98    24069., 24699., 25391., 26185., 27146., 28422., 31416.}, /* En=11460 keV */
 99   {0.000,  2290.,  3391.,  4447.,  5824.,  8351.,  9448.,  10077., 10553., 10952., 11306.,
100    11630., 11935., 12227., 12510., 12789., 13066., 13344., 13628., 13918., 14220., 14538.,
101    14875., 15239., 15637., 16077., 16565., 17101., 17671., 18251., 18820., 19368., 19897.,
102    20413., 20927., 21449., 21989., 22558., 23162., 23802., 24459., 25110., 25731., 26317.,
103    26873., 27410., 27944., 28493., 29091., 29810., 31416.}, /* En=11870 keV */
104   {0.000,  2119.,  3072.,  3874.,  4627.,  5378.,  6155.,  6966.,  7786.,  8568.,  9281.,
105    9921.,  10501., 11034., 11531., 12002., 12452., 12886., 13306., 13715., 14114., 14504.,
106    14888., 15266., 15641., 16013., 16386., 16761., 17142., 17534., 17940., 18367., 18821.,
107    19312., 19849., 20438., 21072., 21724., 22360., 22957., 23511., 24030., 24525., 25008.,
108    25491., 25987., 26514., 27098., 27791., 28730., 31416.}, /* En=12140 keV */
109   {0.000,  2425.,  3423.,  4211.,  4909.,  5563.,  6200.,  6829.,  7455.,  8070.,  8664.,
110    9226.,  9751.,  10239., 10694., 11120., 11522., 11904., 12271., 12625., 12969., 13305.,
111    13636., 13962., 14287., 14610., 14934., 15261., 15591., 15928., 16273., 16630., 17002.,
112    17393., 17810., 18263., 18764., 19330., 19985., 20742., 21575., 22399., 23153., 23836.,
113    24471., 25085., 25704., 26365., 27125., 28138., 31416.}, /* En=12320 keV */
114   {0.000,  2661.,  3692.,  4487.,  5182.,  5829.,  6457.,  7081.,  7706.,  8329.,  8933.,
115    9504.,  10030., 10509., 10946., 11346., 11717., 12063., 12390., 12701., 13000., 13289.,
116    13570., 13846., 14118., 14388., 14658., 14929., 15203., 15482., 15767., 16062., 16370.,
117    16694., 17041., 17417., 17834., 18311., 18879., 19596., 20545., 21666., 22662., 23468.,
118    24158., 24791., 25408., 26049., 26769., 27705., 31416.}, /* En=12570 keV */
119   {0.000,  1941.,  2896.,  3780.,  4713.,  5731.,  6720.,  7532.,  8172.,  8693.,  9137.,
120    9528.,  9880.,  10204., 10506., 10792., 11064., 11325., 11578., 11824., 12064., 12301.,
121    12534., 12766., 12997., 13227., 13458., 13691., 13926., 14165., 14409., 14659., 14916.,
122    15182., 15461., 15753., 16064., 16399., 16765., 17175., 17647., 18217., 18962., 20057.,
123    21618., 23016., 24127., 25137., 26190., 27504., 31416.}, /* En=12940 keV */
124   {0.000,  2028.,  2996.,  3875.,  4801.,  5899.,  7158.,  8182.,  8905.,  9458.,  9914.,
125    10310., 10665., 10992., 11298., 11588., 11868., 12138., 12403., 12664., 12922., 13179.,
126    13436., 13695., 13957., 14223., 14494., 14772., 15057., 15352., 15658., 15977., 16310.,
127    16660., 17031., 17425., 17848., 18305., 18805., 19358., 19977., 20675., 21451., 22281.,
128    23121., 23943., 24754., 25586., 26507., 27691., 31416.}, /* En=13420 keV */
129   {0.000,  2099.,  3016.,  3762.,  4437.,  5088.,  5747.,  6450.,  7233.,  8116.,  9014.,
130    9779.,  10388., 10886., 11312., 11690., 12034., 12354., 12657., 12948., 13231., 13509.,
131    13784., 14060., 14337., 14619., 14908., 15206., 15516., 15842., 16187., 16554., 16950.,
132    17378., 17846., 18361., 18932., 19572., 20296., 21115., 22008., 22896., 23707., 24434.,
133    25101., 25740., 26380., 27056., 27825., 28822., 31416.}, /* En=13760 keV */
134   {0.000,  2216.,  3042.,  3663.,  4192.,  4675.,  5137.,  5595.,  6065.,  6565.,  7121.,
135    7768.,  8547.,  9421.,  10221., 10870., 11403., 11859., 12264., 12634., 12980., 13309.,
136    13626., 13935., 14240., 14543., 14846., 15152., 15463., 15783., 16113., 16456., 16817.,
137    17199., 17610., 18056., 18549., 19104., 19743., 20491., 21350., 22257., 23110., 23873.,
138    24563., 25214., 25858., 26532., 27296., 28301., 31416.}, /* En=14020 keV */
139   {0.000,  2249.,  3077.,  3698.,  4228.,  4711.,  5172.,  5628.,  6094.,  6586.,  7123.,
140    7731.,  8436.,  9232.,  10033., 10745., 11353., 11879., 12348., 12776., 13175., 13553.,
141    13916., 14268., 14613., 14955., 15296., 15639., 15986., 16342., 16709., 17090., 17492.,
142    17919., 18379., 18881., 19436., 20053., 20732., 21456., 22190., 22904., 23585., 24235.,
143    24868., 25501., 26155., 26860., 27672., 28727., 31416.}, /* En=14200 keV */
144   {0.000,  2206.,  3020.,  3628.,  4143.,  4607.,  5043.,  5466.,  5884.,  6307.,  6745.,
145    7205.,  7698.,  8237.,  8833.,  9491.,  10201., 10929., 11634., 12293., 12903., 13472.,
146    14009., 14525., 15026., 15519., 16009., 16500., 16995., 17496., 18005., 18524., 19053.,
147    19590., 20134., 20681., 21225., 21762., 22289., 22805., 23314., 23821., 24333., 24858.,
148    25408., 25997., 26643., 27368., 28202., 29223., 31416.}, /* En=14440 keV */
149   {0.000,  2489.,  3349.,  3983.,  4517.,  4998.,  5447.,  5879.,  6304.,  6729.,  7160.,
150    7600.,  8054.,  8523.,  9008.,  9508.,  10019., 10539., 11065., 11593., 12124., 12657.,
151    13192., 13731., 14275., 14826., 15388., 15964., 16559., 17177., 17821., 18488., 19164.,
152    19832., 20472., 21076., 21644., 22180., 22692., 23187., 23673., 24156., 24643., 25141.,
153    25657., 26203., 26790., 27441., 28196., 29158., 31416.}, /* En=14620 keV */
154   {0.000,  2945.,  3913.,  4642.,  5270.,  5848.,  6400.,  6937.,  7468.,  7992.,  8510.,
155    9019.,  9517.,  10001., 10472., 10930., 11375., 11811., 12238., 12658., 13073., 13486.,
156    13898., 14312., 14731., 15159., 15600., 16059., 16544., 17063., 17627., 18245., 18921.,
157    19641., 20371., 21073., 21726., 22329., 22888., 23413., 23913., 24397., 24871., 25344.,
158    25825., 26323., 26854., 27442., 28132., 29045., 31416.}, /* En=14820 keV */
159   {0.000,  3308.,  4355.,  5183.,  5930.,  6638.,  7311.,  7938.,  8510.,  9030.,  9504.,
160    9941.,  10350., 10738., 11111., 11474., 11831., 12184., 12538., 12894., 13255., 13622.,
161    13996., 14380., 14773., 15177., 15594., 16024., 16472., 16941., 17437., 17970., 18549.,
162    19186., 19878., 20604., 21317., 21982., 22593., 23158., 23691., 24203., 24706., 25211.,
163    25728., 26270., 26850., 27493., 28240., 29191., 31416.}, /* En=15050 keV */
164   {0.000,  2774.,  3718.,  4452.,  5121.,  5796.,  6548.,  7503.,  8803.,  9927.,  10670.,
165    11225., 11684., 12085., 12450., 12790., 13113., 13425., 13729., 14030., 14329., 14629.,
166    14933., 15242., 15558., 15885., 16223., 16575., 16943., 17330., 17738., 18170., 18628.,
167    19115., 19633., 20184., 20765., 21375., 22002., 22634., 23260., 23871., 24467., 25050.,
168    25629., 26216., 26825., 27482., 28229., 29174., 31416.}, /* En=15660 keV */
169   {0.000,  2572.,  3380.,  3971.,  4470.,  4924.,  5357.,  5787.,  6230.,  6709.,  7260.,
170    7971.,  9158.,  10803., 11664., 12234., 12684., 13071., 13419., 13742., 14048., 14344.,
171    14633., 14919., 15207., 15499., 15798., 16110., 16437., 16787., 17167., 17587., 18060.,
172    18603., 19225., 19919., 20643., 21351., 22020., 22648., 23246., 23823., 24387., 24948.,
173    25516., 26099., 26713., 27381., 28144., 29113., 31416.}, /* En=15980 keV */
174   {0.000,  2876.,  3700.,  4293.,  4786.,  5221.,  5622.,  6001.,  6366.,  6726.,  7086.,
175    7452.,  7832.,  8235.,  8674.,  9173.,  9772.,  10542., 11506., 12399., 13077., 13612.,
176    14066., 14471., 14845., 15199., 15542., 15880., 16219., 16564., 16919., 17292., 17689.,
177    18123., 18606., 19163., 19825., 20623., 21523., 22383., 23124., 23768., 24352., 24906.,
178    25454., 26017., 26620., 27294., 28088., 29106., 31416.}, /* En=16470 keV */
179   {0.000,  3534.,  4408.,  5069.,  5638.,  6157.,  6649.,  7126.,  7595.,  8063.,  8537.,
180    9020.,  9518.,  10033., 10565., 11110., 11658., 12195., 12707., 13189., 13639., 14059.,
181    14453., 14827., 15182., 15525., 15856., 16181., 16500., 16817., 17134., 17453., 17779.,
182    18115., 18465., 18836., 19237., 19682., 20194., 20809., 21581., 22515., 23441., 24240.,
183    24941., 25597., 26249., 26940., 27735., 28775., 31416.}, /* En=16940 keV */
184   {0.000,  3063.,  3895.,  4531.,  5092.,  5621.,  6144.,  6674.,  7219.,  7776.,  8334.,
185    8879.,  9398.,  9888.,  10352., 10795., 11221., 11638., 12048., 12455., 12863., 13272.,
186    13685., 14100., 14517., 14935., 15353., 15768., 16181., 16590., 16995., 17397., 17797.,
187    18195., 18594., 18994., 19398., 19810., 20232., 20671., 21133., 21628., 22168., 22773.,
188    23466., 24265., 25162., 26124., 27154., 28382., 31416.}, /* En=17970 keV */
189   {0.000,  2839.,  4027.,  4951.,  5736.,  6437.,  7074.,  7671.,  8230.,  8765.,  9273.,
190    9767.,  10241., 10703., 11152., 11595., 12026., 12453., 12871., 13285., 13697., 14102.,
191    14507., 14909., 15308., 15711., 16110., 16509., 16911., 17316., 17721., 18133., 18547.,
192    18965., 19393., 19823., 20266., 20715., 21177., 21651., 22145., 22654., 23188., 23747.,
193    24344., 24981., 25682., 26467., 27391., 28579., 31416.}, /* En=18000 keV */
194   {0.000,  2839.,  4027.,  4951.,  5736.,  6437.,  7074.,  7671.,  8230.,  8765.,  9273.,
195    9767.,  10241., 10703., 11152., 11595., 12026., 12453., 12871., 13285., 13697., 14102.,
196    14507., 14909., 15308., 15711., 16110., 16509., 16911., 17316., 17721., 18133., 18547.,
197    18965., 19393., 19823., 20266., 20715., 21177., 21651., 22145., 22654., 23188., 23747.,
198    24344., 24981., 25682., 26467., 27391., 28579., 31416.}, /* En=19000 keV */
199   {0.000,  2839.,  4027.,  4951.,  5736.,  6437.,  7074.,  7671.,  8230.,  8765.,  9273.,
200    9767.,  10241., 10703., 11152., 11595., 12026., 12453., 12871., 13285., 13697., 14102.,
201    14507., 14909., 15308., 15711., 16110., 16509., 16911., 17316., 17721., 18133., 18547.,
202    18965., 19393., 19823., 20266., 20715., 21177., 21651., 22145., 22654., 23188., 23747.,
203    24344., 24981., 25682., 26467., 27391., 28579., 31416.} /* En=20000 keV*/
204 };
205 
206 void G4NRESP71M03::DKINMA(G4ReactionProduct* p1, G4ReactionProduct* p2, G4ReactionProduct* p3,
207                           G4ReactionProduct* p4, const G4double Q, const G4double costhcm3)
208 {
209   G4ReactionProduct CM;
210   G4double TotalEnergyCM;
211 
212   if (p2 != nullptr)  // If it is not a decay reaction...
213   {
214     // Calculating (total momentum, energy and mass) of the center of mass.
215     const G4ThreeVector TotalMomentumLAB = p1->GetMomentum() + p2->GetMomentum();
216     CM.SetMomentum(TotalMomentumLAB);
217 
218     const G4double TotalEnergyLAB = p1->GetTotalEnergy() + p2->GetTotalEnergy();
219     CM.SetTotalEnergy(TotalEnergyLAB);
220 
221     CM.SetMass(std::sqrt(TotalEnergyLAB * TotalEnergyLAB - TotalMomentumLAB * TotalMomentumLAB));
222 
223     // Transforming primary particles from laboratory to center of mass system.
224     p1->Lorentz(*p1, CM);
225     p2->Lorentz(*p2, CM);
226 
227     TotalEnergyCM = p1->GetTotalEnergy() + p2->GetTotalEnergy();
228 
229     const G4double m4 =
230       (p1->GetMass() + p2->GetMass())
231       - (p3->GetMass()
232          + Q);  // Mass of the residual nucleus in the excited state (not in the ground state).
233     p4->SetMass(m4);
234   }
235   else  // If it is a decay reaction...
236   {
237     const G4ThreeVector TotalMomentumLAB = p1->GetMomentum();
238     CM.SetMomentum(TotalMomentumLAB);
239 
240     const G4double TotalEnergyLAB = p1->GetTotalEnergy();
241     CM.SetTotalEnergy(TotalEnergyLAB);
242 
243     CM.SetMass(std::sqrt(TotalEnergyLAB * TotalEnergyLAB - TotalMomentumLAB * TotalMomentumLAB));
244 
245     // Transforming primary particles from laboratory to center of mass system (not really necessary
246     // in this case).
247     p1->Lorentz(*p1, CM);
248 
249     const G4double m4 =
250       p1->GetMass()
251       - (p3->GetMass()
252          + Q);  // Mass of the residual nucleus in the excited state (not in the ground state).
253     p4->SetMass(m4);
254 
255     TotalEnergyCM = p1->GetTotalEnergy();
256   }
257 
258   // Calculating momentum and total energy of the reaction products.
259 
260   const G4ThreeVector p1unit = p1->GetMomentum().unit();
261 
262   G4RotationMatrix rot(std::acos(p1unit * G4ThreeVector(0., 1., 0.)),
263                        std::acos(p1unit * G4ThreeVector(0., 0., 1.)), 0.);
264   rot.inverse();
265 
266   const G4double theta = std::acos(costhcm3);
267   const G4double phi = twopi * G4UniformRand();
268 
269   const G4double Energy3CM =
270     (std::pow(TotalEnergyCM, 2.) + std::pow(p3->GetMass(), 2.) - std::pow(p4->GetMass(), 2.))
271     / (2. * TotalEnergyCM);
272   p3->SetTotalEnergy(Energy3CM);
273 
274   const G4double Momentum3CM = std::sqrt(std::pow(Energy3CM, 2.) - std::pow(p3->GetMass(), 2.));
275   p3->SetMomentum(rot
276                   * G4ThreeVector(Momentum3CM * std::sin(theta) * std::cos(phi),
277                                   Momentum3CM * std::sin(theta) * std::sin(phi),
278                                   Momentum3CM * costhcm3));
279 
280   const G4double Energy4CM = TotalEnergyCM - Energy3CM;
281   p4->SetTotalEnergy(Energy4CM);
282 
283   const G4double Momentum4CM = std::sqrt(std::pow(Energy4CM, 2.) - std::pow(p4->GetMass(), 2.));
284   p4->SetMomentum(-Momentum4CM * p3->GetMomentum().unit());
285 
286   // Transforming reaction products from center of mass to laboratory system.
287   p3->Lorentz(*p3, -1. * CM);
288   p4->Lorentz(*p4, -1. * CM);
289 }
290 
291 G4int G4NRESP71M03::ApplyMechanismI_NBeA2A(G4ReactionProduct& neut, G4ReactionProduct& carb,
292                                            G4ReactionProduct* theProds, const G4double QI)
293 {
294   // N+12C --> A+9BE*
295   G4ReactionProduct p4;
296 
297   theProds[0].SetDefinition(G4Alpha::Alpha());
298 
299   DKINMA(&neut, &carb, &(theProds[0]), &p4, QI, 2. * G4UniformRand() - 1.);
300 
301   // 9BE* --> N+8BE
302   G4ReactionProduct p1(p4);
303 
304   theProds[1].SetDefinition(G4Neutron::Neutron());
305 
306   DKINMA(&p1, nullptr, &(theProds[1]), &p4, -QI - 7.369, 2. * G4UniformRand() - 1.);
307 
308   // 8BE --> 2*A
309   p1 = p4;
310 
311   theProds[2].SetDefinition(G4Alpha::Alpha());
312   theProds[3].SetDefinition(G4Alpha::Alpha());
313 
314   DKINMA(&p1, nullptr, &(theProds[2]), &(theProds[3]), 0.09538798439007223351,
315          2. * G4UniformRand() - 1.);
316 
317   return 0;
318 }
319 
320 // Apply kinematics for excited level selected by GEANT4 (it).
321 G4int G4NRESP71M03::ApplyMechanismII_ACN2A(G4ReactionProduct& neut, G4ReactionProduct& carb,
322                                            G4ReactionProduct* theProds, const G4double QI)
323 {
324   // 12C(N,N')12C'
325   G4ReactionProduct p4;
326 
327   theProds[0].SetDefinition(G4Neutron::Neutron());
328 
329   DKINMA(&neut, &carb, &(theProds[0]), &p4, QI, 2. * G4UniformRand() - 1.);
330 
331   // 12C' --> ALPHA+8BE'
332   G4ReactionProduct p1(p4);
333 
334   theProds[1].SetDefinition(G4Alpha::Alpha());
335 
336   DKINMA(&p1, nullptr, &(theProds[1]), &p4, -QI - 7.369, 2. * G4UniformRand() - 1.);
337 
338   // 8BE --> 2*ALPHA
339   p1 = p4;
340 
341   theProds[2].SetDefinition(G4Alpha::Alpha());
342   theProds[3].SetDefinition(G4Alpha::Alpha());
343 
344   DKINMA(&p1, nullptr, &(theProds[2]), &(theProds[3]), 0.09538798439007223351,
345          2. * G4UniformRand() - 1.);
346 
347   return 0;
348 }
349 
350 G4int G4NRESP71M03::ApplyMechanismABE(G4ReactionProduct& neut, G4ReactionProduct& carb,
351                                       G4ReactionProduct* theProds)
352 {
353   G4double CosThetaCMAlpha = 0.;  // Cosine of the angle of emission of the alpha particle (theta).
354 
355   G4double Kn = neut.GetKineticEnergy();  // Neutron energy in the center of mass system.
356   if (Kn > 5.7 * MeV) {
357     // Sorting.
358     for (G4int i = 1; i < ndist; i++) {
359       if (BEN2[i] >= Kn / keV) {
360         // Ok. The neutron energy is between values i-1 and i of BEN2. Taking them.
361         G4double BE1 = BEN2[i - 1];
362         G4double BE2 = BEN2[i];
363 
364         // Performing energy and angle interpolation.
365 
366         G4double FRA =
367           G4UniformRand()
368           * 49.99999999;  // Sorting the bin of the cumulative probability FRA (Rho). There are 51
369                           // values of Rho from 0 to 1 (0; 0.02; 0.04 ... 1).
370         G4double DJTETA =
371           FRA
372           - G4int(FRA);  // Distance in bin units (DJTETA) from the low edge of the bin with Rho.
373         G4int JTETA =
374           G4int(FRA) + 1;  // Getting the bin (JTETA) next to the bin with the value of Rho.
375 
376         // Calculating the angle from the cumulative distribution at energy i.
377 
378         G4double B11 = B2[i - 1][JTETA - 1];
379         G4double B12 = B2[i - 1][JTETA];
380 
381         G4double TCM1 = B11 + (B12 - B11) * DJTETA;  // Angle at energy i.
382 
383         // Calculating the angle from the cumulative distribution at energy i-1.
384 
385         G4double B21 = B2[i][JTETA - 1];
386         G4double B22 = B2[i][JTETA];
387 
388         G4double TCM2 = B21 + (B22 - B21) * DJTETA;  // Angle at energy i-1.
389 
390         // Interpolating the angle between energy values i and i-1.
391         G4double TCM = (TCM1 + (TCM2 - TCM1) * (Kn / keV - BE1) / (BE2 - BE1)) * 1.E-4;
392         CosThetaCMAlpha = std::cos(TCM);
393 
394         break;
395       }
396     }
397   }
398   else {
399     // Isotropic distribution in CM.
400     CosThetaCMAlpha = 1. - 2. * G4UniformRand();
401   }
402 
403   // N+12C --> A+9BE
404   theProds[0].SetDefinition(G4Alpha::Alpha());
405   theProds[1].SetDefinition(G4IonTable::GetIonTable()->GetIon(4, 9, 0.));
406 
407   DKINMA(&neut, &carb, &(theProds[0]), &(theProds[1]), -5.71 * MeV, CosThetaCMAlpha);
408 
409   return 0;
410 }
411