Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RandomDirection.hh File Reference
#include <CLHEP/Units/PhysicalConstants.h>
#include "globals.hh"
#include "Randomize.hh"
#include "G4ThreeVector.hh"

Go to the source code of this file.

Functions

G4ThreeVector G4RandomDirection ()
 

Function Documentation

◆ G4RandomDirection()

G4ThreeVector G4RandomDirection ( )
inline

Definition at line 55 of file G4RandomDirection.hh.

56{
57 G4double cosTheta = 2.*G4UniformRand()-1.;
58 G4double sinTheta2 = 1. - cosTheta*cosTheta;
59 if( sinTheta2 < 0.) sinTheta2 = 0.;
60 G4double sinTheta = std::sqrt(sinTheta2);
61 G4double phi = CLHEP::twopi*G4UniformRand();
62 return G4ThreeVector(sinTheta*std::cos(phi),
63 sinTheta*std::sin(phi), cosTheta).unit();
64}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const

Referenced by G4MuMinusCapturePrecompound::ApplyYourself(), G4MuonMinusBoundDecay::ApplyYourself(), G4QNucleus::ChooseFermiMomenta(), G4QNucleus::ChoosePositions(), G4PionDecayMakeSpin::DaughterPolarization(), G4LambertianRand(), G4QInelastic::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QDiffractionRatio::ProjFragment(), G4HeatedKleinNishinaCompton::SampleSecondaries(), G4eeTo3PiModel::SampleSecondaries(), and G4QDiffractionRatio::TargFragment().