sfc_climo_gen 1.14.0
Loading...
Searching...
No Matches
search_frac_cats.f90
Go to the documentation of this file.
1
4
24 subroutine search_frac_cats (field, mask, idim, jdim, num_categories, tile, field_name)
25
26 use mpi
27 use esmf
28
29 implicit none
30
31 integer, intent(in) :: idim, jdim, tile, num_categories
32 integer(esmf_kind_i4), intent(in) :: mask(idim,jdim)
33
34 real(esmf_kind_r4), intent(inout) :: field(idim,jdim,num_categories)
35
36 character(len=*) :: field_name
37
38 integer :: i, j, krad, ii, jj
39 integer :: istart, iend
40 integer :: jstart, jend
41 integer :: ierr
42 integer :: default_category
43
44 real(esmf_kind_r4), allocatable :: field_save(:,:,:)
45
46!-----------------------------------------------------------------------
47! Set default category.
48!-----------------------------------------------------------------------
49
50 select case (field_name)
51 case ('soil_type') ! soil type
52 default_category = 3
53 case ('vegetation_type') ! vegetation type
54 default_category = 3
55 case default
56 print*,'- FATAL ERROR IN ROUTINE SEARCH. UNIDENTIFIED FIELD : ', field
57 call mpi_abort(mpi_comm_world, 77, ierr)
58 end select
59
60!-----------------------------------------------------------------------
61! Perform search and replace.
62!-----------------------------------------------------------------------
63
64 allocate (field_save(idim,jdim,num_categories))
65 field_save = field
66
67 j_loop : do j = 1, jdim
68 i_loop : do i = 1, idim
69
70 if (mask(i,j) == 1 .and. field_save(i,j,1) < -9999.0) then
71
72 krad_loop : do krad = 1, 100
73
74 istart = i - krad
75 iend = i + krad
76 jstart = j - krad
77 jend = j + krad
78
79 jj_loop : do jj = jstart, jend
80 ii_loop : do ii = istart, iend
81
82!-----------------------------------------------------------------------
83! Search only along outer square.
84!-----------------------------------------------------------------------
85
86 if ((jj == jstart) .or. (jj == jend) .or. &
87 (ii == istart) .or. (ii == iend)) then
88
89 if (jj < 1 .or. jj > jdim) cycle jj_loop
90 if (ii < 1 .or. ii > idim) cycle ii_loop
91
92 if (mask(ii,jj) == 1 .and. maxval(field_save(ii,jj,:)) > 0.0) then
93 field(i,j,:) = field_save(ii,jj,:)
94 write(6,100) tile,i,j,ii,jj
95 cycle i_loop
96 endif
97
98 endif
99
100 enddo ii_loop
101 enddo jj_loop
102
103 enddo krad_loop
104
105 field(i,j,:) = 0.0
106 field(i,j,default_category) = 1.0 ! Search failed. Use 100% of default category.
107
108 write(6,101) tile,i,j,default_category
109
110 endif
111 enddo i_loop
112 enddo j_loop
113
114 deallocate(field_save)
115
116 100 format(1x,"- MISSING POINT TILE: ",i2," I/J: ",i5,i5," SET TO VALUE AT: ",i5,i5)
117 101 format(1x,"- MISSING POINT TILE: ",i2," I/J: ",i5,i5," SET TO DEFAULT VALUE OF: ",i3)
118
119 end subroutine search_frac_cats
subroutine search_frac_cats(field, mask, idim, jdim, num_categories, tile, field_name)
Replace undefined values on the model grid with valid values at a nearby neighbor.