Ignore:
Timestamp:
Jun 18, 2010, 11:42:07 AM (14 years ago)
Author:
garnier
Message:

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File:
1 edited

Legend:

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

    r962 r1315  
    2323// * acceptance of all terms of the Geant4 Software license.          *
    2424// ********************************************************************
     25// $Id: G4InuclSpecialFunctions.cc,v 1.19 2010/03/19 05:03:23 mkelsey Exp $
     26// Geant4 tag: $Name: geant4-09-04-beta-cand-01 $
    2527//
     28// 20100114  M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
     29
    2630#include "G4InuclSpecialFunctions.hh"
     31#include "G4LorentzVector.hh"
     32#include "G4ThreeVector.hh"
    2733#include "Randomize.hh"
     34#include <cmath>
    2835
    2936G4double G4InuclSpecialFunctions::getAL(G4double A) {
     
    3441  }
    3542
    36   return 0.76 + 2.2 / std::pow(A, 0.333333);
     43  return 0.76 + 2.2 / G4cbrt(A);
    3744}
    3845
     
    8491
    8592  const G4double C = 55.4;
    86   G4double Ef;
     93  G4double arg = (ntype==0) ? (A-Z)/A : Z/A;
    8794
    88   if (ntype == 0) {
    89     Ef = C * std::pow((A - Z) / A, 0.666667);
     95  return C * G4cbrt(arg*arg);   // 2/3 power
     96}
    9097
    91   } else {
    92     Ef = C * std::pow(Z / A, 0.666667);
    93   };
    94 
    95   return Ef;
     98G4double G4InuclSpecialFunctions::G4cbrt(G4double x) {
     99  return x==0 ? 0. : (x<0?-1.:1.)*std::exp(std::log(std::fabs(x))/3.);
    96100}
    97101
     
    150154}
    151155
    152 G4CascadeMomentum G4InuclSpecialFunctions::generateWithFixedTheta(G4double ct,
    153                                                                         G4double p) {
    154   G4int verboseLevel = 2;
    155 
     156G4LorentzVector
     157G4InuclSpecialFunctions::generateWithFixedTheta(G4double ct, G4double p,
     158                                                G4double m) {
     159  const G4int verboseLevel = 0;
    156160  if (verboseLevel > 3) {
    157161    G4cout << " >>> G4InuclSpecialFunctions::generateWithFixedTheta" << G4endl;
    158162  }
    159163
    160   G4CascadeMomentum momr;
    161164  G4double phi = randomPHI();
    162165  G4double pt = p * std::sqrt(std::fabs(1.0 - ct * ct));
    163   //  not used:  G4CascadeMomentum mom1;
    164   momr[1] = pt * std::cos(phi);
    165   momr[2] = pt * std::sin(phi);
    166   momr[3] = p * ct;
     166
     167  static G4ThreeVector pvec;    // Buffers to avoid memory thrashing
     168  static G4LorentzVector momr;
     169
     170  pvec.set(pt*std::cos(phi), pt*std::sin(phi), p*ct);
     171  momr.setVectM(pvec, m);
    167172
    168173  return momr;
    169174}
     175
     176G4LorentzVector
     177G4InuclSpecialFunctions::generateWithRandomAngles(G4double p, G4double m) {
     178  std::pair<G4double, G4double> COS_SIN = randomCOS_SIN();
     179  G4double phi = randomPHI();
     180  G4double pt = p * COS_SIN.second;
     181 
     182  static G4ThreeVector pvec;    // Buffers to avoid memory thrashing
     183  static G4LorentzVector momr;
     184
     185  pvec.set(pt*std::cos(phi), pt*std::sin(phi), p*COS_SIN.first);
     186  momr.setVectM(pvec, m);
     187
     188  return momr;
     189}
Note: See TracChangeset for help on using the changeset viewer.