source: trunk/source/global/HEPRandom/test/G4RandomDirectionTest.cc @ 1347

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

File size: 7.7 KB
Line 
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//
27// $Id: G4RandomDirectionTest.cc,v 1.1 2008/03/19 16:53:55 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-04-ref-00 $
29//
30//
31// Test for random direction unit vector algorithm
32//
33// Author: V. Grichine
34//
35//
36// History:
37//
38// 19.03.08 - First implementation
39// ----------------------------------------------------------------------
40
41#include "G4ios.hh"
42#include <fstream>
43#include <cmath>
44#include <iomanip>
45
46#include "globals.hh"
47#include "Randomize.hh"
48#include "G4RandomDirection.hh"
49#include "G4UnitsTable.hh"
50#include "G4Timer.hh"
51
52///////////////////////////////////////////////////////////////////////
53//
54// Random algorithm from Geant3 RAND3
55
56G4ThreeVector IsotropicCubeRand() 
57{
58  /* Returns a random isotropic unit vector. */
59
60  G4ThreeVector vect;
61  G4double len2;
62
63  do {
64
65    vect.setX(G4UniformRand() - 0.5);
66    vect.setY(G4UniformRand() - 0.5);
67    vect.setZ(G4UniformRand() - 0.5);
68
69    len2 = vect.mag2();
70
71  } while (len2 < 0.01 || len2 > 0.25);
72
73  return vect.unit();
74}
75
76/////////////////////////////////////////////////////////////////////////////
77//
78// Random distribution over unit radius sphere
79
80G4ThreeVector IsotropicSphereRand() 
81{
82  G4double cosTheta  = 2*G4UniformRand()-1.;
83  G4double sinTheta2 = 1. - cosTheta*cosTheta;
84  if( sinTheta2 < 0.)  sinTheta2 = 0.;
85  G4double sinTheta  = std::sqrt(sinTheta2); 
86  G4double phi       = twopi*G4UniformRand();
87  return G4ThreeVector(sinTheta*std::cos(phi),
88                       sinTheta*std::sin(phi), cosTheta).unit(); 
89}
90
91///////////////////////////////////////////////////////////////////////////////
92//
93// 8 quadrants algorithm
94
95G4ThreeVector MKRandomDirection()
96{
97  // Randomization in one of 8 Quadrants (x>0, y>0, z>0)
98  //
99  G4double x=G4UniformRand(), y=G4UniformRand(), z=G4UniformRand();
100  G4double r2= x*x+y*y+z*z;
101  while(r2>1.||r2<.000001)
102  {
103    x = G4UniformRand(); y = G4UniformRand(); z = G4UniformRand();
104    r2=x*x+y*y+z*z;
105  }
106  G4double r=std::sqrt(r2), quad=G4UniformRand();
107
108  if(quad>0.5)
109  {
110    if(quad>0.75)
111    {
112      if(quad>0.875)    return G4ThreeVector(-x/r,-y/r,-z/r);
113      else              return G4ThreeVector(-x/r,-y/r, z/r);
114    }
115    else
116    {
117      if(quad>0.625)    return G4ThreeVector(-x/r, y/r,-z/r);
118      else              return G4ThreeVector(-x/r, y/r, z/r);
119    }
120  }
121  else
122  {
123    if(quad>0.25)
124    {
125      if(quad>0.375)    return G4ThreeVector( x/r,-y/r,-z/r);
126      else              return G4ThreeVector( x/r,-y/r, z/r);
127    }
128    else if(quad>0.125) return G4ThreeVector( x/r, y/r,-z/r);
129  }
130  return                       G4ThreeVector( x/r, y/r, z/r);
131}
132
133
134//////////////////////////////////////////////////////////////////////////////
135//
136// Test main program
137
138int main()
139{
140  G4int i, iMax = 20;
141
142  G4Timer timer;
143
144  iMax = 1000000;
145
146  timer.Start();
147  for( i = 0; i < iMax; i++ )
148  {
149    G4ThreeVector isoVectr = IsotropicCubeRand(); 
150  }
151  timer.Stop();
152  G4cout<<"Total time of volume "<<iMax<<" calls = "
153        <<timer.GetUserElapsed()<<" s"<<G4endl<<G4endl;
154
155  timer.Start();
156  for( i = 0; i < iMax; i++ )
157  {
158    G4ThreeVector isoVectr = IsotropicSphereRand(); 
159  }
160  timer.Stop();
161  G4cout<<"Total time of surface "<<iMax<<" calls = "
162        <<timer.GetUserElapsed()<<" s"<<G4endl<<G4endl;
163 
164  timer.Start();
165  for( i = 0; i < iMax; i++ )
166  {
167    G4ThreeVector isoVectr = G4RandomDirection(); 
168  }
169  timer.Stop();
170  G4cout<<"Total time of G4RandomDirection() "<<iMax
171        <<" calls = "<<timer.GetUserElapsed()
172        <<" s"<<G4endl<<G4endl;
173
174  timer.Start();
175  for( i = 0; i < iMax; i++ )
176  {
177    G4ThreeVector isoVectr = MKRandomDirection(); 
178  }
179  timer.Stop();
180  G4cout<<"Total time of MKRandomDirection() "<<iMax
181        <<" calls = "<<timer.GetUserElapsed()
182        <<" s"<<G4endl<<G4endl;
183
184  iMax = 1000000;
185  G4int j, jMax = 100;
186  G4double cosThetaDistr[100], phi[100];
187
188  for( j = 0; j < jMax; j++ )
189  {
190    cosThetaDistr[j] = 0.; 
191    phi[j]           = 0.; 
192  }   
193  G4double xyPlane, phiNow, cosThetaNow, cosThetaTmp, phiTmp;
194 
195  for( i = 0; i < iMax; i++ )
196  {
197    // G4ThreeVector isoVectr = IsotropicCubeRand();
198    // G4ThreeVector isoVectr = IsotropicSphereRand();
199    G4ThreeVector isoVectr = G4RandomDirection(); 
200
201    xyPlane = std::sqrt( isoVectr.x()*isoVectr.x()
202                       + isoVectr.y()*isoVectr.y() );
203    if ( xyPlane ) 
204    {
205      phiNow  = std::atan2(isoVectr.y(),isoVectr.x());
206      phiNow += pi;                                     // 0-twopi range
207
208      cosThetaNow = isoVectr.z();
209    }
210    else
211    {
212      if ( isoVectr.z() >= 0. ) cosThetaNow =  1.;
213      else                      cosThetaNow = -1.;
214      phiNow = twopi*G4UniformRand();     
215    }
216    for( j = 0; j < jMax; j++ )
217    {
218      cosThetaTmp = -1. + 2.*j/jMax;
219
220      if( cosThetaTmp >= cosThetaNow )
221      {
222        cosThetaDistr[j] += 1.;
223        break;
224      }
225    }   
226    for( j = 0; j < jMax; j++ )
227    {
228      phiTmp      = twopi*j/jMax;
229      if( phiTmp >= phiNow )
230      {
231        phi[j] += 1.;
232        break;
233      }
234    }     
235  }
236  G4cout << G4endl;
237  G4cout <<"cosThetaTmp"<<"\t"<<"cosThetaDistr[j]"<<"\t"
238         <<"phi/degree"<<"\t"<<"phi[j]"<<G4endl;
239  G4cout << G4endl;
240
241  for( j = 0; j < jMax; j++ )
242  {
243    cosThetaTmp = -1. + 2.*j/jMax;
244    phiTmp      = twopi*j/jMax;
245    G4cout <<cosThetaTmp<<"\t"<<cosThetaDistr[j]<<"\t"
246           <<phiTmp/degree<<"\t"<<phi[j]<<G4endl;
247  }
248   
249  G4double mean=0., rms2=0., rms, relDelta;
250
251  for( j = 0; j < jMax; j++ )
252  {
253    mean += cosThetaDistr[j];
254  }
255  mean /= jMax;
256
257  for( j = 0; j < jMax; j++ )
258  {
259    rms2 += (mean-cosThetaDistr[j])*(mean-cosThetaDistr[j]);
260  }
261  rms2 /= jMax-1;
262
263  rms = std::sqrt(rms2);
264
265  relDelta = rms/mean;
266  G4cout << G4endl;
267  G4cout << "meanCosThetaDistr = "<<mean<<" +- "<<rms
268         <<" with relative error = "<<relDelta <<G4endl;
269
270  mean = 0;
271  rms2 = 0.;
272
273  for( j = 0; j < jMax; j++ )
274  {
275    mean += phi[j];
276  }
277  mean /= jMax;
278
279  for( j = 0; j < jMax; j++ )
280  {
281    rms2 += (mean-phi[j])*(mean-phi[j]);
282  }
283  rms2 /= jMax-1;
284
285  rms = std::sqrt(rms2);
286
287  relDelta = rms/mean;
288  G4cout << G4endl;
289  G4cout << "meanPhiDistr = "<<mean<<" +- "<<rms
290         <<" with relative error = "<<relDelta <<G4endl;
291
292  return 1 ;
293}
Note: See TracBrowser for help on using the repository browser.