da_verif_obs_init.f90

References to this file elsewhere.
1 MODULE da_verif_obs_init
2 !----------------------------------------------------------------------------   
3 ! History:
4 !
5 !  Abstract:  
6 !   Main module for 
7 !   initializing various arrays                   
8 !
9 !  Author:   Syed RH Rizvi     NCAR/MMM         05/30/2006
10 !----------------------------------------------------------------------------   
11    USE da_verif_obs_control
12 
13 CONTAINS
14 
15    subroutine da_advance_cymdh( start_date, dh, end_date )
16 
17    implicit none
18 
19    character (len=10), intent(in)  :: start_date ! In date (ccyymmddhh).
20    integer, intent(in)             :: dh         ! Period to advance (-ve for past).
21    character (len=10), intent(out) :: end_date   ! Out date (ccyymmddhh).
22 
23    integer :: ccyy, mm, dd, hh
24 
25    read(start_date(1:10), fmt='(i4, 3i2)')  ccyy, mm, dd, hh
26 
27    hh = hh + dh
28 
29    do while (hh < 0)
30       hh = hh + 24
31       call da_change_date ( ccyy, mm, dd, -1 )
32    end do
33 
34    do while (hh > 23)
35       hh = hh - 24
36       call da_change_date ( ccyy, mm, dd, 1 )
37    end do
38 
39    write(UNIT=end_date(1:10), fmt='(i4, 3i2.2)')  ccyy, mm, dd, hh
40 
41 end subroutine da_advance_cymdh
42 subroutine da_change_date( ccyy, mm, dd, delta )
43 
44    implicit none
45 
46    integer, intent(inout) :: ccyy, mm, dd
47    integer, intent(in)    :: delta
48 
49    integer, dimension(12) :: mmday
50 
51    mmday = (/31,28,31,30,31,30,31,31,30,31,30,31/)
52 
53    mmday(2) = 28
54    if (mod(ccyy,4) == 0) then
55       mmday(2) = 29
56 
57       if ( mod(ccyy,100) == 0) then
58          mmday(2) = 28
59       endif
60 
61       if(mod(ccyy,400) == 0) then
62          mmday(2) = 29
63       end if
64    endif
65 
66    dd = dd + delta
67 
68    if(dd == 0) then
69       mm = mm - 1
70 
71       if(mm == 0) then
72          mm = 12
73          ccyy = ccyy - 1
74       endif
75 
76       dd = mmday(mm)
77    elseif ( dd .gt. mmday(mm) ) then
78       dd = 1
79       mm = mm + 1
80       if(mm > 12 ) then
81          mm = 1
82          ccyy = ccyy + 1
83       end if
84    end if
85 end subroutine da_change_date
86    subroutine initialize_surface_type(surface)
87    type(surface_type), intent(inout)    :: surface
88 !
89    call initialize_stats_type(surface%uomb,surface%uoma)
90    call initialize_stats_type(surface%vomb,surface%voma)
91    call initialize_stats_type(surface%tomb,surface%toma)
92    call initialize_stats_type(surface%pomb,surface%poma)
93    call initialize_stats_type(surface%qomb,surface%qoma)
94 
95    end subroutine initialize_surface_type
96 !
97    subroutine initialize_upr_type(upr)
98    type(upr_type), intent(inout)    :: upr
99 !
100    integer    :: k
101 !
102    do k = 1, nstd
103    call initialize_stats_type(upr%uomb(k),upr%uoma(k))
104    call initialize_stats_type(upr%vomb(k),upr%voma(k))
105    call initialize_stats_type(upr%tomb(k),upr%toma(k))
106    call initialize_stats_type(upr%qomb(k),upr%qoma(k))
107    enddo
108 
109    end subroutine initialize_upr_type
110 !
111    subroutine initialize_gpspw_type(gpspw)
112    type(gpspw_type), intent(inout)    ::  gpspw
113 !
114       call initialize_stats_type(gpspw%tpwomb,gpspw%tpwoma)
115    end subroutine initialize_gpspw_type
116 !  
117    subroutine initialize_gpsref_type(gpsref)
118    type(gpsref_type), intent(inout)    ::  gpsref
119    integer     :: k
120 
121 !  
122    do k = 1, nstdh
123    call initialize_stats_type(gpsref%refomb(k),gpsref%refoma(k))
124    enddo
125    end subroutine initialize_gpsref_type
126 !  
127    subroutine initialize_stats_type(omb, oma)
128    type(stats_value), intent(inout)    :: omb, oma
129    omb%num   = 0 ; oma%num   = 0
130    omb%bias  = 0 ; oma%bias  = 0
131    omb%abias = 0 ; oma%abias = 0
132    omb%rmse  = 0 ; oma%rmse  = 0
133    end subroutine initialize_stats_type
134 !
135    subroutine initialize_t_tab
136    implicit none
137    integer :: i
138 
139  !
140   ! Initalize Student t table for alpha=0.025
141   !
142 
143   ! Degrees of freedom
144    do i=1,30
145     alpha(i,1) = i
146    end do
147    alpha(31,1) = 40.
148    alpha(32,1) = 60.
149    alpha(33,1) = 120.
150    alpha(34,1) = 400.
151  
152    ! Critical t value
153    alpha(1,2) = 12.706
154    alpha(2,2) = 4.303
155    alpha(3,2) = 3.182
156    alpha(4,2) = 2.776
157    alpha(5,2) = 2.571
158    alpha(6,2) = 2.447
159    alpha(7,2) = 2.365
160    alpha(8,2) = 2.306
161    alpha(9,2) = 2.262
162    alpha(10,2) = 2.228
163    alpha(11,2) = 2.201
164    alpha(12,2) = 2.179
165    alpha(13,2) = 2.160
166    alpha(14,2) = 2.145
167    alpha(15,2) = 2.131
168    alpha(16,2) = 2.120
169    alpha(17,2) = 2.110
170    alpha(18,2) = 2.101
171    alpha(19,2) = 2.093
172    alpha(20,2) = 2.086
173    alpha(21,2) = 2.080
174    alpha(22,2) = 2.074
175    alpha(23,2) = 2.069
176    alpha(24,2) = 2.064
177    alpha(25,2) = 2.060
178    alpha(26,2) = 2.056
179    alpha(27,2) = 2.052
180    alpha(28,2) = 2.048
181    alpha(29,2) = 2.045
182    alpha(30,2) = 2.042
183    alpha(31,2) = 2.021
184    alpha(32,2) = 2.000
185    alpha(33,2) = 1.980
186    alpha(34,2) = 1.960
187 
188    end subroutine initialize_t_tab
189 
190 !-------------------------------------------------
191 end MODULE da_verif_obs_init