Ignore:
Timestamp:
Nov 5, 2010, 3:45:55 PM (14 years ago)
Author:
garnier
Message:

update ti head

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/source/processes/hadronic/models/cascade/cascade/src/G4InuclSpecialFunctions.cc

    r1337 r1340  
    2424// ********************************************************************
    2525//
    26 // $Id: G4InuclSpecialFunctions.cc,v 1.20 2010/06/25 09:44:46 gunter Exp $
    27 // Geant4 tag: $Name: geant4-09-04-beta-01 $
     26// $Id: G4InuclSpecialFunctions.cc,v 1.21 2010/09/14 17:51:36 mkelsey Exp $
     27// Geant4 tag: $Name: hadr-casc-V09-03-85 $
    2828//
    2929// 20100114  M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
     30// 20100914  M. Kelsey -- Migrate to integer A and Z.  Discard pointless
     31//              verbosity.
    3032
    3133#include "G4InuclSpecialFunctions.hh"
     
    3537#include <cmath>
    3638
    37 G4double G4InuclSpecialFunctions::getAL(G4double A) {
    38   G4int verboseLevel = 2;
    3939
    40   if (verboseLevel > 3) {
    41     G4cout << " >>> G4InuclSpecialFunctions::getAL" << G4endl;
    42   }
    43 
     40G4double G4InuclSpecialFunctions::getAL(G4int A) {
    4441  return 0.76 + 2.2 / G4cbrt(A);
    4542}
    4643
    4744G4double G4InuclSpecialFunctions::csNN(G4double e) {
    48   G4int verboseLevel = 2;
    49  
    50   if (verboseLevel > 3) {
    51     G4cout << " >>> G4InuclSpecialFunctions::csNN" << G4endl;
    52   }
    53 
    5445  G4double snn;
    5546
    5647  if (e < 40.0) {
    5748    snn = -1174.8 / (e * e) + 3088.5 / e + 5.3107;
    58 
    5949  } else {
    6050    snn = 93074.0 / (e * e) - 11.148 / e + 22.429;
    61   };
     51  }
    6252
    6353  return snn;
     
    6555
    6656G4double G4InuclSpecialFunctions::csPN(G4double e) {
    67   G4int verboseLevel = 2;
    68 
    69   if (verboseLevel > 3) {
    70     G4cout << " >>> G4InuclSpecialFunctions::csPN" << G4endl;
    71   }
    72 
    7357  G4double spn;
    7458
    7559  if (e < 40.0) {
    7660    spn = -5057.4 / (e * e) + 9069.2 / e + 6.9466;
    77 
    7861  } else {
    7962    spn = 239380.0 / (e * e) + 1802.0 / e + 27.147;
    80   };
     63  }
    8164
    8265  return spn;
    8366}
    8467
    85 G4double G4InuclSpecialFunctions::FermiEnergy(G4double A, G4double Z, G4int ntype) {
    86   // calculates the nuclei Fermi energy for 0 - neutron and 1 - proton
    87   G4int verboseLevel = 2;
     68// calculates the nuclei Fermi energy for 0 - neutron and 1 - proton
    8869
    89   if (verboseLevel > 3) {
    90     G4cout << " >>> G4InuclSpecialFunctions::FermiEnergy" << G4endl;
    91   }
    92 
     70G4double G4InuclSpecialFunctions::FermiEnergy(G4int A, G4int Z, G4int ntype) {
    9371  const G4double C = 55.4;
    94   G4double arg = (ntype==0) ? (A-Z)/A : Z/A;
     72  G4double arg = (ntype==0) ? G4double(A-Z)/A : G4double(Z)/A;
    9573
    9674  return C * G4cbrt(arg*arg);   // 2/3 power
     
    10280
    10381G4double G4InuclSpecialFunctions::inuclRndm() {
    104   G4int verboseLevel = 2;
    105 
    106   if (verboseLevel > 3) {
    107     G4cout << " >>> G4InuclSpecialFunctions::inuclRndm" << G4endl;
    108   }
    109 
    11082  G4double rnd = G4UniformRand();
    111 
    11283  return rnd;
    11384}
    11485
    11586G4double G4InuclSpecialFunctions::randomGauss(G4double sigma) {
    116   G4int verboseLevel = 2;
    117 
    118   if (verboseLevel > 3) {
    119     G4cout << " >>> G4InuclSpecialFunctions::randomGauss" << G4endl;
    120   }
    121 
    12287  const G4double eps = 1.0e-6;
    12388  const G4double twopi = 6.2831854;
     
    13297
    13398G4double G4InuclSpecialFunctions::randomPHI() {
    134   G4int verboseLevel = 2;
    135 
    136   if (verboseLevel > 3) {
    137     G4cout << " >>> G4InuclSpecialFunctions::randomPHI" << G4endl;
    138   }
    139 
    14099  const G4double twopi = 6.2831853;
    141 
    142100  return twopi * inuclRndm();
    143101}
    144102
    145103std::pair<G4double, G4double> G4InuclSpecialFunctions::randomCOS_SIN() {
    146   G4int verboseLevel = 2;
    147 
    148   if (verboseLevel > 3) {
    149     G4cout << " >>> G4InuclSpecialFunctions::randomCOS_SIN" << G4endl;
    150   }
    151 
    152104  G4double CT = 1.0 - 2.0 * inuclRndm();
    153105
    154   return std::pair<G4double, G4double>(CT, std::sqrt(1.0 - CT * CT));
     106  return std::pair<G4double, G4double>(CT, std::sqrt(1.0 - CT*CT));
    155107}
    156108
     
    158110G4InuclSpecialFunctions::generateWithFixedTheta(G4double ct, G4double p,
    159111                                                G4double m) {
    160   const G4int verboseLevel = 0;
    161   if (verboseLevel > 3) {
    162     G4cout << " >>> G4InuclSpecialFunctions::generateWithFixedTheta" << G4endl;
    163   }
    164 
    165112  G4double phi = randomPHI();
    166113  G4double pt = p * std::sqrt(std::fabs(1.0 - ct * ct));
Note: See TracChangeset for help on using the changeset viewer.