COSMOS v7.655  COSMOSv7655
(AirShowerMC)
ksampAF.f
Go to the documentation of this file.
1 !c teste ksampAF
2 !c (x,y) is wave length vs light yield for Li2B4O7
3 ! integer n
4 ! parameter (n=33)
5 ! real*8 x(n), y(n)
6 !c data x/300., 320., 340., 360., 380., 390., 400.,
7 ! data x/ 320., 340., 360., 380., 390., 400.,
8 ! * 410., 420., 430., 440.,450.,460.,470., 480., 490, 500.,
9 ! * 510., 520., 530., 540.,550.,560.,570., 580., 590, 600.,
10 ! * 610., 620., 630., 640., 650., 660., 670./
11 !c data y/0., 0.1, 1.,4., 15.,20.,32.,
12 ! data y/ 0.0, 1., 4., 15.,20.,32.,
13 ! * 39., 56., 70., 78., 90., 95., 99., 105.,111.,118.,
14 ! * 117., 99., 79., 63.,56., 47., 43., 40.,37.,31.,23.,
15 ! * 20., 17., 16., 10., 5., 0./
16 ! integer nc
17 ! parameter (nc=n+1)
18 ! real*8 coef(nc, 3)
19 !
20 ! real*8 yi(n)
21 ! real*8 coef2(nc,3)
22 ! real*8 a, b, dx, xx, u, yy
23 !
24 ! real*8 xs
25 ! integer i
26 ! call ksampAF0(x, y, n, coef, nc, yi, coef2)
27 ! to see how well interpolation was done
28 ! do i = 1, n
29 ! write(*,'(a, 3g15.4)') 'o ', x(i), y(i), yi(i)
30 ! enddo
31 ! a= x(1)
32 ! b= x(n)
33 ! dx = (b-a)/(5*n)
34 ! do i = 1, 5*n+1
35 ! xx = a + (i-1)*dx
36 ! u = (i-1)/(5.0*n)
37 ! call kcsplIntp(x, y, n, coef, nc, xx, xs)
38 ! call kcsplIntp(yi, x, n, coef2, nc, u, yy)
39 ! write(*,'(a, 4e15.4)') 'i ', xx, xs, u, yy
40 ! enddo
41 !
42 !c generate random number
43 ! do i = 1, 100000
44 ! call ksampAF(x, yi, n, coef2, nc, xs)
45 ! write(*,*) xs
46 ! enddo
47 ! end
48 !
49 ! sample a random number following an arbitrary function
50 ! which is expressed as a table of n-point (x,y)
51 ! This is initialization routine
52 !
53  subroutine ksampaf0(x, y, n, coef, nc, yi, total,
54  * coef2)
55  implicit none
56  integer n ! input number of data points
57  real*8 x(n), y(n) ! input data points (x,y)
58 ! zero is assumed at x<x(1) or x>x(n)
59  integer nc ! input. must be >=n-1
60  real*8 coef(nc, 3) ! output. to keep the spline coefficents
61  real*8 yi(n) ! output. yi(i)= Normalized integral from x=x(1) to x(i)
62 ! (i = 1, n) (yia(1)=0. yi(n)= 1.0)
63  real*8 total ! integral valeu of the function
64 !
65  real*8 coef2(nc, 3) !ouptut. coef2 to be used for (yi, x)
66 ! spline interploation
67 !
68 
69 
70 
71  real*8 dx, a, b, sum
72  integer i
73  real*8 x1, x2
74 
75  if( nc .lt. n-1 ) then
76  write(0,*) ' nc< n-1 in ksamplAF0'
77  stop 1234
78  endif
79 
80  call kcsplcoef(x, y, n, coef, nc)
81  a = x(1)
82  b = x(n)
83  call kcsplinteg(x, y, n, coef, nc, a, b, total)
84 
85  yi(1) = 0.
86  do i = 2, n
87  b = x(i)
88  call kcsplinteg(x, y, n, coef, nc, a, b, sum)
89  yi(i) = sum/total
90  enddo
91 ! to assure yi(n) = 1.
92  yi(n) = 1.
93  call kcsplcoef(yi, x, n, coef2, nc)
94  end
95 
96 
97  subroutine ksampaf0_b(x, y, n, coef, nc)
98  implicit none
99 ! this is only for getting interplation
100  integer n ! input number of data points
101  real*8 x(n), y(n) ! input data points (x,y)
102 ! zero is assumed at x<x(1) or x>x(n)
103  integer nc ! input. must be >=n-1
104  real*8 coef(nc, 3) ! output. to keep the spline coefficents
105 
106 
107  real*8 dx, a, b
108  integer i
109  real*8 x1, x2
110 
111  if( nc .lt. n-1 ) then
112  write(0,*) ' nc< n-1 in ksamplAF0'
113  stop 1234
114  endif
115 
116  call kcsplcoef(x, y, n, coef, nc)
117  end
118 
119  subroutine ksampaf(x, yi, n, coef2, nc, xs)
120  implicit none
121  integer n
122  real*8 x(n), yi(n)
123  integer nc
124  real*8 coef2(nc, 3)
125  real*8 xs !output. sampled x
126  real*8 u
127 
128  call rndc(u)
129  call kcsplintp(yi, x, n, coef2, nc, u, xs)
130  end
subroutine kcsplinteg(x, y, n, coef, nc, a, b, s)
Definition: kcsplInteg.f:2
subroutine kcsplcoef(x, y, n, coef, nc)
Definition: kcsplCoef.f:2
subroutine ksampaf(x, yi, n, coef2, nc, xs)
Definition: ksampAF.f:120
subroutine rndc(u)
Definition: rnd.f:91
subroutine ksampaf0(x, y, n, coef, nc, yi, total, coef2)
Definition: ksampAF.f:55
subroutine ksampaf0_b(x, y, n, coef, nc)
Definition: ksampAF.f:98
subroutine kcsplintp(x, y, n, coef, nc, v, f)
Definition: kcsplIntp.f:2