ocnice_prep 1.14.0
Loading...
Searching...
No Matches
init_mod.F90
Go to the documentation of this file.
1
7module init_mod
8
9 implicit none
10
11 public
12
13 integer, parameter :: maxvars = 60
14 character(len=10) :: maskvar = 'h'
16 type :: vardefs
17 character(len= 20) :: var_name
18 character(len=120) :: long_name
19 character(len= 20) :: units
20 character(len= 20) :: var_remapmethod
21 integer :: var_dimen
22 character(len= 4) :: var_grid
23 character(len= 20) :: var_pair
24 character(len= 4) :: var_pair_grid
25 end type vardefs
26
27 type(vardefs) :: outvars(maxvars)
29
30 character(len=10) :: ftype
31 character(len=10) :: fsrc
32 character(len=10) :: fdst
33 character(len=120) :: wgtsdir
34 character(len=120) :: griddir
35 character(len=20) :: input_file
36
37 integer :: nxt
38 integer :: nyt
39 integer :: nlevs
40
41 integer :: nxr
42 integer :: nyr
43
44 integer :: logunit
45 logical :: debug
46 logical :: do_ocnprep
47
48contains
56 subroutine readnml(fname,errmsg,rc)
57
58 character(len=*), intent(in) :: fname
59 character(len=*), intent(out) :: errmsg
60 integer, intent(out) :: rc
61
62 ! local variable
63 logical :: fexist
64 integer :: ierr, iounit
65 integer :: srcdims(2), dstdims(2)
66 !----------------------------------------------------------------------------
67
68 namelist /ocniceprep_nml/ ftype, wgtsdir, griddir, srcdims, dstdims, debug
69
70 srcdims = 0; dstdims = 0
71 errmsg='' ! for successful return
72 rc = 0 ! for successful retun
73
74 inquire(file=trim(fname), exist=fexist)
75 if (.not. fexist) then
76 write (errmsg, '(a)') 'FATAL ERROR: input file '//trim(fname)//' does not exist.'
77 rc = 1
78 return
79 else
80 ! Open and read namelist file.
81 open (action='read', file=trim(fname), iostat=ierr, newunit=iounit)
82 read (nml=ocniceprep_nml, iostat=ierr, unit=iounit)
83 if (ierr /= 0) then
84 rc = 1
85 write (errmsg, '(a)') 'FATAL ERROR: invalid namelist format.'
86 return
87 end if
88 close (iounit)
89 end if
90
91 ! append slash to wgtsdir and griddir
92 wgtsdir = trim(wgtsdir)//'/'
93 griddir = trim(griddir)//'/'
94
95 ! check that model is either ocean or ice
96 if (trim(ftype) /= 'ocean' .and. trim(ftype) /= 'ice') then
97 rc = 1
98 write (errmsg, '(a)') 'FATAL ERROR: ftype must be ocean or ice'
99 return
100 end if
101
102 ! set grid dimensions and names
103 nxt = srcdims(1); nyt = srcdims(2)
104 nxr = dstdims(1); nyr = dstdims(2)
105 fsrc = '' ; fdst = ''
106 if (nxt == 1440 .and. nyt == 1080) fsrc = 'mx025' ! 1/4deg tripole
107 if (nxt == 720 .and. nyt == 576) fsrc = 'mx050' ! 1/2deg tripole
108 if (nxt == 360 .and. nyt == 320) fsrc = 'mx100' ! 1deg tripole
109 if (len_trim(fsrc) == 0) then
110 rc = 1
111 write(errmsg,'(a)')'FATAL ERROR: source grid dimensions invalid'
112 return
113 end if
114
115 if (nxr == 1440 .and. nyr == 1080) fdst = 'mx025' ! 1/4deg tripole
116 if (nxr == 720 .and. nyr == 576) fdst = 'mx050' ! 1/2deg tripole
117 if (nxr == 360 .and. nyr == 320) fdst = 'mx100' ! 1deg tripole
118 if (nxr == 72 .and. nyr == 35) fdst = 'mx500' ! 5deg tripole
119 if (len_trim(fdst) == 0) then
120 rc = 1
121 write(errmsg,'(a)')'FATAL ERROR: destination grid dimensions invalid'
122 return
123 end if
124
125 if (trim(fsrc) .eq. trim(fdst)) then
126 rc = 1
127 write(errmsg,'(a)')'FATAL ERROR: Source and destination grids must differ'
128 return
129 end if
130
131 if (trim(fdst) == 'mx500') then
132 rc = 1
133 write(errmsg,'(a)')'FATAL ERROR: 5deg destination grid is not supported'
134 return
135 end if
136
137 ! initialize the source file types
138 if (trim(ftype) == 'ocean') then
139 do_ocnprep = .true.
140 else
141 do_ocnprep = .false.
142 end if
143
144 input_file = trim(ftype)//'.nc'
145 inquire (file=trim(input_file), exist=fexist)
146 if (.not. fexist) then
147 write (errmsg, '(a)') 'FATAL ERROR: input file '//trim(input_file)//' does not exist.'
148 rc=1
149 return
150 end if
151
152 ! log file
153 open(newunit=logunit, file=trim(ftype)//'.prep.log',form='formatted')
154 if (debug) write(logunit, '(a)')'input file: '//trim(input_file)
155
156 ! all checks pass, continue
157 write(errmsg,'(a)') 'Namelist successfully read, continue'
158 rc = 0
159
160 end subroutine readnml
161
170 subroutine readcsv(fname,errmsg,rc,nvalid)
171
172 character(len=*), intent(in) :: fname
173 character(len=*), intent(out) :: errmsg
174 integer, intent(out) :: rc
175 integer, intent(out) :: nvalid
176
177 ! local variables
178 character(len=100) :: chead
179 character(len= 20) :: c1,c3,c4,c5,c6
180 integer :: i2, idx1,idx2
181 integer :: nn,n,ierr,iounit
182 !----------------------------------------------------------------------------
183
184 errmsg='' ! for successful return
185 rc = 0 ! for successful retun
186
187 open(newunit=iounit, file=trim(fname), status='old', iostat=ierr)
188 if (ierr /= 0) then
189 rc = 1
190 write (errmsg, '(a)') 'FATAL ERROR: input file '//trim(fname)//' does not exist.'
191 return
192 end if
193
194 read(iounit,*)chead
195 nn=0
196 do n = 1,maxvars
197 read(iounit,*,iostat=ierr)c1,i2,c3,c4,c5,c6
198 if (ierr .ne. 0) exit
199 if (len_trim(c1) > 0) then
200 nn = nn+1
201 outvars(nn)%var_name = trim(c1)
202 outvars(nn)%var_dimen = i2
203 outvars(nn)%var_grid = trim(c3)
204 outvars(nn)%var_remapmethod = trim(c4)
205 outvars(nn)%var_pair = trim(c5)
206 outvars(nn)%var_pair_grid = trim(c6)
207 end if
208 end do
209 close(iounit)
210 nvalid = nn
211 ! check for u,v pairs, these should be listed in csv file in ordered pairs
212 idx1 = 0; idx2 = 0
213 do n = 1,nvalid
214 if (len_trim(outvars(n)%var_pair) > 0 .and. idx1 .eq. 0) then
215 idx1 = n
216 idx2 = n+1
217 end if
218 end do
219
220 if (idx1*idx2 > 0) then
221 if (trim(outvars(idx1)%var_pair) /= trim(outvars(idx2)%var_name)) then
222 rc = 1
223 write(errmsg,'(a)')'FATAL ERROR: vector pair for '//trim(outvars(idx1)%var_name) &
224 //' is not set correctly'
225 return
226 end if
227 if (trim(outvars(idx2)%var_pair) /= trim(outvars(idx1)%var_name)) then
228 rc = 1
229 write(errmsg,'(a)')'FATAL ERROR: vector pair for '//trim(outvars(idx2)%var_name) &
230 //' is not set correctly'
231 return
232 end if
233
234 ! check for u velocities on u-staggers and v-velocities on v-staggers
235 if (outvars(idx1)%var_name(1:1) == 'u') then
236 if ((outvars(idx1)%var_grid(1:2) /= 'Cu') .and. outvars(idx1)%var_grid(1:2) /= 'Bu') then
237 rc = 1
238 write(errmsg,'(a)')'FATAL ERROR: u-vector has wrong grid '
239 return
240 end if
241 end if
242 if (outvars(idx2)%var_name(1:1) == 'v') then
243 if ((outvars(idx2)%var_grid(1:2) /= 'Cv') .and. outvars(idx2)%var_grid(1:2) /= 'Bu') then
244 rc = 1
245 write(errmsg,'(a)')'FATAL ERROR: v-vector has wrong grid '
246 return
247 end if
248 end if
249 end if
250
251 ! all checks pass, continue
252 write(errmsg,'(a)')'CSV successfully read, continue'
253 rc = 0
254
255 end subroutine readcsv
256end module init_mod