fvcom_tools 1.14.0
Loading...
Searching...
No Matches
process_FVCOM.f90
Go to the documentation of this file.
1
6
33
34 use mpi
35 use kinds, only: r_kind, i_kind, r_single
36 use module_ncio, only: ncio
37 use module_nwp, only: fcst_nwp
38
39 implicit none
40
41! MPI variables
42 integer :: npe, mype, mypelocal,ierror
43 integer,allocatable :: mpi_layout_begin(:),mpi_layout_end(:)
44!
45
46! New object-oriented declarations
47
48 type(ncio) :: geo
49 type(fcst_nwp) :: fcst
50
51! Grid variables
52
53 character*180 :: thisfv3file
54 character*2 :: workpath
55 character*1 :: char1
56
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
63
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
70 real :: zero
71
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(:)
79
80 real(r_single), allocatable :: fv3ice(:,:), fv3sst(:,:)
81 real(r_single), allocatable :: fv3sfct(:,:), fv3mask(:,:)
82 real(r_single), allocatable :: fv3icet(:,:), fv3sfctl(:,:)
83 real(r_single), allocatable :: fv3zorl(:,:), fv3hice(:,:)
84 real(r_single), allocatable :: lbcice(:,:), lbcsst(:,:)
85 real(r_single), allocatable :: lbcsfct(:,:), lbcmask(:,:)
86 real(r_single), allocatable :: lbcicet(:,:), lbczorl(:,:)
87 real(r_single), allocatable :: lbchice(:,:)
88
89! Declare namelists
90! SETUP (general control namelist) :
91
92 integer :: update_type
93 character(len=4) :: clayout_y
94 integer :: iy,iblock
95
96! namelist/setup/update_type, t2
97
98! MPI setup
99 call mpi_init(ierror)
100 call mpi_comm_size(mpi_comm_world,npe,ierror)
101 call mpi_comm_rank(mpi_comm_world,mype,ierror)
102!
103! NCEP LSF has to use all cores allocated to run this application
104! but this if check can make sure only one core run through the real code.
105!
106
107 write(*,*) 'number of cores=',npe
108
109 zero = 0.0
110! Get file names from command line arguements
111 num_args = command_argument_count()
112 allocate(args(num_args))
113
114 do ix = 1, num_args
115 call get_command_argument(ix,args(ix))
116 end do
117! fv3file: location of UFS grid
118! fvcomfile: location of FVCOM data file
119! wcstart: warm (restart) or cold start
120! inputFVCOMtimes: string of time to use
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))
128! read(inputFVCOMselStr,*) inputFVCOMsel
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
133 fv3_io_layout_y=1
134 endif
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
140 endif
141!
142! figure out subdomain grid begin and end grid index in Y for each core
143! The match the subdomains to each core
144!
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))
149 fv3_layout_begin=0
150 fv3_layout_end=0
151
152 workpath='./'
153 iy=0
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
157 else
158 write(thisfv3file,'(a,a)') trim(workpath), trim(fv3file)
159 endif
160 call geo%open(trim(thisfv3file),'r',200)
161 call geo%get_dim("yaxis_1",nlat)
162 call geo%close
163 fv3_layout_begin(ix)=iy+1
164 fv3_layout_end(ix)=iy+nlat
165 iy=fv3_layout_end(ix)
166 enddo ! find dimension
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
169
170 mypelocal=mype+1
171 if(npe==1) then
172 mpi_layout_begin(mypelocal)=1
173 mpi_layout_end(mypelocal)=fv3_io_layout_y
174 else
175 mpi_layout_begin=0
176 mpi_layout_end=0
177 do ix=1,fv3_io_layout_y
178 iy=mod(ix-1,npe)+1
179 mpi_layout_end(iy)=mpi_layout_end(iy)+1
180 enddo
181 iy=0
182 do ix=1,npe
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)
187 else
188 mpi_layout_begin(ix)=0
189 mpi_layout_end(ix)=0
190 endif
191 enddo
192 endif
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
196
197if(mypelocal <= fv3_io_layout_y) then
198
199 do ix=mpi_layout_begin(mypelocal),mpi_layout_end(mypelocal)
200
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
204 else
205 write(thisfv3file,'(a,a)') trim(workpath), trim(fv3file)
206 endif
207 write(6,*) 'sfc data file', trim(thisfv3file)
208
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
213 call geo%close
214
215 write(6,*) 'Finished reading sfc_data grid information.'
216 write(6,*) ' '
217!
218! figure out subdomain grid dimension in Y
219!
220 write(6,*) 'the Y dimensions=',fv3_layout_begin(ix),fv3_layout_end(ix)
221
222! Allocate variables for I/O
223
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))
232
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))
240! Read fv3 sfc_data.nc before update
241
242! fv3file='sfc_data.nc'
243! fv3times: length of time dimension of UFS atmospheric grid (should be 1)
244! t1: index of time dimension to pull (should be 1)
245 fv3times=1
246 t1=1
247
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, &
252 1,nlat)
253 call fcst%finish('FV3LAM',wcstart)
254
255
256 write(6,*) 'fv3times: ', fv3times
257 write(6,*) 'time to use: ', t1
258
259! Read FVCOM input datasets
260
261! fvcomfile='fvcom.nc'
262! lbctimes: length of time dimension of FVCOM input data (command line input)
263! Space infront of ' FVCOM' below is important!!
264 call fcst%initial(' FVCOM',wcstart)
265 call fcst%list_initial
266 call fcst%get_time_ind(trim(fvcomfile),inputfvcomselstr,indexfvcomsel)
267! t2: index of time dimension to pull from FVCOM
268 t2=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)
275
276! Check that the dimensions match
277
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
284 stop 'error'
285 endif
286
287 write(6,*) 'lbctimes: ', lbctimes
288 write(6,*) 'time to use: ', t2
289
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
294 stop 'error'
295 endif
296
297! Update with FVCOM fields and process
298! ice cover data. Ice fraction is set
299! to a minimum of 15% due to FV3-LAM
300! raising any value below 15% to 15%.
301 if (wcstart == 'warm') then
302 do j=1,nlat
303 do i=1,nlon
304 if (lbcmask(i,j) > 0. .and. lbcsst(i,j) .ge. -90.0) then !GL Points
305 !If ice fraction below 15%, set to 0
306 if (lbcice(i,j) < 0.15) then
307 lbcice(i,j) = 0.0
308 lbchice(i,j) = 0.0 !remove ice thickness
309 endif
310 fv3ice(i,j) = lbcice(i,j)
311 fv3hice(i,j) = lbchice(i,j)
312
313 !If ice in FVCOM, but not in FV3-LAM, change to ice
314 if (lbcice(i,j) > 0. .and. fv3mask(i,j) == 0.) then
315 fv3mask(i,j) = 2.
316 fv3zorl(i,j) = 1.1
317 endif
318 !If ice in FV3-LAM and not FVCOM, remove it from FV3-LAM
319 if (fv3mask(i,j) == 2. .and. lbcice(i,j) == 0.) then
320 fv3mask(i,j) = 0.
321 fv3zorl(i,j) = zero / zero !Use Fill_Value
322 endif
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
327 !If ice exists in FVCOM, change ice surface temp
328 if (lbcice(i,j) > 0.) then
329 fv3icet(i,j) = lbcicet(i,j) + 273.15
330 end if
331 end if
332 enddo
333 enddo
334 else if (wcstart == 'cold') then
335 do j=1,nlat
336 do i=1,nlon
337 if (lbcmask(i,j) > 0. .and. lbcsst(i,j) .ge. -90.0) then
338 !If ice fraction below 15%, set to 0
339 if (lbcice(i,j) < 0.15) then
340 lbcice(i,j) = 0.0
341 lbchice(i,j) = 0.0 !remove ice thickness
342 endif
343 fv3ice(i,j) = lbcice(i,j)
344 fv3hice(i,j) = lbchice(i,j)
345 !If ice in FVCOM, but not in FV3-LAM, change to ice
346 if (lbcice(i,j) > 0. .and. fv3mask(i,j) == 0.) then
347 fv3mask(i,j) = 2.
348 fv3zorl(i,j) = 1.1
349 endif
350 !If ice in FV3-LAM and not FVCOM, remove it from FV3-LAM
351 if (fv3mask(i,j) == 2. .and. lbcice(i,j) == 0.) then
352 fv3mask(i,j) = 0.
353 fv3zorl(i,j) = zero / zero !Use Fill_Value
354 endif
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
358 !If ice exists in FVCOM, change ice surface temp
359 if (lbcice(i,j) > 0.) then
360 fv3icet(i,j) = lbcicet(i,j) + 273.15
361 end if
362 end if
363 enddo
364 enddo
365 else
366 write(6,*) 'Variable wcstart is not set to either warm or cold'
367 end if
368
369! Write out sfc file again
370
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
379! Add_New_Var takes names of (Variable,Dim1,Dim2,Dim3,Long_Name,Units)
380 call geo%replace_var("zorl",nlon,nlat,fv3zorl)
381 call geo%add_new_var('glmsk','xaxis_1','yaxis_1','Time','glmsk','none','float')
382 call geo%replace_var('glmsk',nlon,nlat,lbcmask)
383 end if
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','float')
389 call geo%replace_var('glmsk',nlon,nlat,lbcmask)
390 end if
391 call geo%close
392
393 deallocate(fv3ice)
394 deallocate(fv3sfct)
395 deallocate(fv3sst)
396 deallocate(fv3mask)
397 deallocate(fv3icet)
398 deallocate(fv3sfctl)
399 deallocate(fv3zorl)
400 deallocate(fv3hice)
401
402 deallocate(lbcice)
403 deallocate(lbcsfct)
404 deallocate(lbcsst)
405 deallocate(lbcmask)
406 deallocate(lbcicet)
407 deallocate(lbczorl)
408 deallocate(lbchice)
409
410 write(6,*) "=== LOWBC UPDATE SUCCESS ===",ix
411
412 enddo ! loop through subdomain
413
414endif ! mype==0
415
416call mpi_finalize(ierror)
417
418
419end program process_fvcom
Module to hold specification kinds for variable declaration.
Definition kinds.f90:11
integer, parameter, public r_single
specification kind for single precision (4-byte) real variable.
Definition kinds.f90:25
integer, parameter, public i_kind
generic specification kind for default integer.
Definition kinds.f90:21
integer, parameter, public r_kind
generic specification kind for default floating point
Definition kinds.f90:26
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...