UPP (develop)
Loading...
Searching...
No Matches
IFI.F
Go to the documentation of this file.
1
3module upp_ifi_mod
4
5#ifdef USE_IFI
6 ! Use the libIFI API modules.
7 use ifi_type_mod
8 use ifi_mod
9#endif
10
11 use iso_c_binding, only: c_bool, c_int64_t, c_ptr, c_double, c_float
12 implicit none
13
14 private
15
16 public run_ifi, set_ifi_dims, ifi_real_t, write_ifi_debug_files, &
17 first_supported_ifi_hour, last_supported_ifi_hour
18
19 logical :: write_ifi_debug_files = .false.
20 real, parameter :: first_supported_ifi_hour = 1
21 real, parameter :: last_supported_ifi_hour = 21
22
23#ifndef USE_IFI
24 ! Stubs need their own version of the ifi_real_t since it is defined in ifi_type_mod
25 integer, parameter :: ifi_real_t = c_float
26
27#else
28 ! Actual IFI code needs a few more module-scope variables.
29
30 real, parameter :: feet2meters = 0.3048
31
32 type(ificonfig) :: ifi_config
33 logical :: have_read_ifi_config = .false.
34
35 ! Communication-related variables for gathers and
36 ! scatters. Gathering is done in the i direction first, all to one
37 ! rank in each row. Then that rank in each row gathers to the
38 ! destination rank. Scatter follows the same process, in reverse.
39
40 integer :: ifi_mpi_real_kind=-999
41 integer :: row_comm = -999 ! communicator for i-direction messages
42 integer :: col_comm = -999 ! communicator for j-direction messages
43 integer :: row_comm_rank=-999, col_comm_rank=-999 ! Rank in the i & j comm communicators
44 integer :: row_comm_size=-999, col_comm_size=-999 ! Size of the i & j comm communicators
45 integer :: grid_rank ! rank within mpi_comm_comp
46
47 integer, allocatable :: ista_grid(:) ! ista for each rank in mpi_comm_comp in i direction
48 integer, allocatable :: jsta_grid(:) ! jsta for each rank in mpi_comm_comp in j direction
49
50 real, allocatable :: rearrange(:,:) ! for rearranging data on local processor
51 integer, allocatable :: row_comm_count(:) ! scatterv&gatherv receive counts for i-direction collectives
52 integer, allocatable :: row_comm_displ(:) ! scatterv&gatherv displacements for i-direction collectives
53 integer, allocatable :: row_ista(:) ! i starts for rearranging data (icomm_size)
54 integer, allocatable :: row_iend(:) ! i ends for rearranging data (icomm_size)
55 real, allocatable :: rearrange_row1d(:) ! for rearranging data for an entire row (roots of icomm communicators)
56
57 real, allocatable :: rearrange_row2d(:,:) ! for rearranging data for an entire row (roots of icomm communicators)
58 integer, allocatable :: col_comm_count(:) ! scatterv&gatherv receive counts for j-direction collectives
59 integer, allocatable :: col_comm_displ(:) ! scatterv&gatherv displacements for j-direction collectives
60 integer, allocatable :: col_jsta(:) ! j starts for rearranging data (jcomm_size)
61 integer, allocatable :: col_jend(:) ! j ends for rearranging data (jcomm_size)
62
63 integer, allocatable :: isize_grid(:) ! iend-ista+1 for each gridpoint
64#endif
65
66contains !==============================================================
67
68 subroutine send_missing_data(ient)
69 ! Fills IFI fields in index ient (of IGET) with missing values.
70 !
71 ! This subroutine is both in the stub and actual IFI code.
72 ! The stubs always send missing data, while the actual code
73 ! sends missing data for unsupported hours.
74 use ctlblk_mod, only: ifi_nflight, ifi_flight_levels
75 use ctlblk_mod, only: spval, ista, iend, jsta, jend, lm, im, cfld, datapd, fld_info, ifi_flight_levels, jsta_2l, jend_2u,me
76 use rqstfld_mod, only: iget, iavblfld, lvlsxml, lvls
77 implicit none
78
79 integer, intent(in) :: ient
80
81 ! Locals
82 integer :: i,j,k
83
84 logical, save :: wrote_message = .false. ! guarded by an OMP CRITICAL below
85
86 if(size(iget)<ient) then
87 return
88 endif
89
90 if(iget(ient)<1 .or. iget(ient)>size(lvls,2)) then
91 return
92 endif
93
94 ! Go level-by-level writing grib2 if requested
95 do k=1,ifi_nflight
96 if(k>size(lvls,1)) then
97 write(0,*) 'ERROR: LVLS is too small to handle ',ifi_nflight,' vertical levels'
98 return
99 endif
100 if(lvls(k,iget(ient))>0) then
101 if(.not.wrote_message) then
102 !$OMP CRITICAL
103 if(.not.wrote_message.and.me==0) then
104 write(*,*) 'This post cannot produce IFI icing products because it was not compiled with libIFI.'
105 wrote_message = .true.
106 endif
107 !$OMP END CRITICAL
108 endif
109
110 cfld = cfld+1
111 fld_info(cfld)%ifld = iavblfld(iget(ient))
112 fld_info(cfld)%lvl = k
113
114 !$OMP PARALLEL DO PRIVATE(i,j) COLLAPSE(2)
115 do j=jsta,jend
116 do i=ista,iend
117 datapd(i-ista+1,j-jsta+1,cfld) = spval
118 enddo
119 enddo
120 endif
121 enddo
122 end subroutine send_missing_data
123
124#ifndef USE_IFI
125! IFI stubs
127 subroutine set_ifi_dims()
128 use ctlblk_mod, only: ifi_nflight, ifi_flight_levels
129 implicit none
130 integer :: i
131
132 ! Bogus fill value for flight levels to prevent a crash
133 ifi_nflight = 60
134 if(allocated(ifi_flight_levels)) then
135 deallocate(ifi_flight_levels)
136 endif
137 allocate(ifi_flight_levels(ifi_nflight))
138 do i=1,ifi_nflight
139 ifi_flight_levels(i) = 500*i
140 enddo
141 end subroutine set_ifi_dims
142
143 ! --------------------------------------------------------------------
145 subroutine run_ifi()
146 call send_missing_data(1007) ! ICE_PROB = missing
147 call send_missing_data(1008) ! SLD = missing
148 call send_missing_data(1009) ! ICE_SEV_CAT = missing
149 call send_missing_data(1010) ! ICE_SEV_CAT = missing
150 end subroutine run_ifi
151
152!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
153
154#else
155
156!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
157
158! Actual IFI code
159
160 subroutine make_communicators
161 use ctlblk_mod, only: jsta, jend, ista, iend, mpi_comm_comp, im, jm
162 use mpi
163 use iso_c_binding, only: c_sizeof
164 implicit none
165
166 integer, allocatable :: key(:) ! ensure ranks are in a consistent order
167 integer :: size, ierr, local_count, accum, irank, i
168 real(kind=ifi_real_t) :: testreal
169
170 if(allocated(jsta_grid)) then
171 return ! Already made communicators
172 endif
173
174 ! Use the right MPI type for the ifi real kind
175 if(c_sizeof(testreal)==4) then
176 ifi_mpi_real_kind=mpi_real4
177 else
178 ifi_mpi_real_kind=mpi_real8
179 endif
180
181 ! Get the correct size of the communicator before we try to split it.
182 call mpi_comm_size(mpi_comm_comp,size,ierr)
183 call mpi_comm_rank(mpi_comm_comp,grid_rank,ierr)
184
185! 929 format('Rank ',I0,' ista=',I0,' jsta=',I0)
186! print 929,grid_rank,ista,jsta
187
188 if(allocated(ista_grid)) then
189 deallocate(ista_grid)
190 deallocate(jsta_grid)
191 deallocate(rearrange)
192 deallocate(rearrange_row1d)
193 deallocate(rearrange_row2d)
194 deallocate(row_ista)
195 deallocate(row_iend)
196 deallocate(row_comm_count)
197 deallocate(row_comm_displ)
198 deallocate(col_jsta)
199 deallocate(col_jend)
200 deallocate(col_comm_count)
201 deallocate(col_comm_displ)
202 endif
203
204 ! Get the start locations on every rank. We need this for the key,
205 ! and for determining root ranks.
206 allocate(ista_grid(size))
207 allocate(jsta_grid(size))
208 ista_grid=-1
209 jsta_grid=-1
210 call mpi_allgather(ista,1,mpi_integer,ista_grid,1,mpi_integer,mpi_comm_comp,ierr)
211 call mpi_allgather(jsta,1,mpi_integer,jsta_grid,1,mpi_integer,mpi_comm_comp,ierr)
212
213 ! Ensure ranks are arranged in a consistent order.
214 allocate(key(size))
215 do i=1,size
216 if(ista_grid(i)==-1) then
217 write(0,*) 'invalid ista_grid ',ista_grid(i)
218 call mpi_abort(mpi_comm_world,1,ierr)
219 endif
220 if(jsta_grid(i)==-1) then
221 write(0,*) 'invalid jsta_grid ',jsta_grid(i)
222 call mpi_abort(mpi_comm_world,1,ierr)
223 endif
224 key(i) = ista_grid(i) + im*(jsta_grid(i)-1)
225 enddo
226
227 ! Create a communicator for row gather/scatter (all ranks that share a jsta)
228 if(all(jsta_grid==jsta_grid(1))) then
229 row_comm = mpi_comm_comp
230 else
231 call mpi_comm_split(mpi_comm_comp,jsta_grid(grid_rank+1),key(grid_rank+1),row_comm,ierr)
232 endif
233 if(row_comm==mpi_comm_null .or. row_comm==0) then
234 write(0,*) 'MPI_Comm_split gave MPI_COMM_NULL for row_comm'
235 call mpi_abort(mpi_comm_world,1,ierr)
236 endif
237 call mpi_comm_rank(row_comm,row_comm_rank,ierr)
238 call mpi_comm_size(row_comm,row_comm_size,ierr)
239
240 ! Create a communicator for column gather/scatter (all ranks that share an ista)
241 if(all(ista_grid==ista_grid(1))) then
242 col_comm = mpi_comm_comp
243 else
244 call mpi_comm_split(mpi_comm_comp,ista_grid(grid_rank+1),key(grid_rank+1),col_comm,ierr)
245 endif
246 if(col_comm==mpi_comm_null .or. col_comm==0) then
247 write(0,*) 'MPI_Comm_split gave MPI_COMM_NULL for col_comm'
248 call mpi_abort(mpi_comm_world,1,ierr)
249 endif
250 call mpi_comm_rank(col_comm,col_comm_rank,ierr)
251 call mpi_comm_size(col_comm,col_comm_size,ierr)
252
253 ! Done with the key.
254 deallocate(key)
255
256 ! Allocate the arrays used to rearrange data.
257 allocate(rearrange(ista:iend,jsta:jend))
258 allocate(rearrange_row1d(im*(jend-jsta+1)))
259 allocate(rearrange_row2d(1:im,jsta:jend))
260
261 ! Calculate information for row MPI_Gatherv and MPI_Scatterv calls.
262 allocate(row_ista(row_comm_size))
263 call mpi_allgather(ista,1,mpi_integer,row_ista,1,mpi_integer,row_comm,ierr)
264 allocate(row_iend(row_comm_size))
265 call mpi_allgather(iend,1,mpi_integer,row_iend,1,mpi_integer,row_comm,ierr)
266 allocate(row_comm_count(row_comm_size))
267 allocate(row_comm_displ(row_comm_size))
268 call get_count_and_displ(row_comm,row_comm_size,(iend-ista+1)*(jend-jsta+1), &
269 row_comm_count,row_comm_displ)
270
271 ! Calculate information for j-direction MPI_Gatherv and MPI_Scatterv calls.
272 allocate(col_jsta(col_comm_size))
273 call mpi_allgather(jsta,1,mpi_integer,col_jsta,1,mpi_integer,col_comm,ierr)
274 allocate(col_jend(col_comm_size))
275 call mpi_allgather(jend,1,mpi_integer,col_jend,1,mpi_integer,col_comm,ierr)
276 allocate(col_comm_count(col_comm_size))
277 allocate(col_comm_displ(col_comm_size))
278 call get_count_and_displ(col_comm,col_comm_size,im*(jend-jsta+1), &
279 col_comm_count,col_comm_displ)
280
281 end subroutine make_communicators
282
283 ! --------------------------------------------------------------------
284
285 subroutine get_count_and_displ(comm,comm_size,here,counts,displacements)
286 use mpi
287 implicit none
288 integer, intent(in) :: comm,comm_size,here
289 integer, intent(inout) :: counts(comm_size),displacements(comm_size)
290 integer :: accum, ierr, i
291
292 ! Get counts from all ranks
293 call mpi_allgather(here,1,mpi_integer,counts,1,mpi_integer,comm,ierr)
294
295 ! Displacements are a cumulative sum starting at 0 for rank 0
296 accum=0
297 do i=1,comm_size
298 displacements(i) = accum
299 accum = accum+counts(i)
300 end do
301 end subroutine get_count_and_displ
302
303 ! --------------------------------------------------------------------
304
305 subroutine find_comm_roots(global_rank,row_root,col_root)
306 ! Find the roots of the row and column communicators for a gather
307 ! or scatter with the specified rank as the location of the global
308 ! array.
309 use mpi
310 implicit none
311
312 integer, intent(in) :: global_rank
313 integer, intent(inout) :: row_root,col_root
314 integer :: r, destination_ista, destination_jsta, ierr
315
316 ! print *,'find roots for rank ',global_rank
317
318 destination_ista=ista_grid(global_rank+1)
319 destination_jsta=jsta_grid(global_rank+1)
320
321 row_root=-1
322 do r=1,row_comm_size
323 if(row_ista(r)==destination_ista) then
324 row_root=r-1
325 exit
326 endif
327 enddo
328
329 col_root=-1
330 do r=1,col_comm_size
331 if(col_jsta(r)==destination_jsta) then
332 col_root=r-1
333 exit
334 endif
335 enddo
336
337! 201 format('root is row_root=',I0,' col_root=',I0,' for rank ',I0)
338! print 201,row_root,col_root,global_rank
339
340 if(col_root==-1) then
341 write(0,'(A,I0)') 'ABORT: Could not find column root rank for rank ',global_rank
342 call mpi_abort(mpi_comm_world,1,ierr)
343 endif
344 if(row_root==-1) then
345 write(0,'(A,I0)') 'ABORT: Could not find row root rank for rank ',global_rank
346 call mpi_abort(mpi_comm_world,1,ierr)
347 endif
348 end subroutine find_comm_roots
349
350 ! --------------------------------------------------------------------
351
352 subroutine global_gather(local,global,global_rank,ista_local,iend_local,jsta_local,jend_local)
353 ! Gather from local arrays on all ranks to global array on specified rank
354 use mpi
355 use ctlblk_mod, only: jsta, jend, ista, iend, im, jm
356 use iso_c_binding, only: c_sizeof
357 implicit none
358
359 real(kind=ifi_real_t), intent(in) :: local(ista_local:iend_local,jsta_local:jend_local)
360 real(kind=ifi_real_t), intent(out) :: global(im,jm)
361 integer, intent(in) :: global_rank,ista_local,iend_local,jsta_local,jend_local
362
363 integer :: i,j,r,inindex,row_root,col_root, destination_ista, destination_jsta,ni_rank,nj_rank,idxstart, ierr
364
365 if(col_comm==0 .or. col_comm==mpi_comm_null) then
366 write(0,*) 'somehow, col_comm became invalid ',col_comm
367 call mpi_abort(mpi_comm_world,1,ierr)
368 endif
369
370 if(row_comm==0 .or. row_comm==mpi_comm_null) then
371 write(0,*) 'somehow, row_comm became invalid ',row_comm
372 call mpi_abort(mpi_comm_world,1,ierr)
373 endif
374
375 ! Find out who is the root of the i and j direction communications.
376 call find_comm_roots(global_rank,row_root,col_root)
377
378 ! Store the data in a contiguous local array
379 !$OMP PARALLEL DO COLLAPSE(2) PRIVATE(i,j)
380 do j=jsta,jend
381 do i=ista,iend
382 rearrange(i,j) = local(i,j)
383 enddo
384 enddo
385
386 ! Gather across all ranks in row
387 call mpi_gatherv(rearrange,(jend-jsta+1)*(iend-ista+1),ifi_mpi_real_kind, &
388 rearrange_row1d,row_comm_count,row_comm_displ,ifi_mpi_real_kind,row_root,row_comm,ierr)
389
390 if(row_comm_rank==row_root) then
391 ! Rearrange from i-j-rank dimensions to i-j dimensions.
392 idxstart=1
393 do r=1,row_comm_size
394
395 ni_rank=row_iend(r)-row_ista(r)+1
396 nj_rank=jend-jsta+1
397 ! print *,'r,ista,iend=',r,row_ista(r),row_iend(r)
398
399 !$OMP PARALLEL DO COLLAPSE(2) PRIVATE(inindex)
400 do j=0,nj_rank-1
401 do i=0,ni_rank-1
402 inindex = idxstart + i + ni_rank*j
403 rearrange_row2d(i+row_ista(r),j+jsta) = rearrange_row1d(inindex)
404 enddo
405 enddo
406
407 idxstart = idxstart + ni_rank*nj_rank
408 enddo
409
410 ! Gather along columns
411 ! print *,'global gatherv to rank ',grid_rank
412 call mpi_gatherv(rearrange_row2d,im*(jend-jsta+1),ifi_mpi_real_kind, &
413 global,col_comm_count,col_comm_displ,ifi_mpi_real_kind,col_root,col_comm,ierr)
414 endif
415
416 ! print *,'global_gather success'
417 end subroutine global_gather
418
419 ! --------------------------------------------------------------------
420
421 subroutine global_scatter(local,global,global_rank,ista_local,iend_local,jsta_local,jend_local)
422 ! Scatter from global array at specified rank to local arrays on all ranks.
423 ! NOTE: Does NOT update halo regions!
424 use mpi
425 use ctlblk_mod, only: jsta, jend, ista, iend, im, jm
426 use iso_c_binding, only: c_sizeof
427 implicit none
428
429 real(kind=ifi_real_t), intent(out) :: local(ista_local:iend_local,jsta_local:jend_local)
430 real(kind=ifi_real_t), intent(in) :: global(im,jm)
431 integer, intent(in) :: global_rank,ista_local,iend_local,jsta_local,jend_local
432
433 integer :: i,j,r,outindex,idxstart,ni_rank,nj_rank,col_root,ierr,row_root
434
435 if(col_comm==0 .or. col_comm==mpi_comm_null) then
436 write(0,*) 'somehow, col_comm became invalid ',col_comm
437 call mpi_abort(mpi_comm_world,1,ierr)
438 endif
439
440 if(row_comm==0 .or. row_comm==mpi_comm_null) then
441 write(0,*) 'somehow, row_comm became invalid ',row_comm
442 call mpi_abort(mpi_comm_world,1,ierr)
443 endif
444
445 ! Find out who is the root of the i and j direction communications.
446 call find_comm_roots(global_rank,row_root,col_root)
447
448 if(row_comm_rank==row_root) then
449 ! Distribute in J direction.
450 ! print *,'global scatterv from rank ',grid_rank
451 call mpi_scatterv(global, col_comm_count, col_comm_displ, ifi_mpi_real_kind, &
452 rearrange_row2d ,im*(jend-jsta+1), ifi_mpi_real_kind, col_root, col_comm, ierr)
453 ! Rearrange from i-j dimensions to i-j-rank dimensions
454 idxstart=1
455 do r=1,row_comm_size
456
457 ni_rank=row_iend(r)-row_ista(r)+1
458 nj_rank=jend-jsta+1
459
460 !$OMP PARALLEL DO COLLAPSE(2) PRIVATE(outindex)
461 do j=0,nj_rank-1
462 do i=0,ni_rank-1
463 outindex = idxstart + i + ni_rank*j
464 rearrange_row1d(outindex) = rearrange_row2d(i+row_ista(r),j+jsta)
465 enddo
466 enddo
467
468 idxstart = idxstart + ni_rank*nj_rank
469 enddo
470 endif
471
472 ! Distribute across row
473 call mpi_scatterv(rearrange_row1d,row_comm_count,row_comm_displ,ifi_mpi_real_kind, &
474 rearrange,(iend-ista+1)*(jend-jsta+1),ifi_mpi_real_kind,row_root,row_comm,ierr)
475
476 ! Copy back to contiguous local array.
477 ! NOTE: Does NOT update the halo!
478 !$OMP PARALLEL DO COLLAPSE(2) PRIVATE(i,j)
479 do j=jsta,jend
480 do i=ista,iend
481 local(i,j) = rearrange(i,j)
482 enddo
483 enddo
484
485 ! print *,'global_scatter success'
486
487 end subroutine global_scatter
488
489 ! --------------------------------------------------------------------
490
491 subroutine ifi_check(status,error_message)
492 use mpi, only: mpi_abort, mpi_comm_world
493 implicit none
494 integer(c_int64_t), intent(in) :: status
495 character(*), intent(in) :: error_message
496 integer :: ierr
497
498 ! Exit program if status is non-zero
499 if(status/=0) then
500 write(0,'("IFI Failed: ",A)') trim(error_message)
501 call mpi_abort(mpi_comm_world,1,ierr)
502 endif
503 end subroutine ifi_check
504
505 ! --------------------------------------------------------------------
506
507 subroutine read_ifi_config()
508 implicit none
509 ! Only read the config if we haven't already:
510 if(.not.have_read_ifi_config) then
511 call ifi_check(ifi_config%init(),'read IFI config')
512 have_read_ifi_config = .true.
513 endif
514 end subroutine read_ifi_config
515
516 ! --------------------------------------------------------------------
517
518 subroutine set_ifi_dims()
519 use ctlblk_mod, only: ifi_nflight, ifi_flight_levels
520 implicit none
521
522 ! Warning: do not deallocate config_flight_levels_feet.
523 ! The IFI C++ library manages that array.
524 integer(kind=c_int64_t), pointer :: config_flight_levels_feet(:)
525 integer :: i
526
527 ! Make sure the config was read in
528 call read_ifi_config()
529
530 ! Get the flight levels
531 call ifi_check(ifi_config%get_flight_levels_feet(config_flight_levels_feet), &
532 'cannot get flight levels in feet from IFI config')
533
534 ! Convert from integer to real:
535 ifi_nflight = size(config_flight_levels_feet)
536 if(allocated(ifi_flight_levels)) then
537 deallocate(ifi_flight_levels)
538 endif
539 allocate(ifi_flight_levels(ifi_nflight))
540 ifi_flight_levels = config_flight_levels_feet
541
542! print '(A)','IFI flight levels:'
543! 38 format(' ifi_flight_level[',I0,'] = ',F15.5)
544! do i=1,ifi_nflight
545! print 38,i,ifi_flight_levels(i)
546! enddo
547
548 end subroutine set_ifi_dims
549
550 ! --------------------------------------------------------------------
551
552 subroutine send_data(vars,name,ient)
553 use ctlblk_mod, only: spval, jsta, jend, lm, cfld, datapd, fld_info, ifi_flight_levels, jsta_2l, jend_2u, ista, iend, ista_2l, iend_2u
554 use rqstfld_mod, only: iget, iavblfld, lvlsxml, lvls
555 implicit none
556
557 class(ifidata) :: vars
558 integer, intent(in) :: ient
559 character(*), intent(in) :: name
560
561 ! Locals
562 integer(c_int64_t) :: &
563 ims,ime,jms,jme,kms,kme, ids,ide,jds,jde,kds,kde, ips,ipe,jps,jpe,kps,kpe
564 real(kind=ifi_real_t), pointer :: data(:)
565 logical(c_bool) :: missing_value_is_set
566 real(kind=ifi_real_t) :: missing_value
567 integer :: i,j,k,nj_local,jpad,ilen,jlen,kstartm1,jstartm1,iloc,ndata,ji,ni_local,ipad
568
569 if(.not. iget(ient)>0) then
570 return
571 endif
572
573 ! WARNING: do not deallocate the data pointer. It is managed by the IFI C++ library.
574 data => vars%get_data(trim(name),missing_value_is_set,missing_value, &
575 ims,ime,jms,jme,kms,kme, ids,ide,jds,jde,kds,kde, ips,ipe,jps,jpe,kps,kpe)
576 if(.not.missing_value_is_set) then
577 missing_value = missing
578 endif
579
580 ! Get dimensions and do some sanity checks
581 nj_local = jpe-jps+1
582 ni_local = ipe-ips+1
583 jpad = jps-jms
584 ipad = ips-ims
585 if(nj_local/=jend-jsta+1 .or. jpad/=jsta-jsta_2l .or. ni_local/=iend-ista+1 .or. ipad/=ista-ista_2l) then
58694 format(' ',a,' = ',i0,', ',i0)
58783 format('Warning: ',a,': IFI output j bounds do not match UPP j bounds')
588 print 83,trim(name)
589 print 94,'jps,jsta',jps,jsta
590 print 94,'jpe,jend',jpe,jend
591 print 94,'jms,jsta_2l',jms,jsta_2l
592 print 94,'jme,jend_2u',jme,jend_2u
593 print 94,'jpad,jsta-jsta_2l',jpad,jsta-jsta_2l
594 print 94,'nj_local,jend-jsta+1',nj_local,jend-jsta+1
595 print 94,'ips,ista',ips,ista
596 print 94,'ipe,iend',ipe,iend
597 print 94,'ims,ista_2l',ims,ista_2l
598 print 94,'ime,iend_2u',ime,iend_2u
599 print 94,'ipad,ista-ista_2l',ipad,ista-ista_2l
600 print 94,'ni_local,iend-ista+1',ni_local,iend-ista+1
601 !call ifi_check(-1,'Internal error: IFI output j bounds do not match UPP j bounds')
602 end if
603 ilen=ime-ims+1
604 jlen=jme-jms+1
605 ndata=ilen*jlen*(kme-kms+1)
606
607 ! Go level-by-level writing grib2 if requested
608 do k=kps,kpe
609 kstartm1=ilen*jlen*(k-kms)
610 if(lvls(k,iget(ient))>0) then
611 cfld = cfld+1
612 fld_info(cfld)%ifld = iavblfld(iget(ient))
613 fld_info(cfld)%lvl = k ! ifi_flight_levels(k)*feet2meters
614 !$OMP PARALLEL DO PRIVATE(i,j,iloc,jstartm1) COLLAPSE(2)
615 do j=jps,jpe
616 do i=ips,ipe
617 jstartm1 = kstartm1 + (j-jms)*ilen
618 iloc = jstartm1+(i-ims)+1
619 if(iloc<1 .or. iloc>size(data)) then
620 call ifi_check(-1,'Internal error: out of bounds in send_data.')
621 endif
622 if(data(iloc)==missing_value) then
623 datapd(i-ips+1,j-jps+1,cfld) = spval
624 else
625 datapd(i-ips+1,j-jps+1,cfld) = data(iloc)
626 endif
627 enddo
628 enddo
629 endif
630 enddo
631 end subroutine send_data
632
633 ! --------------------------------------------------------------------
634
635 subroutine gather_ifi(local_buf_c,global_buf_c,receiving_rank) bind(C)
636 use iso_c_binding, only: c_float, c_double,c_int64_t,c_f_pointer
637 use mpi
638 use ctlblk_mod, only: me, mpi_comm_comp, jsta, jend, im, jm, icnt, idsp, ista, iend
639 implicit none
640
641 integer(kind=c_int64_t), value :: receiving_rank
642 type(c_ptr), value :: global_buf_c, local_buf_c
643 real(kind=ifi_real_t), pointer :: global_buf(:,:), local_buf(:,:)
644
645 integer :: type, iret
646
647 call c_f_pointer(local_buf_c,local_buf,(/ iend-ista+1, jend-jsta+1 /))
648 if(me==receiving_rank) then
649 call c_f_pointer(global_buf_c,global_buf,(/ im, jm /))
650 else
651 call c_f_pointer(global_buf_c,global_buf,(/ 1,1 /))
652 endif
653
654 call gather_for_write(local_buf,global_buf,receiving_rank)
655
656 end subroutine gather_ifi
657
658 ! --------------------------------------------------------------------
659
660 subroutine gather_for_write(local_buf,global_buf,receiving_rank_c) bind(C)
661 use iso_c_binding, only: c_float, c_double,c_int64_t,c_f_pointer
662 use mpi
663 use ctlblk_mod, only: me, mpi_comm_comp, jsta, jend, im, jm, ista, ista, iend
664 implicit none
665
666 integer(kind=c_int64_t), value :: receiving_rank_c
667 real(kind=ifi_real_t) :: global_buf(:,:), local_buf(:,:)
668 integer :: receiving_rank
669 receiving_rank = receiving_rank_c
670 call global_gather(local_buf,global_buf,receiving_rank,ista,iend,jsta,jend)
671 end subroutine gather_for_write
672
673 ! --------------------------------------------------------------------
674
675 subroutine scatter_ifi(local_buf_c,global_buf_c,sending_rank_c)bind(C)
676 use iso_c_binding, only: c_float, c_double, c_int64_t, c_f_pointer
677 use mpi
678 use ctlblk_mod, only: jsta, jend, im, jm, ista, iend, me
679 implicit none
680
681 integer(kind=c_int64_t), value :: sending_rank_c
682 type(c_ptr), value :: global_buf_c, local_buf_c
683 real(kind=ifi_real_t), pointer :: global_buf(:,:), local_buf(:,:)
684 integer :: sending_rank
685
686 sending_rank=sending_rank_c
687
688 call c_f_pointer(local_buf_c,local_buf,(/ iend-ista+1, jend-jsta+1 /))
689 if(me==sending_rank) then
690 call c_f_pointer(global_buf_c,global_buf,(/ im, jm /))
691 else
692 call c_f_pointer(global_buf_c,global_buf,(/ 1,1 /))
693 endif
694
695 call global_scatter(local_buf,global_buf,sending_rank,ista,iend,jsta,jend)
696 end subroutine scatter_ifi
697
698 ! --------------------------------------------------------------------
700 subroutine run_ifi()
701 use ctlblk_mod, only: spval, lm, lp1, im, jsta_2l,jend_2u, itprec, ifhr, ifmin, grib, &
702 jm, jsta,jend, me, num_procs, mpi_comm_comp, ista, iend, ista_2l, iend_2u, dtq2
703 use vrbls3d, only: &
704 zmid, & ! = IFI "HGT"
705 zint, & ! zint(:,:,LM+1) = IFI "HGT_surface"
706 qqi, & ! = IFI "CIMIXR"
707 qqg, & ! = IFI "GRLE"
708 qqw, & ! = IFI "CLMR"
709 pmid, & ! = IFI "PRES"
710 qqr, & ! "RWMR"
711 q, & ! "SPFH"
712 t, & ! "TMP"
713 qqs, & ! "SNMR"
714 omga ! "VVEL"
715 use vrbls2d, only : &
716 ifi_apcp, & ! IFI "APCP_surface" over ITPREC bucket time
717 cape, cin, &
718 avgprec_cont ! Backup plan for points where IFI_APCP is invalid or 0
719 use rqstfld_mod, only: iget
720 use iso_c_binding, only: c_bool, c_int64_t
721
722 implicit none
723
724 INTERFACE ! implemented in EXCH_c_float.f
725 SUBROUTINE exch_c_float(A) BIND(C)
726 use ifi_type_mod, only: ifi_real_t
727 implicit none
728 real(kind=ifi_real_t) :: a(*)
729 END SUBROUTINE exch_c_float
730 END INTERFACE
731
732 type(ifidata) :: hybr_vars, pres_vars, derived_vars, fip_algo_vars, flight_vars, cat_vars
733 type(ifialgo) :: algo
734 character(88) :: outfile
735 real :: to_hourly
736 real(c_double) :: fcst_lead_sec
737 integer :: i, j
738
739 if(ifhr < first_supported_ifi_hour .or. ifhr > last_supported_ifi_hour) then
7401838 format('You requested IFI fields for hour ',i0)
741 print 1838, ifhr
7421839 format('IFI fields are only available from hours ',i0,' through ',i0,'.')
743 print 1839, first_supported_ifi_hour, last_supported_ifi_hour
7441840 format(a)
745 print 1840, 'Will output missing data for IFI fields for this time.'
746 call send_missing_data(1007)
747 call send_missing_data(1008)
748 call send_missing_data(1009)
749 call send_missing_data(1010)
750 return
751 endif
752
753 ! FIXME: Once the post is i-decomposed, replace the appropriate 1..im in this routine.
754
755 if(grib/='grib2' .or. (iget(1007)<=0 .and. iget(1008)<=0 .and. iget(1009)<=0 .and. iget(1010)<=0)) then
756 if(me==0) then
757 print '(A)','IFI fields were not requested; skipping IFI'
758 endif
759 return ! nothing to do
760 endif
761
762 if(me==0) then
763 print '(A)','Running libIFI to get icing products.'
764 call ifi_print_copyright()
765
766 if(write_ifi_debug_files) then
767 print '(A)','Will output IFI debug files.'
768 else
769 print '(A)','Will NOT output IFI debug files.'
770 endif
771 endif
772
773 if(write_ifi_debug_files) then
774 call make_communicators()
775 endif
776
777 ! Read config and initialize input data structures:
778
779 call read_ifi_config()
780 call ifi_check(hybr_vars%init(int(1,c_int64_t),int(im,c_int64_t),int(1,c_int64_t),int(jm,c_int64_t),&
781 int(1,c_int64_t),int(lm,c_int64_t),int( ista,c_int64_t),int(iend,c_int64_t),int(jsta,c_int64_t),&
782 int(jend,c_int64_t),int(1,c_int64_t),int(lm,c_int64_t)),&
783 'could not initialize IFI input data structures')
784
785 ! Copy 2D vars to IFI internal structures:
786
787 call ifi_check(hybr_vars%add_ij_var('topography',int(ista_2l,c_int64_t),int(iend_2u,c_int64_t),&
788 int(jsta_2l,c_int64_t),int(jend_2u,c_int64_t), &
789 zint(:,:,lp1)),'could not send topography to IFI')
790
791 to_hourly=3600*1000.0/dtq2
792!$OMP PARALLEL DO COLLAPSE(2)
793 do j=jsta,jend
794 do i=ista,iend
795 if(ifi_apcp(i,j)>9e9 .or. ifi_apcp(i,j)<-9e9 .or. ifi_apcp(i,j)==0) then
796 ifi_apcp(i,j) = avgprec_cont(i,j)*to_hourly
797 else if(itprec>1e-5) then
798 ifi_apcp(i,j) = ifi_apcp(i,j)/itprec
799 endif
800 enddo
801 enddo
802
803 call ifi_check(hybr_vars%add_ij_var('APCP_surface',int(ista_2l,c_int64_t),int(iend_2u,c_int64_t),&
804 int(jsta_2l,c_int64_t),int(jend_2u,c_int64_t), &
805 ifi_apcp),'could not send IFI_APCP to IFI')
806 call ifi_check(hybr_vars%add_ij_var('CAPE_surface',int(ista_2l,c_int64_t),int(iend_2u,c_int64_t), &
807 int(jsta_2l,c_int64_t),int(jend_2u,c_int64_t),cape), &
808 'could not send CAPE to IFI')
809 call ifi_check(hybr_vars%add_ij_var('CIN_surface',int(ista_2l,c_int64_t),int(iend_2u,c_int64_t), &
810 int(jsta_2l,c_int64_t),int(jend_2u,c_int64_t),cin), &
811 'could not send CIN to IFI')
812
813 ! Copy 3D vars to IFI internal structures, inverting K dimension
814
815 call add_ijk_var('HGT','zmid',zmid)
816 call add_ijk_var('CIMIXR','QQI',qqi)
817 call add_ijk_var('CLMR','QQW',qqw)
818 call add_ijk_var('GRLE','QQG',qqg)
819 call add_ijk_var('PRES','pmid',pmid)
820 call add_ijk_var('RWMR','QQR',qqr)
821 call add_ijk_var('SPFH','Q',q)
822 call add_ijk_var('TMP','T',t)
823 call add_ijk_var('SNMR','QQS',qqs)
824 call add_ijk_var('VVEL','OMGA',omga)
825
826308 format(a,'_',i0,'.nc')
827
828 fcst_lead_sec=ifhr*3600.0+ifmin*60.0
829
830 if(write_ifi_debug_files) then
831 write(outfile,308) 'hybr_vars',me
832 call write_fip_output(hybr_vars,fcst_lead_sec,trim(outfile),.true.,'z0','z1')
833 endif
834
835 ! Initialize the IFI algorithm
836
837 call ifi_check(algo%init(ifi_config,fcst_lead_sec,hybr_vars,me,num_procs,mpi_comm_comp), &
838 'could not initialize IFI algorithm')
839
840 ! Run the IFI algorithm
841
842 call ifi_check(algo%calc_exner(),'calc_exner() failed')
843 call ifi_check(algo%hybrid_to_pressure(),'hybrid_to_pressure() failed')
844
845 call ifi_check(algo%get_pres_vars(pres_vars),'get_pres_vars()')
846 if(write_ifi_debug_files) then
847 write(outfile,308) 'pres_vars',me
848 call write_fip_output(pres_vars,fcst_lead_sec,trim(outfile),.true.,'z1','z0')
849 endif
850
851 call ifi_check(algo%discard_hybrid_level_vars(),'discard_hybrid_level_vars() failed')
852
853 call ifi_check(algo%derive_fields(),'derive_fields() failed')
854
855 call ifi_check(algo%get_derived_vars(derived_vars),'get_derived_vars()')
856 if(write_ifi_debug_files) then
857 write(outfile,308) 'derived_vars',me
858 call write_fip_output(derived_vars,fcst_lead_sec,trim(outfile),.false.,'z1','z0')
859 endif
860
861 call ifi_check(algo%run_fip_algo(),'run_fip_algo() failed')
862
863 call ifi_check(algo%get_fip_algo_vars(pres_vars),'get_fip_algo_vars()')
864 if(write_ifi_debug_files) then
865 write(outfile,308) 'fip_algo_vars',me
866 call write_fip_output(pres_vars,fcst_lead_sec,trim(outfile),.false.,'z1','z0')
867 endif
868
869 call ifi_check(algo%discard_pressure_level_vars(),'discard_pressure_level_vars() failed')
870 call ifi_check(algo%discard_derived_vars(),'discard_derived_vars() failed')
871 call ifi_check(algo%pressure_to_flight(),'pressure_to_flight() failed')
872
873 call ifi_check(algo%get_flight_vars(flight_vars),'get_flight_vars()')
874 if(write_ifi_debug_files) then
875 write(outfile,308) 'flight_vars',me
876 call write_fip_output(flight_vars,fcst_lead_sec,trim(outfile),.false.,'z1','z0')
877 endif
878
879 call ifi_check(algo%discard_fip_algo_vars(),'discard_fip_algo_vars() failed')
880 call ifi_check(algo%make_icing_category(),'make_icing_category() failed')
881
882 ! Get the final output fields:
883
884 call ifi_check(algo%get_cat_vars(cat_vars),'get_cat_vars()')
885
886 if(write_ifi_debug_files) then
887 write(outfile,308) 'cat_vars',me
888 call write_fip_output(cat_vars,fcst_lead_sec,trim(outfile),.false.,'z1','z0')
889 endif
890
891 call send_data(cat_vars,'ICE_PROB',1007)
892 call send_data(cat_vars,'SLD',1008)
893 call send_data(cat_vars,'ICE_SEV_CAT',1009)
894 call send_data(cat_vars,'WMO_ICE_SEV_CAT',1010)
895
896 ! When this subroutine ends, a Fortran-2003-compliant compiler
897 ! will free all memory IFI uses, by calling the destructors (final
898 ! routines) for cat_vars, hybr_vars, algo, and config.
899
900 contains
901
902 subroutine add_ijk_var(ifi_name,upp_name,upp_var)
903 implicit none
904 character(len=*), intent(in) :: ifi_name,upp_name
905 real, intent(in) :: upp_var(ista_2l:iend_2u,jsta_2l:jend_2u,lm)
906 real(kind=ifi_real_t) :: ifi_var(ista_2l:iend_2u,jsta_2l:jend_2u,lm)
907 integer i,j,k
908
909!$OMP PARALLEL DO COLLAPSE(2)
910 do k=1,lm
911 do j=jsta_2l,jend_2u
912 do i=ista_2l,iend_2u
913 ifi_var(i,j,k) = upp_var(i,j,lm-k+1)
914 enddo
915 enddo
916 enddo
917
918 call ifi_check(hybr_vars%add_ijk_var(ifi_name,int(ista_2l,c_int64_t),int(iend_2u,c_int64_t),&
919 int(jsta_2l,c_int64_t),int(jend_2u,c_int64_t),&
920 int(1,c_int64_t),int(lm,c_int64_t),ifi_var), &
921 'could not send '//ifi_name//' ('//upp_name//') to IFI')
922 end subroutine add_ijk_var
923
924 end subroutine run_ifi
925
926!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
927
928 subroutine find_range(var,count,missing_value_is_set,missing_value,min_not_miss,max_not_miss,all_missing)
929 USE ieee_arithmetic
930 use mpi
931 use ctlblk_mod, only: mpi_comm_comp
932 implicit none
933 real(kind=ifi_real_t), intent(in) :: var(*)
934 real(kind=ifi_real_t), intent(in) :: missing_value
935 real(kind=ifi_real_t), intent(out) :: min_not_miss,max_not_miss
936 integer, intent(in) :: count
937 logical, intent(out) :: all_missing
938 logical(c_bool), intent(in) :: missing_value_is_set
939
940 real(kind=ifi_real_t) :: minv,maxv,epsilon, global_minv,global_maxv
941 integer :: i,type,iret
942 real(kind=ifi_real_t), parameter :: zero = 0
943
944 print *,'find range begin'
945
946 if(missing_value_is_set) then
947 epsilon = abs(missing_value)*1e-4
948 if(.not. epsilon+1>epsilon) then
949 epsilon=1
950 endif
951 endif
952
953 ! Initialize to out-of-bounds values so they'll stay that way if no valid values are found:
954 minv = 1e30
955 maxv = -1e30
956
957 ! Find the min and max values in the array:
958 if(missing_value_is_set) then
959 !$OMP PARALLEL DO REDUCTION(min:minv) REDUCTION(max:maxv)
960 do i=1,count
961 if(.not. abs(var(i)-missing_value)>epsilon .and. var(i)<9e9 .and. var(i)>-9e9) then
962 minv=min(minv,var(i))
963 maxv=max(maxv,var(i))
964 endif
965 end do
966 else
967 !$OMP PARALLEL DO REDUCTION(min:minv) REDUCTION(max:maxv)
968 do i=1,count
969 if(var(i)+1>var(i)) then
970 if(var(i)<9e9 .and. var(i)>-9e9) then
971 minv=min(minv,var(i))
972 maxv=max(maxv,var(i))
973 endif
974 endif
975 end do
976 endif
977
978 if(ifi_real_t==c_double) then
979 type=mpi_real8
980 else
981 type=mpi_real4
982 endif
983
984 ! call MPI_Allreduce(minv,global_minv,1,type,MPI_MIN,mpi_comm_comp,iret)
985 ! call MPI_Allreduce(maxv,global_maxv,1,type,MPI_MAX,mpi_comm_comp,iret)
986 global_minv=minv
987 global_maxv=maxv
988
989 ! If min or max are inf, -inf, or NaN, assume all values are missing.
990 all_missing = global_minv<-9e9 .or. global_minv>9e9 .or. global_maxv<-9e9 .or. global_maxv>9e9
991
992 print *,'range min,max = ',global_minv,global_maxv
993
994 if(all_missing) then
995 min_not_miss=-1
996 max_not_miss=1
997 else
998 min_not_miss=global_minv
999 max_not_miss=global_maxv
1000 endif
1001
1002 print *,'find range end'
1003 end subroutine find_range
1004
1005!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1006
1007 subroutine nc_check(code,file,message)
1008 use netcdf
1009 use mpi
1010 implicit none
1011 integer, intent(in) :: code
1012 character(*), intent(in) :: file, message
1013 integer :: i, ierr
1014 if(code/=0) then
1015 write(0,93) trim(file),trim(message),trim(nf90_strerror(code)),code
1016 call mpi_abort(mpi_comm_world,1,ierr)
1017 endif
1018 93 format(a,': ',a,': ',a,'(',i0,')')
1019 end subroutine nc_check
1020
1021 subroutine write_fip_output(ifi_data,fcst_lead_sec,output_file,rename,z2dname,z3dname)
1022 use ifi_mod
1023 use iso_c_binding
1024 use netcdf
1025 use ctlblk_mod, only: me
1026 implicit none
1027 integer, parameter :: maxname=80
1028 integer, parameter :: maxvars=999
1029 logical, intent(in) :: rename
1030 double precision, intent(in) :: fcst_lead_sec
1031 character(len=*), intent(in) :: z2dname,z3dname
1032 type var_info
1033 character(len=maxname) :: varname, outname
1034 integer :: varid, ndims, dims(4), dimids(4), count
1035 real(ifi_real_t), pointer :: data(:)
1036 integer :: ims,ime,jms,jme,kms,kme
1037 integer :: ids,ide,jds,jde,kds,kde
1038 integer :: ips,ipe,jps,jpe,kps,kpe
1039 logical :: should_dealloc
1040 real(kind=ifi_real_t) :: missing_value
1041 real(kind=ifi_real_t) :: min_not_miss,max_not_miss
1042 logical(c_bool) :: missing_value_is_set
1043 logical :: all_missing
1044 end type var_info
1045 character(len=*), intent(in) :: output_file
1046
1047 type(ifidata) :: ifi_data
1048 integer(c_int64_t) :: ids,ide, jds,jde, kds,kde
1049 integer(c_int64_t) :: ips,ipe, jps,jpe, kps,kpe
1050 integer(c_int64_t) :: ims,ime, jms,jme, kms,kme
1051 integer :: count
1052
1053 character(len=200) :: varname,nextname
1054
1055 integer :: dimids(4), ncid, ivar, id, nvars, nx0, ny0, nz0, ntime, dims(4)
1056 type(var_info),target :: var_data(maxvars)
1057
1058 if(me==0) then
1059 write(*,'(A,A)') trim(output_file),': writing IFI debug data'
1060 endif
1061
1062 call ifi_check(ifi_data%get_dims(ids,ide, jds,jde, kds,kde, ips,ipe, jps,jpe, kps,kpe),"get_dims")
1063
1064 nx0=ide-ids+1
1065 ny0=jde-jds+1
1066 nz0=kde-kds+1
1067 ntime=1
1068
1069 dims = (/nx0,ny0,nz0,ntime/)
1070
1071 varname=' ' ! special value indicating "start at the beginning"
1072
1073 ivar=0
1074 do while(ivar<maxvars)
1075 nextname='internal_error_nextname_was_not_updated'
1076 call ifi_check(ifi_data%next_varname(varname,nextname),"next_varname")
1077 varname=nextname
1078
1079 if(len_trim(nextname)<1) then
1080 exit
1081 else
1082 ivar=ivar+1
1083 endif
1084
1085 if(ivar==1 .and. me==0) then
108659 format(a,': write to file')
1087 print 59,output_file
1088
1089 call nc_check(nf90_create(output_file,ior(nf90_noclobber,nf90_64bit_offset),ncid),&
1090 output_file,'nf90_create')
1091
1092 call nc_check(nf90_put_att(ncid,nf90_global,"fcst_lead_sec",fcst_lead_sec),&
1093 output_file,"nf90_put_att fcst_lead_sec")
1094
1095 call nc_check(nf90_def_dim(ncid,"x0",nx0,dimids(1)),output_file,"nf90_def_dim x0")
1096 call nc_check(nf90_def_dim(ncid,"y0",ny0,dimids(2)),output_file,"nf90_def_dim y0")
1097 call nc_check(nf90_def_dim(ncid,z3dname,nz0,dimids(3)),output_file,"nf90_def_dim "//z3dname)
1098 call nc_check(nf90_def_dim(ncid,"time",ntime,dimids(4)),output_file,"nf90_def_dim time")
1099 endif
1100
1101 var_data(ivar)%varname = varname
1102 var_data(ivar)%outname = varname
1103 var_data(ivar)%should_dealloc = .false.
1104 var_data(ivar)%data => read_var(var_data(ivar),trim(varname),var_data(ivar)%should_dealloc)
1105 if(me==0) then
1106 call set_dims(var_data(ivar),dims,dimids)
1107 var_data(ivar)%varid = def_var(var_data(ivar),rename)
1108 endif
1109 enddo
1110
1111 if(me/=0) then
1112 ! Ranks that are not writing are done now
1113 return
1114 endif
1115
1116 nvars=ivar
1117
1118 if(nvars<1) then
1119 write(*,"(A,A)") output_file,': no data to write!'
1120 return
1121 endif
1122
1123 call nc_check(nf90_enddef(ncid),output_file,'nf90_enddef')
1124
1125 write(*,*) 'before write loop ',me
1126 do ivar=1,nvars
112724 format(' put var ',a)
1128 print 24,trim(var_data(ivar)%varname)
1129
1130 call find_range(var_data(ivar)%data,var_data(ivar)%count, &
1131 var_data(ivar)%missing_value_is_set,var_data(ivar)%missing_value,&
1132 var_data(ivar)%min_not_miss,var_data(ivar)%max_not_miss, &
1133 var_data(ivar)%all_missing)
1134 call write_var(var_data(ivar))
1135 if(var_data(ivar)%should_dealloc) then
1136 deallocate(var_data(ivar)%data)
1137 endif
1138 nullify(var_data(ivar)%data)
1139 enddo
1140
1141 call nc_check(nf90_close(ncid),output_file,"nf90_close")
1142 contains
1143
1144 subroutine set_dims(var,dims,dimids)
1145 implicit none
1146 type(var_info), intent(inout), target :: var
1147 integer, intent(in) :: dims(:), dimids(:)
1148 character(len=:), pointer :: varname
1149
1150 varname=>var%varname(1:len_trim(var%varname))
1151
1152 if(me/=0) then
1153 ! Ranks that are not writing are done now
1154 return
1155 endif
1156
1157 ! These must also be in IFITest.cc IFITest::write_netcdf
1158
1159 if(varname=='x' .or. varname=='x0') then
1160 var%dims=(/nx0,1,1,1/)
1161 var%dimids=(/dimids(1),-1,-1,-1/)
1162 var%ndims=1
1163 else if(varname=='y' .or. varname=='y0') then
1164 var%dims=(/ny0,1,1,1/)
1165 var%dimids=(/dimids(2),-1,-1,-1/)
1166 var%ndims=1
1167 else if(varname=='z' .or. varname=='z0' .or. varname=='z1' .or. &
1168 varname=='pressure_levels' .or. varname=='exner_levels') then
1169 var%dims=(/nz0,1,1,1/)
1170 var%dimids=(/dimids(3),-1,-1,-1/)
1171 var%ndims=1
1172 else if(varname=='latitude' .or. varname=='longitude') then
1173 var%dims=(/nx0,ny0,1,1/)
1174 var%dimids=(/dimids(1),dimids(2),-1,-1/)
1175 var%ndims=2
1176 else if(varname=='time') then
1177 var%dims=(/ntime,1,1,1/)
1178 var%dimids=(/dimids(4),-1,-1,-1/)
1179 var%ndims=1
1180 else if(kme<=kms) then
1181 var%dims=(/nx0,ny0,ntime,1/)
1182 var%dimids=(/dimids(1),dimids(2),dimids(4),-1/)
1183 var%ndims=3
1184 else
1185 var%dims=(/nx0,ny0,nz0,ntime/)
1186 var%dimids=dimids
1187 var%ndims=4
1188 endif
1189
1190 var%count=var%dims(1)*var%dims(2)*var%dims(3)*var%dims(4)
1191 end subroutine set_dims
1192
1193 subroutine write_var(var)
1194 implicit none
1195 type(var_info), intent(inout) :: var
1196
1197 integer :: ones(var%ndims), dims(var%ndims)
1198 real(ifi_real_t) :: put(var%count)
1199 integer :: i,j,k,n,ilen,jlen,m,imlen,jmlen
1200
1201 ones = 1
1202 dims = var%dims(1:var%ndims)
1203
1204 if(me/=0) then
1205 ! Ranks that are not writing are done now
1206 return
1207 endif
1208
1209 call nc_check(nf90_put_var(ncid=ncid,varid=var%varid,values=var%data, &
1210 start=ones,count=dims), &
1211 output_file,"nf90_put_var "//trim(var%outname))
1212
1213 call nc_check(nf90_put_att(ncid,var%varid,"min_value",var%min_not_miss), &
1214 output_file,"nf90_put_att "//trim(var%outname)//" min_value")
1215
1216 call nc_check(nf90_put_att(ncid,var%varid,"max_value",var%max_not_miss), &
1217 output_file,"nf90_put_att "//trim(var%outname)//" max_value")
1218
1219 end subroutine write_var
1220
1221 integer function def_var(var,rename)
1222 use iso_c_binding, only: c_float
1223 implicit none
1224 logical :: rename
1225 type(var_info), intent(inout) :: var
1226
1227 integer :: varid, xtype, dimids(var%ndims)
1228 character(len=100) :: outname
1229
1230 if(rename) then
1231 var%outname=var%varname
1232 select case(trim(var%varname))
1233 case('CIMIXR')
1234 var%outname = 'ICMR'
1235 case('SPFH')
1236 var%outname = 'MIXR'
1237 case('GRLE')
1238 var%outname = 'GRMR'
1239 case('CAPE_surface')
1240 var%outname = 'CAPE'
1241 case('CIN_surface')
1242 var%outname = 'CIN'
1243 case('CLMR')
1244 var%outname = 'CLWMR'
1245 case('APCP_surface')
1246 var%outname = 'APCP1Hr'
1247 end select
1248 endif
1249
1250 var%ims=ims ; var%ime=ime ; var%jms=jms ; var%jme=jme ; var%kms=kms ; var%kme=kme
1251 var%ids=ids ; var%ide=ide ; var%jds=jds ; var%jde=jde ; var%kds=kds ; var%kde=kde
1252 var%ips=ips ; var%ipe=ipe ; var%jps=jps ; var%jpe=jpe ; var%kps=kps ; var%kpe=kpe
1253
1254 dimids = var%dimids(1:var%ndims)
1255
1256 if(ifi_real_t==c_float) then
1257 xtype = nf90_float
1258 else
1259 xtype = nf90_double
1260 endif
1261
1262 call nc_check(nf90_def_var(ncid,trim(var%outname),xtype,dimids,def_var), &
1263 output_file,"nf90_def_var "//trim(var%outname))
1264
1265 call nc_check(nf90_put_att(ncid,def_var,"min_value",var%min_not_miss), &
1266 output_file,"nf90_put_att "//trim(var%outname)//" min_value")
1267
1268 call nc_check(nf90_put_att(ncid,def_var,"max_value",var%max_not_miss), &
1269 output_file,"nf90_put_att "//trim(var%outname)//" max_value")
1270
1271 if(var%missing_value_is_set) then
1272 call nc_check(nf90_put_att(ncid,def_var,"_FillValue",-9999.0), &
1273 output_file,"nf90_put_att "//trim(var%outname)//" _FillValue")
1274 end if
1275
1276 end function def_var
1277
1278 function read_var(var,varname,should_dealloc)
1279 use mpi
1280 implicit none
1281 type(var_info), intent(inout) :: var
1282 real(kind=ifi_real_t), pointer :: read_var(:)
1283 real(kind=ifi_real_t), pointer :: local_data_1d(:)
1284 character(len=*), intent(in) :: varname
1285 logical, intent(out) :: should_dealloc
1286
1287 real(kind=ifi_real_t), allocatable :: local_data(:,:),global_data(:,:)
1288 real(kind=ifi_real_t), pointer :: global_data_1d(:)
1289
1290 integer :: nxny_local,nxny_global,nz,count,i,j,k
1291 integer :: local_ilen,local_jlen,local_klen,local_index
1292 integer :: global_ilen,global_jlen,global_klen,global_index, ierr
1293 type(c_ptr) :: global_cptr,local_cptr
1294
1295 local_data_1d => ifi_data%get_data(trim(varname),var%missing_value_is_set,var%missing_value, &
1296 ims,ime,jms,jme,kms,kme, ids,ide,jds,jde,kds,kde, ips,ipe,jps,jpe,kps,kpe)
1297
1298 global_ilen=ide-ids+1
1299 global_jlen=jde-jds+1
1300 global_klen=kde-kds+1
1301 local_ilen=ipe-ips+1
1302 local_jlen=jpe-jps+1
1303 local_klen=kpe-kps+1
1304
1305 if(.not.associated(local_data_1d)) then
1306 write(0,38) trim(varname)
130738 format("IFI did not produce ",a," variable.")
1308 call mpi_abort(mpi_comm_world,1,ierr)
1309 end if
1310
1311 count=global_ilen*global_jlen*global_klen
1312 if(count<=0) then
1313 write(0,39) trim(varname),count
131439 format("IFI variable ",a," had no data (size=",i0,")")
1315 call mpi_abort(mpi_comm_world,1,ierr)
1316 endif
1317
1318 if(ids==ide .or. jds==jde) then
1319 ! This is a 1D variable so we write it as is from rank 0
1320 read_var => local_data_1d
1321 should_dealloc = .false.
1322 return
1323 endif
1324
132540 format('var ',a,' has bad ',a,'. Got: ',i0,' but expected: ',i0)
1326
1327 if(global_ilen /= nx0) then
1328 write(0,40) trim(varname),'global_ilen',global_ilen,nx0
1329 call mpi_abort(mpi_comm_world,1,ierr)
1330 endif
1331
1332 if(global_jlen /= ny0) then
1333 write(0,40) trim(varname),'global_jlen',global_jlen,ny0
1334 call mpi_abort(mpi_comm_world,1,ierr)
1335 endif
1336
1337 if(global_klen/=1 .and. global_klen /= nz0) then
1338 write(0,40) trim(varname),'global_klen',global_klen,nz0
1339 call mpi_abort(mpi_comm_world,1,ierr)
1340 endif
1341
1342
1343 allocate(local_data(local_ilen,local_jlen))
1344 if(me==0) then
1345 allocate(global_data(global_ilen,global_jlen))
1346 allocate(global_data_1d(global_ilen*global_jlen*global_klen))
1347 else
1348 allocate(global_data(1,1))
1349 allocate(global_data_1d(1))
1350 endif
1351
1352 do k=kds,kde
1353 do j=jps,jpe
1354 do i=ips,ipe
1355 local_data(i-ips+1,j-jps+1) = local_data_1d( 1 + (i-ims) + ((j-jms) + (k-kms)*(jme-jms+1))*(ime-ims+1) )
1356 enddo
1357 enddo
1358 call gather_for_write(local_data,global_data,0)
1359 if(me==0) then
1360 do j=1,global_jlen
1361 do i=1,global_ilen
1362 global_data_1d(1 + (i-1) + ((j-1) + (k-kds)*global_jlen)*global_ilen) = &
1363 global_data(i,j)
1364 enddo
1365 enddo
1366 endif
1367 enddo
1368
1369 ! Ranks that are not writing data are done.
1370 if(me/=0) then
1371 ! This rank does not have the global data, and will not use the data anyway,
1372 ! but it wants an array, so we'll send it the local data that is managed
1373 ! internally by libIFI.
1374 read_var => local_data_1d
1375 should_dealloc = .false.
1376 deallocate(global_data_1d)
1377 else
1378 read_var => global_data_1d
1379 should_dealloc = .true.
1380 endif
1381
1382 deallocate(global_data)
1383 deallocate(local_data)
1384 end function read_var
1385 end subroutine write_fip_output
1386
1387#endif
1388
1389end module upp_ifi_mod
subroutine exch_c_float(a)
EXCH_c_float Subroutines that exchange one halo row.