35 use kinds, only: r_kind, i_kind, r_single
42 integer :: npe, mype, mypelocal,ierror
43 integer,
allocatable :: mpi_layout_begin(:),mpi_layout_end(:)
53 character*180 :: thisfv3file
54 character*2 :: workpath
57 integer :: map_proj, nlon, nlat
58 integer :: fv3lon, fv3lat, fv3times
59 integer :: lbclon, lbclat, lbctimes
60 integer :: i, j, t1, t2
61 integer :: num_args, ix
62 integer :: indexfvcomsel
64 real :: rad2deg = 180.0/3.1415926
65 real :: userdx, userdy, cen_lat, cen_lon
66 real :: usertruelat1, usertruelat2, moad_cen_lat, stand_lon
67 real :: truelat1, truelat2, stdlon, lat1, lon1, r_earth
68 real :: knowni, knownj, dx
69 real :: one, pi, deg2rad
72 character(len=180) :: fv3file
73 character(len=180) :: fvcomfile
74 character(len=180) :: wcstart
75 character(len=180) :: inputfvcomselstr
76 character(len=180),
dimension(:),
allocatable :: args
77 integer :: fv3_io_layout_y
78 integer,
allocatable :: fv3_layout_begin(:),fv3_layout_end(:)
80 real(r_kind),
allocatable :: fv3ice(:,:), fv3sst(:,:)
81 real(r_kind),
allocatable :: fv3sfct(:,:), fv3mask(:,:)
82 real(r_kind),
allocatable :: fv3icet(:,:), fv3sfctl(:,:)
83 real(r_kind),
allocatable :: fv3zorl(:,:), fv3hice(:,:)
84 real(r_kind),
allocatable :: lbcice(:,:), lbcsst(:,:)
85 real(r_kind),
allocatable :: lbcsfct(:,:), lbcmask(:,:)
86 real(r_kind),
allocatable :: lbcicet(:,:), lbczorl(:,:)
87 real(r_kind),
allocatable :: lbchice(:,:)
92 integer :: update_type
93 character(len=4) :: clayout_y
100 call mpi_comm_size(mpi_comm_world,npe,ierror)
101 call mpi_comm_rank(mpi_comm_world,mype,ierror)
107 write(*,*)
'number of cores=',npe
111 num_args = command_argument_count()
112 allocate(args(num_args))
115 call get_command_argument(ix,args(ix))
121 fv3file=trim(args(1))
122 write(*,*)
'surface file=',trim(fv3file)
123 fvcomfile=trim(args(2))
124 write(*,*)
'fvcom file=',trim(fvcomfile)
125 wcstart=trim(args(3))
126 write(*,*)
'warm or cold start = ', trim(wcstart)
127 inputfvcomselstr=trim(args(4))
129 write(*,*)
'select time = ', trim(inputfvcomselstr)
130 clayout_y=trim(args(5))
131 read(clayout_y,
'(I2)') fv3_io_layout_y
132 if (wcstart ==
'cold')
then
135 write(*,*)
'subdomain number=',fv3_io_layout_y
136 if(mype < fv3_io_layout_y)
then
137 write(thisfv3file,
'(a,I3.3)')
'stdout_fvcom.',mype
138 open(6, file=trim(thisfv3file),form=
'formatted',status=
'unknown')
139 write(6,*)
'===> start process for core = ', mype
145 allocate(fv3_layout_begin(fv3_io_layout_y))
146 allocate(fv3_layout_end(fv3_io_layout_y))
147 allocate(mpi_layout_begin(npe))
148 allocate(mpi_layout_end(npe))
154 do ix=1,fv3_io_layout_y
155 if(fv3_io_layout_y > 1)
then
156 write(thisfv3file,
'(a,a,a1,I4.4)') trim(workpath), trim(fv3file),
'.',ix-1
158 write(thisfv3file,
'(a,a)') trim(workpath), trim(fv3file)
160 call geo%open(trim(thisfv3file),
'r',200)
161 call geo%get_dim(
"yaxis_1",nlat)
163 fv3_layout_begin(ix)=iy+1
164 fv3_layout_end(ix)=iy+nlat
165 iy=fv3_layout_end(ix)
167 write(6,
'(a,20I5)')
'begin index for each subdomain',fv3_layout_begin
168 write(6,
'(a,20I5)')
' end index for each subdomain',fv3_layout_end
172 mpi_layout_begin(mypelocal)=1
173 mpi_layout_end(mypelocal)=fv3_io_layout_y
177 do ix=1,fv3_io_layout_y
179 mpi_layout_end(iy)=mpi_layout_end(iy)+1
183 if(mpi_layout_end(ix) > 0)
then
184 mpi_layout_begin(ix)=iy+1
185 mpi_layout_end(ix)=iy+mpi_layout_end(ix)
186 iy=mpi_layout_end(ix)
188 mpi_layout_begin(ix)=0
193 write(6,
'(a)')
'begin and end domain index for each core:'
194 write(6,
'(20I5)') mpi_layout_begin
195 write(6,
'(20I5)') mpi_layout_end
197 if(mypelocal <= fv3_io_layout_y)
then
199 do ix=mpi_layout_begin(mypelocal),mpi_layout_end(mypelocal)
201 write(6,*)
'process subdomain ',ix,
' with core ',mype
202 if(fv3_io_layout_y > 1)
then
203 write(thisfv3file,
'(a,a,a1,I4.4)') trim(workpath), trim(fv3file),
'.',ix-1
205 write(thisfv3file,
'(a,a)') trim(workpath), trim(fv3file)
207 write(6,*)
'sfc data file', trim(thisfv3file)
209 call geo%open(trim(thisfv3file),
'r',200)
210 call geo%get_dim(
"xaxis_1",nlon)
211 call geo%get_dim(
"yaxis_1",nlat)
212 write(6,*)
'NLON,NLAT:', nlon, nlat
215 write(6,*)
'Finished reading sfc_data grid information.'
220 write(6,*)
'the Y dimensions=',fv3_layout_begin(ix),fv3_layout_end(ix)
224 allocate(fv3ice(nlon,nlat))
225 allocate(fv3sfct(nlon,nlat))
226 allocate(fv3sst(nlon,nlat))
227 allocate(fv3mask(nlon,nlat))
228 allocate(fv3icet(nlon,nlat))
229 allocate(fv3sfctl(nlon,nlat))
230 allocate(fv3zorl(nlon,nlat))
231 allocate(fv3hice(nlon,nlat))
233 allocate(lbcice(nlon,nlat))
234 allocate(lbcsfct(nlon,nlat))
235 allocate(lbcsst(nlon,nlat))
236 allocate(lbcmask(nlon,nlat))
237 allocate(lbcicet(nlon,nlat))
238 allocate(lbczorl(nlon,nlat))
239 allocate(lbchice(nlon,nlat))
248 call fcst%initial(
'FV3LAM',wcstart)
249 call fcst%list_initial
250 call fcst%read_n(trim(thisfv3file),
'FV3LAM',wcstart,fv3lon,fv3lat,fv3times,t1,&
251 fv3mask,fv3sst,fv3ice,fv3sfct,fv3icet,fv3sfctl,fv3zorl,fv3hice, &
253 call fcst%finish(
'FV3LAM',wcstart)
256 write(6,*)
'fv3times: ', fv3times
257 write(6,*)
'time to use: ', t1
264 call fcst%initial(
' FVCOM',wcstart)
265 call fcst%list_initial
266 call fcst%get_time_ind(trim(fvcomfile),inputfvcomselstr,indexfvcomsel)
269 write(6,*)
'time asked for =', trim(inputfvcomselstr)
270 write(6,*)
'time index selected = ', t2
271 call fcst%read_n(trim(fvcomfile),
' FVCOM',wcstart,lbclon,lbclat,lbctimes,t2, &
272 lbcmask,lbcsst,lbcice,lbcsfct,lbcicet,fv3sfctl,lbczorl,lbchice, &
273 fv3_layout_begin(ix),fv3_layout_end(ix))
274 call fcst%finish(
' FVCOM',wcstart)
278 if (lbclon .ne. nlon .or. lbclat .ne. nlat)
then
279 write(6,*)
'ERROR: FVCOM/FV3 dimensions do not match:'
280 write(6,*)
'lbclon: ', lbclon
281 write(6,*)
'nlon: ', nlon
282 write(6,*)
'lbclat: ', lbclat
283 write(6,*)
'nlat: ', nlat
287 write(6,*)
'lbctimes: ', lbctimes
288 write(6,*)
'time to use: ', t2
290 if (t2 .gt. lbctimes)
then
291 write(6,*)
'ERROR: Requested time dimension out of range'
292 write(6,*)
'Length of time dimension: ', lbctimes
293 write(6,*)
'Time index to use: ', t2
301 if (wcstart ==
'warm')
then
304 if (lbcmask(i,j) > 0. .and. lbcsst(i,j) .ge. -90.0)
then
306 if (lbcice(i,j) < 0.15)
then
310 fv3ice(i,j) = lbcice(i,j)
311 fv3hice(i,j) = lbchice(i,j)
314 if (lbcice(i,j) > 0. .and. fv3mask(i,j) == 0.)
then
319 if (fv3mask(i,j) == 2. .and. lbcice(i,j) == 0.)
then
321 fv3zorl(i,j) = zero / zero
323 fv3sst(i,j) = lbcsst(i,j) + 273.15
324 fv3sfct(i,j) = lbcsst(i,j) + 273.15
325 fv3icet(i,j) = lbcsst(i,j) + 273.15
326 fv3sfctl(i,j)= lbcsst(i,j) + 273.15
328 if (lbcice(i,j) > 0.)
then
329 fv3icet(i,j) = lbcicet(i,j) + 273.15
334 else if (wcstart ==
'cold')
then
337 if (lbcmask(i,j) > 0. .and. lbcsst(i,j) .ge. -90.0)
then
339 if (lbcice(i,j) < 0.15)
then
343 fv3ice(i,j) = lbcice(i,j)
344 fv3hice(i,j) = lbchice(i,j)
346 if (lbcice(i,j) > 0. .and. fv3mask(i,j) == 0.)
then
351 if (fv3mask(i,j) == 2. .and. lbcice(i,j) == 0.)
then
353 fv3zorl(i,j) = zero / zero
355 fv3sst(i,j) = lbcsst(i,j) + 273.15
356 fv3sfct(i,j) = lbcsst(i,j) + 273.15
357 fv3icet(i,j) = lbcsst(i,j) + 273.15
359 if (lbcice(i,j) > 0.)
then
360 fv3icet(i,j) = lbcicet(i,j) + 273.15
366 write(6,*)
'Variable wcstart is not set to either warm or cold'
371 write(6,*)
"udpate file=",trim(thisfv3file)
372 call geo%open(trim(thisfv3file),
'w',300)
373 call geo%replace_var(
"tsea",nlon,nlat,fv3sst)
374 call geo%replace_var(
"fice",nlon,nlat,fv3ice)
375 call geo%replace_var(
"slmsk",nlon,nlat,fv3mask)
376 call geo%replace_var(
"tisfc",nlon,nlat,fv3icet)
377 call geo%replace_var(
"hice",nlon,nlat,fv3hice)
378 if (wcstart ==
'cold')
then
380 call geo%replace_var(
"zorl",nlon,nlat,fv3zorl)
381 call geo%add_new_var(
'glmsk',
'xaxis_1',
'yaxis_1',
'Time',
'glmsk',
'none')
382 call geo%replace_var(
'glmsk',nlon,nlat,lbcmask)
384 if (wcstart ==
'warm')
then
385 call geo%replace_var(
"zorli",nlon,nlat,fv3zorl)
386 call geo%replace_var(
"tsfc",nlon,nlat,fv3sfct)
387 call geo%replace_var(
"tsfcl",nlon,nlat,fv3sfctl)
388 call geo%add_new_var(
'glmsk',
'xaxis_1',
'yaxis_1',
'Time',
'glmsk',
'none')
389 call geo%replace_var(
'glmsk',nlon,nlat,lbcmask)
410 write(6,*)
"=== LOWBC UPDATE SUCCESS ===",ix
416 call mpi_finalize(ierror)
Module to hold specification kinds for variable declaration.
Functions to read and write netcdf files.
This module defines FV3LAM and FVCOM forecast data structure and the method to read and write observa...
program process_fvcom
Put lake surface temperature and aerial ice concentration from GLERL-provided FVCOM forecast files (w...