FV3 Bundle
fv_surf_map_nlm.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU General Public License *
3 !* This file is a part of fvGFS. *
4 !* *
5 !* fvGFS is free software; you can redistribute it and/or modify it *
6 !* and are expected to follow the terms of the GNU General Public *
7 !* License as published by the Free Software Foundation; either *
8 !* version 2 of the License, or (at your option) any later version. *
9 !* *
10 !* fvGFS is distributed in the hope that it will be useful, but *
11 !* WITHOUT ANY WARRANTY; without even the implied warranty of *
12 !* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
13 !* General Public License for more details. *
14 !* *
15 !* For the full text of the GNU General Public License, *
16 !* write to: Free Software Foundation, Inc., *
17 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
18 !* or see: http://www.gnu.org/licenses/gpl.html *
19 !***********************************************************************
21 
22  use fms_mod, only: file_exist, check_nml_error, &
23  open_namelist_file, close_file, stdlog, &
24  mpp_pe, mpp_root_pe, fatal, error_mesg
25  use mpp_mod, only: get_unit, input_nml_file, mpp_error
27  use constants_mod, only: grav, radius, pi=>pi_8
28 
31  use fv_mp_nlm_mod, only: ng
32  use fv_mp_nlm_mod, only: mp_stop, mp_reduce_min, mp_reduce_max, is_master
35 
36  implicit none
37 
38  private
39 !-----------------------------------------------------------------------
40 ! NAMELIST
41 ! Name, resolution, and format of XXmin USGS datafile
42 ! 1min ---------> 1.85 km
43 ! nlon = 10800 * 2
44 ! nlat = 5400 * 2
45 ! 2min ---------> 3.7 km
46 ! nlon = 10800
47 ! nlat = 5400
48 ! 5min
49 ! nlon = 4320
50 ! nlat = 2160
51 ! surf_format: netcdf (default)
52 ! binary
53 ! New NASA SRTM30 data: SRTM30.nc
54 ! nlon = 43200
55 ! nlat = 21600
56  logical:: zs_filter = .true.
57  logical:: zero_ocean = .true. ! if true, no diffusive flux into water/ocean area
58  integer :: nlon = 21600
59  integer :: nlat = 10800
60  real:: cd4 = 0.15 ! Dimensionless coeff for del-4 diffusion (with FCT)
61  real:: cd2 = -1. ! Dimensionless coeff for del-2 diffusion (-1 gives resolution-determined value)
62  real:: peak_fac = 1.05 ! overshoot factor for the mountain peak
63  real:: max_slope = 0.15 ! max allowable terrain slope: 1 --> 45 deg
64  ! 0.15 for C768 or lower; 0.25 C1536; 0.3 for C3072
65  integer:: n_del2_weak = 12
66  integer:: n_del2_strong = -1
67  integer:: n_del4 = -1
68 
69 
70  character(len=128):: surf_file = "INPUT/topo1min.nc"
71  character(len=6) :: surf_format = 'netcdf'
72  logical :: namelist_read = .false.
73 
74  real(kind=R_GRID) da_min
75  real cos_grid
76  character(len=3) :: grid_string = ''
77 
78  namelist /surf_map_nml/ surf_file,surf_format,nlon,nlat, zero_ocean, zs_filter, &
80 !
81  real, allocatable:: zs_g(:,:), sgh_g(:,:), oro_g(:,:)
82 
83  public sgh_g, oro_g, zs_g
84  public surfdrv
86 
87  contains
88 
89  subroutine surfdrv(npx, npy, grid, agrid, area, dx, dy, dxa, dya, dxc, dyc, sin_sg, phis, &
90  stretch_fac, nested, npx_global, domain,grid_number, bd)
91 
92  implicit none
93 #include <netcdf.inc>
94  integer, intent(in):: npx, npy
95 
96  ! INPUT arrays
97  type(fv_grid_bounds_type), intent(IN) :: bd
98  real(kind=R_GRID), intent(in)::area(bd%is-ng:bd%ie+ng, bd%js-ng:bd%je+ng)
99  real, intent(in):: dx(bd%is-ng:bd%ie+ng, bd%js-ng:bd%je+ng+1)
100  real, intent(in):: dy(bd%is-ng:bd%ie+ng+1, bd%js-ng:bd%je+ng)
101  real, intent(in), dimension(bd%is-ng:bd%ie+ng, bd%js-ng:bd%je+ng)::dxa, dya
102  real, intent(in)::dxc(bd%is-ng:bd%ie+ng+1, bd%js-ng:bd%je+ng)
103  real, intent(in)::dyc(bd%is-ng:bd%ie+ng, bd%js-ng:bd%je+ng+1)
104 
105  real(kind=R_GRID), intent(in):: grid(bd%is-ng:bd%ie+ng+1, bd%js-ng:bd%je+ng+1,2)
106  real(kind=R_GRID), intent(in):: agrid(bd%is-ng:bd%ie+ng, bd%js-ng:bd%je+ng,2)
107  real, intent(IN):: sin_sg(bd%isd:bd%ied,bd%jsd:bd%jed,9)
108  real(kind=R_GRID), intent(IN):: stretch_fac
109  logical, intent(IN) :: nested
110  integer, intent(IN) :: npx_global
111  type(domain2d), intent(INOUT) :: domain
112  integer, intent(IN) :: grid_number
113 
114  ! OUTPUT arrays
115  real, intent(out):: phis(bd%is-ng:bd%ie+ng, bd%js-ng:bd%je+ng)
116 ! Local:
117  real, allocatable :: z2(:,:)
118 ! Position of edges of the box containing the original data point:
119  integer londim
120  integer latdim
121  character(len=80) :: topoflnm
122  real(kind=4), allocatable :: ft(:,:), zs(:,:)
123  real, allocatable :: lon1(:), lat1(:)
124  real dx1, dx2, dy1, dy2, lats, latn, r2d
125  real(kind=R_GRID) da_max
126  real zmean, z2mean, delg, rgrav
127 ! real z_sp, f_sp, z_np, f_np
128  integer i, j, n, mdim
129  integer igh, jt
130  integer ncid, lonid, latid, ftopoid, htopoid
131  integer jstart, jend, start(4), nread(4)
132  integer status
133 
134  integer :: is, ie, js, je
135  integer :: isd, ied, jsd, jed
136  real phis_coarse(bd%isd:bd%ied, bd%jsd:bd%jed)
137  real wt
138 
139  is = bd%is
140  ie = bd%ie
141  js = bd%js
142  je = bd%je
143  isd = bd%isd
144  ied = bd%ied
145  jsd = bd%jsd
146  jed = bd%jed
147  if (nested) then
148  !Divide all by grav
149  rgrav = 1./grav
150  do j=jsd,jed
151  do i=isd,ied
152  phis(i,j) = phis(i,j)*rgrav
153  enddo
154  enddo
155  !Save interpolated coarse-grid data for blending
156  do j=jsd,jed
157  do i=isd,ied
158  phis_coarse(i,j) = phis(i,j)
159  enddo
160  enddo
161  endif
162 
163  do j=js,je
164  do i=is,ie
165  phis(i,j) = 0.0
166  enddo
167  enddo
168 
169  call read_namelist
170 
171  if (grid_number > 1) write(grid_string, '(A, I1)') ' g', grid_number
172 
173 
174 !
175 ! surface file must be in NetCDF format
176 !
177  if ( file_exist(surf_file) ) then
178  if (surf_format == "netcdf") then
179 
180  status = nf_open(surf_file, nf_nowrite, ncid)
181  if (status .ne. nf_noerr) call handle_err(status)
182 
183  status = nf_inq_dimid(ncid, 'lon', lonid)
184  if (status .ne. nf_noerr) call handle_err(status)
185  status = nf_inq_dimlen(ncid, lonid, londim)
186  if (status .ne. nf_noerr) call handle_err(status)
187  nlon = londim
188 
189  status = nf_inq_dimid(ncid, 'lat', latid)
190  if (status .ne. nf_noerr) call handle_err(status)
191  status = nf_inq_dimlen(ncid, latid, latdim)
192  if (status .ne. nf_noerr) call handle_err(status)
193  nlat = latdim
194 
195  if ( is_master() ) then
196  if ( nlon==43200 ) then
197  write(*,*) 'Opening NASA datset file:', surf_file, surf_format, nlon, nlat
198  else
199  write(*,*) 'Opening USGS datset file:', surf_file, surf_format, nlon, nlat
200  endif
201  endif
202 
203  else
204  call error_mesg ( 'surfdrv','Raw IEEE data format no longer supported !!!', fatal )
205  endif
206  else
207  call error_mesg ( 'surfdrv','surface file '//trim(surf_file)//' not found !', fatal )
208  endif
209 
210  allocate ( lat1(nlat+1) )
211  allocate ( lon1(nlon+1) )
212 
213  r2d = 180./pi
214 
215  cos_grid = cos( 2.*pi/real(nlat) ) ! two-data_grid distance
216 
217  dx1 = 2.*pi/real(nlon)
218  dy1 = pi/real(nlat)
219 
220  do i=1,nlon+1
221  lon1(i) = dx1 * real(i-1) ! between 0 2pi
222  enddo
223 
224  lat1(1) = - 0.5*pi
225  lat1(nlat+1) = 0.5*pi
226  do j=2,nlat
227  lat1(j) = -0.5*pi + dy1*(j-1)
228  enddo
229 
230 !-------------------------------------
231 ! Compute raw phis and oro
232 !-------------------------------------
233  call timing_on('map_to_cubed')
234 
235  if (surf_format == "netcdf") then
236 
237 ! Find latitude strips reading data
238  lats = pi/2.
239  latn = -pi/2.
240  do j=js,je
241  do i=is,ie
242  lats = min( lats, grid(i,j,2), grid(i+1,j,2), grid(i,j+1,2), grid(i+1,j+1,2), agrid(i,j,2) )
243  latn = max( latn, grid(i,j,2), grid(i+1,j,2), grid(i,j+1,2), grid(i+1,j+1,2), agrid(i,j,2) )
244  enddo
245  enddo
246 
247 ! Enlarge the search zone:
248 ! To account for the curvature of the coordinates:
249 
250  !I have had trouble running c90 with 600 pes unless the search region is expanded
251  ! due to failures in finding latlon points in the source data.
252  !This sets a larger search region if the number of cells on a PE is too small.
253  !(Alternately you can just cold start the topography using a smaller number of PEs)
254  if (min(je-js+1,ie-is+1) < 15) then
255  delg = max( 0.4*(latn-lats), pi/real(npx_global-1), 2.*pi/real(nlat) )
256  else
257  delg = max( 0.2*(latn-lats), pi/real(npx_global-1), 2.*pi/real(nlat) )
258  endif
259  lats = max( -0.5*pi, lats - delg )
260  latn = min( 0.5*pi, latn + delg )
261 
262  jstart = 1
263  do j=2,nlat
264  if ( lats < lat1(j) ) then
265  jstart = j-1
266  exit
267  endif
268  enddo
269  jstart = max(jstart-1, 1)
270 
271  jend = nlat
272  do j=2,nlat
273  if ( latn < lat1(j+1) ) then
274  jend = j+1
275  exit
276  endif
277  enddo
278  jend = min(jend+1, nlat)
279 
280  jt = jend - jstart + 1
281  igh = nlon/8 + nlon/(2*(npx_global-1))
282 
283  if (is_master()) write(*,*) 'Terrain dataset =', nlon, 'jt=', jt
284  if (is_master()) write(*,*) 'igh (terrain ghosting)=', igh
285 
286  status = nf_inq_varid(ncid, 'ftopo', ftopoid)
287  if (status .ne. nf_noerr) call handle_err(status)
288  nread = 1; start = 1
289  nread(1) = nlon
290  start(2) = jstart; nread(2) = jend - jstart + 1
291 
292  allocate ( ft(-igh:nlon+igh,jt) )
293  status = nf_get_vara_real(ncid, ftopoid, start, nread, ft(1:nlon,1:jt))
294  if (status .ne. nf_noerr) call handle_err(status)
295 
296  do j=1,jt
297  do i=-igh,0
298  ft(i,j) = ft(i+nlon,j)
299  enddo
300  do i=nlon+1,nlon+igh
301  ft(i,j) = ft(i-nlon,j)
302  enddo
303  enddo
304 
305  status = nf_inq_varid(ncid, 'htopo', htopoid)
306  if (status .ne. nf_noerr) call handle_err(status)
307  allocate ( zs(-igh:nlon+igh,jt) )
308  status = nf_get_vara_real(ncid, htopoid, start, nread, zs(1:nlon,1:jt))
309  if (status .ne. nf_noerr) call handle_err(status)
310  status = nf_close(ncid)
311  if (status .ne. nf_noerr) call handle_err(status)
312 ! Ghost Data
313  do j=1,jt
314  do i=-igh,0
315  zs(i,j) = zs(i+nlon,j)
316  enddo
317  do i=nlon+1,nlon+igh
318  zs(i,j) = zs(i-nlon,j)
319  enddo
320  enddo
321 
322  endif
323 
324 ! special SP treatment:
325 ! if ( jstart == 1 ) then
326 ! call zonal_mean(nlon, zs(1,1), z_sp)
327 ! call zonal_mean(nlon, ft(1,1), f_sp)
328 ! endif
329 
330  allocate ( oro_g(isd:ied, jsd:jed) )
331  allocate ( sgh_g(isd:ied, jsd:jed) )
332  call timing_on('map_to_cubed')
333  call map_to_cubed_raw(igh, nlon, jt, lat1(jstart:jend+1), lon1, zs, ft, grid, agrid, &
334  phis, oro_g, sgh_g, npx, npy, jstart, jend, stretch_fac, nested, npx_global, bd)
335  if (is_master()) write(*,*) 'map_to_cubed_raw: master PE done'
336  call timing_off('map_to_cubed')
337 
338  deallocate ( zs )
339  deallocate ( ft )
340  deallocate ( lon1 )
341  deallocate ( lat1 )
342 
343  allocate ( zs_g(is:ie, js:je) )
344  allocate ( z2(is:ie,js:je) )
345  do j=js,je
346  do i=is,ie
347  zs_g(i,j) = phis(i,j)
348  z2(i,j) = phis(i,j)**2
349  enddo
350  enddo
351 !--------
352 ! Filter:
353 !--------
354  call global_mx(real(phis,kind=R_GRID), ng, da_min, da_max, bd)
355  zmean = g_sum(domain, zs_g(is:ie,js:je), is, ie, js, je, ng, area, 1)
356  z2mean = g_sum(domain, z2(is:ie,js:je) , is, ie, js, je, ng, area, 1)
357  if ( is_master() ) then
358  write(*,*) 'Before filter ZS', trim(grid_string), ' min=', da_min, ' Max=', da_max,' Mean=',zmean
359  write(*,*) '*** Mean variance', trim(grid_string), ' *** =', z2mean
360  endif
361 
362  !On a nested grid blend coarse-grid and nested-grid
363  ! orography near boundary.
364  ! This works only on the height of topography; assume
365  ! land fraction and sub-grid variance unchanged
366 
367  ! Here, we blend in the four cells nearest to the boundary.
368  ! In the halo we set the value to that interpolated from the coarse
369  ! grid. (Previously this was erroneously not being done, which was causing
370  ! the halo to be filled with unfiltered nested-grid terrain, creating
371  ! ugly edge artifacts.)
372  if (nested) then
373  if (is_master()) write(*,*) 'Blending nested and coarse grid topography'
374  do j=jsd,jed
375  do i=isd,ied
376  wt = max(0.,min(1.,real(5 - min(i,j,npx-i,npy-j,5))/5. ))
377  phis(i,j) = (1.-wt)*phis(i,j) + wt*phis_coarse(i,j)
378 
379  enddo
380  enddo
381  endif
382 
383  call global_mx(real(oro_g,kind=R_GRID), ng, da_min, da_max, bd)
384  if ( is_master() ) write(*,*) 'ORO', trim(grid_string), ' min=', da_min, ' Max=', da_max
385 
386  call global_mx(area, ng, da_min, da_max, bd)
387  call timing_on('Terrain_filter')
388 ! Del-2: high resolution only
389  if ( zs_filter ) then
390  if(is_master()) then
391  write(*,*) 'Applying terrain filters. zero_ocean is', zero_ocean
392  endif
393  call fv3_zs_filter (bd, isd, ied, jsd, jed, npx, npy, npx_global, &
394  stretch_fac, nested, domain, area, dxa, dya, dx, dy, dxc, dyc, grid, &
395  agrid, sin_sg, phis, oro_g)
396  call mpp_update_domains(phis, domain)
397  endif ! end terrain filter
398  call timing_off('Terrain_filter')
399 
400  do j=js,je
401  do i=is,ie
402  z2(i,j) = phis(i,j)**2
403  end do
404  end do
405 
406  call global_mx(real(phis,kind=R_GRID), ng, da_min, da_max, bd)
407  zmean = g_sum(domain, phis(is:ie,js:je), is, ie, js, je, ng, area, 1)
408  z2mean = g_sum(domain, z2, is, ie, js, je, ng, area, 1)
409  deallocate ( z2 )
410 
411  if ( is_master() ) then
412  write(*,*) 'After filter Phis', trim(grid_string), ' min=', da_min, ' Max=', da_max, 'Mean=', zmean
413  write(*,*) '*** Mean variance', trim(grid_string), ' *** =', z2mean
414  endif
415 
416  !FOR NESTING: Unless we fill the outermost halo with topography
417  ! interpolated from the coarse grid, the values of phi there
418  ! will be wrong because they are just z instead of g*z; they
419  ! have not had gravity properly multiplied in so we can get phi
420 
421  ! For now we compute phis and sgh_g on the full data domain; for
422  ! nested grids this allows us to do the smoothing near the boundary
423  ! without having to fill the boundary halo from the coarse grid
424 
425  !ALSO for nesting: note that we are smoothing the terrain using
426  ! the nested-grid's outer halo filled with the terrain computed
427  ! directly from the input file computed here, and then
428  ! replacing it with interpolated topography in fv_restart, so
429  ! as to be consistent when doing the boundary condition
430  ! interpolation. We would ideally replace the nested-grid
431  ! halo topography BEFORE smoothing, which could be more
432  ! consistent, but this would require moving calls to
433  ! nested_grid_BC in this routine.
434 
435  do j=jsd,jed
436  do i=isd,ied
437  phis(i,j) = grav * phis(i,j)
438  if ( sgh_g(i,j) <= 0. ) then
439  sgh_g(i,j) = 0.
440  else
441  sgh_g(i,j) = sqrt(sgh_g(i,j))
442  endif
443  end do
444  end do
445 
446  call global_mx(real(sgh_g,kind=R_GRID), ng, da_min, da_max, bd)
447  if ( is_master() ) write(*,*) 'Before filter SGH', trim(grid_string), ' min=', da_min, ' Max=', da_max
448 
449 
450 !-----------------------------------------------
451 ! Filter the standard deviation of mean terrain:
452 !-----------------------------------------------
453  call global_mx(area, ng, da_min, da_max, bd)
454 
455  if(zs_filter) call del4_cubed_sphere(npx, npy, sgh_g, area, dx, dy, dxc, dyc, sin_sg, 1, zero_ocean, oro_g, nested, domain, bd)
456 
457  call global_mx(real(sgh_g,kind=R_GRID), ng, da_min, da_max, bd)
458  if ( is_master() ) write(*,*) 'After filter SGH', trim(grid_string), ' min=', da_min, ' Max=', da_max
459  do j=js,je
460  do i=is,ie
461  sgh_g(i,j) = max(0., sgh_g(i,j))
462  enddo
463  enddo
464 
465  end subroutine surfdrv
466 
467  subroutine fv3_zs_filter (bd, isd, ied, jsd, jed, npx, npy, npx_global, &
468  stretch_fac, nested, domain, area, dxa, dya, dx, dy, dxc, dyc, grid, &
469  agrid, sin_sg, phis, oro )
470  integer, intent(in):: isd, ied, jsd, jed, npx, npy, npx_global
471  type(fv_grid_bounds_type), intent(IN) :: bd
472  real(kind=R_GRID), intent(in), dimension(isd:ied,jsd:jed)::area
473  real, intent(in), dimension(isd:ied,jsd:jed)::dxa, dya
474  real, intent(in), dimension(isd:ied, jsd:jed+1):: dx, dyc
475  real, intent(in), dimension(isd:ied+1,jsd:jed):: dy, dxc
476 
477  real(kind=R_GRID), intent(in):: grid(isd:ied+1, jsd:jed+1,2)
478  real(kind=R_GRID), intent(in):: agrid(isd:ied, jsd:jed, 2)
479  real, intent(IN):: sin_sg(9,isd:ied,jsd:jed)
480  real(kind=R_GRID), intent(IN):: stretch_fac
481  logical, intent(IN) :: nested
482  real, intent(inout):: phis(isd:ied,jsd,jed)
483  real, intent(inout):: oro(isd:ied,jsd,jed)
484  type(domain2d), intent(INOUT) :: domain
485  integer mdim
486  real(kind=R_GRID) da_max
487 
488  if (is_master()) print*, ' Calling FV3_zs_filter...'
489 
490  if (.not. namelist_read) call read_namelist !when calling from external_ic
491  call global_mx(area, ng, da_min, da_max, bd)
492 
493  mdim = nint( real(npx_global) * min(10., stretch_fac) )
494 
495 ! Del-2: high resolution only
496 ! call del2_cubed_sphere(npx, npy, phis, area, dx, dy, dxc, dyc, sin_sg, n_del2, cd2, zero_ocean, oro, nested, domain, bd)
497  if (n_del2_strong < 0) then
498  if ( npx_global<=97) then
499  n_del2_strong = 0
500  elseif ( npx_global<=193 ) then
501  n_del2_strong = 1
502  else
503  n_del2_strong = 2
504  endif
505  endif
506  if (cd2 < 0.) cd2 = 0.16*da_min
507 ! Applying strong 2-delta-filter:
508  if ( n_del2_strong > 0 ) &
509  call two_delta_filter(npx, npy, phis, area, dx, dy, dxa, dya, dxc, dyc, sin_sg, cd2, zero_ocean, &
510  .true., 0, oro, nested, domain, bd, n_del2_strong)
511 
512 ! MFCT Del-4:
513  if (n_del4 < 0) then
514  if ( mdim<=193 ) then
515  n_del4 = 1
516  elseif ( mdim<=1537 ) then
517  n_del4 = 2
518  else
519  n_del4 = 3
520  endif
521  endif
522  call del4_cubed_sphere(npx, npy, phis, area, dx, dy, dxc, dyc, sin_sg, n_del4, zero_ocean, oro, nested, domain, bd)
523 ! Applying weak 2-delta-filter:
524  cd2 = 0.12*da_min
525  call two_delta_filter(npx, npy, phis, area, dx, dy, dxa, dya, dxc, dyc, sin_sg, cd2, zero_ocean, &
526  .true., 1, oro, nested, domain, bd, n_del2_weak)
527 
528 
529  end subroutine fv3_zs_filter
530 
531 
532  subroutine two_delta_filter(npx, npy, q, area, dx, dy, dxa, dya, dxc, dyc, sin_sg, cd, zero_ocean, &
533  check_slope, filter_type, oro, nested, domain, bd, ntmax)
534  type(fv_grid_bounds_type), intent(IN) :: bd
535  integer, intent(in):: npx, npy
536  integer, intent(in):: ntmax
537  integer, intent(in):: filter_type ! 0: strong, 1: weak
538  real, intent(in):: cd
539 ! INPUT arrays
540  real(kind=R_GRID), intent(in)::area(bd%isd:bd%ied, bd%jsd:bd%jed)
541  real, intent(in):: dx(bd%isd:bd%ied, bd%jsd:bd%jed+1)
542  real, intent(in):: dy(bd%isd:bd%ied+1,bd%jsd:bd%jed)
543  real, intent(in):: dxa(bd%isd:bd%ied, bd%jsd:bd%jed)
544  real, intent(in):: dya(bd%isd:bd%ied, bd%jsd:bd%jed)
545  real, intent(in):: dxc(bd%isd:bd%ied+1,bd%jsd:bd%jed)
546  real, intent(in):: dyc(bd%isd:bd%ied, bd%jsd:bd%jed+1)
547  real, intent(in):: sin_sg(9,bd%isd:bd%ied,bd%jsd:bd%jed)
548  real, intent(in):: oro(bd%isd:bd%ied, bd%jsd:bd%jed) ! 0==water, 1==land
549  logical, intent(in):: zero_ocean, check_slope
550  logical, intent(in):: nested
551  type(domain2d), intent(inout) :: domain
552 ! OUTPUT arrays
553  real, intent(inout):: q(bd%isd:bd%ied, bd%jsd:bd%jed)
554 ! Local:
555  real, parameter:: p1 = 7./12.
556  real, parameter:: p2 = -1./12.
557  real, parameter:: c1 = -2./14.
558  real, parameter:: c2 = 11./14.
559  real, parameter:: c3 = 5./14.
560  real:: ddx(bd%is:bd%ie+1,bd%js:bd%je), ddy(bd%is:bd%ie,bd%js:bd%je+1)
561  logical:: extm(bd%is-1:bd%ie+1)
562  logical:: ext2(bd%is:bd%ie,bd%js-1:bd%je+1)
563  real:: a1(bd%is-1:bd%ie+2)
564  real:: a2(bd%is:bd%ie,bd%js-1:bd%je+2)
565  real(kind=R_GRID):: a3(bd%is:bd%ie,bd%js:bd%je)
566  real(kind=R_GRID):: smax, smin
567  real:: m_slope, fac
568  integer:: i,j, nt
569  integer:: is, ie, js, je
570  integer:: isd, ied, jsd, jed
571  integer:: is1, ie2, js1, je2
572 
573  is = bd%is
574  ie = bd%ie
575  js = bd%js
576  je = bd%je
577  isd = bd%isd
578  ied = bd%ied
579  jsd = bd%jsd
580  jed = bd%jed
581 
582  if ( nested ) then
583  is1 = is-1; ie2 = ie+2
584  js1 = js-1; je2 = je+2
585  else
586  is1 = max(3,is-1); ie2 = min(npx-2,ie+2)
587  js1 = max(3,js-1); je2 = min(npy-2,je+2)
588  end if
589 
590  if ( check_slope ) then
591  m_slope = max_slope
592  else
593  m_slope = 10.
594  endif
595 
596 
597  do 777 nt=1, ntmax
598  call mpp_update_domains(q, domain)
599 
600 ! Check slope
601  if ( nt==1 .and. check_slope ) then
602  do j=js,je
603  do i=is,ie+1
604  ddx(i,j) = (q(i,j) - q(i-1,j))/dxc(i,j)
605  ddx(i,j) = abs(ddx(i,j))
606  enddo
607  enddo
608  do j=js,je+1
609  do i=is,ie
610  ddy(i,j) = (q(i,j) - q(i,j-1))/dyc(i,j)
611  ddy(i,j) = abs(ddy(i,j))
612  enddo
613  enddo
614  do j=js,je
615  do i=is,ie
616  a3(i,j) = max( ddx(i,j), ddx(i+1,j), ddy(i,j), ddy(i,j+1) )
617  enddo
618  enddo
619  call global_mx(a3, 0, smin, smax, bd)
620  if ( is_master() ) write(*,*) 'Before filter: Max_slope=', smax
621  endif
622 
623 ! First step: average the corners:
624  if ( .not. nested .and. nt==1 ) then
625  if ( is==1 .and. js==1 ) then
626  q(1,1) = (q(1,1)*area(1,1)+q(0,1)*area(0,1)+q(1,0)*area(1,0)) &
627  / ( area(1,1)+ area(0,1)+ area(1,0) )
628  q(0,1) = q(1,1)
629  q(1,0) = q(1,1)
630  endif
631  if ( (ie+1)==npx .and. js==1 ) then
632  q(ie, 1) = (q(ie,1)*area(ie,1)+q(npx,1)*area(npx,1)+q(ie,0)*area(ie,0)) &
633  / ( area(ie,1)+ area(npx,1)+ area(ie,0))
634  q(npx,1) = q(ie,1)
635  q(ie, 0) = q(ie,1)
636  endif
637  if ( is==1 .and. (je+1)==npy ) then
638  q(1, je) = (q(1,je)*area(1,je)+q(0,je)*area(0,je)+q(1,npy)*area(1,npy)) &
639  / ( area(1,je)+ area(0,je)+ area(1,npy))
640  q(0, je) = q(1,je)
641  q(1,npy) = q(1,je)
642  endif
643  if ( (ie+1)==npx .and. (je+1)==npy ) then
644  q(ie, je) = (q(ie,je)*area(ie,je)+q(npx,je)*area(npx,je)+q(ie,npy)*area(ie,npy)) &
645  / ( area(ie,je)+ area(npx,je)+ area(ie,npy))
646  q(npx,je) = q(ie,je)
647  q(ie,npy) = q(ie,je)
648  endif
649  call mpp_update_domains(q, domain)
650  endif
651 
652 ! x-diffusive flux:
653  do 333 j=js,je
654 
655  do i=is1, ie2
656  a1(i) = p1*(q(i-1,j)+q(i,j)) + p2*(q(i-2,j)+q(i+1,j))
657  enddo
658 
659  if ( .not. nested ) then
660  if ( is==1 ) then
661  a1(0) = c1*q(-2,j) + c2*q(-1,j) + c3*q(0,j)
662  a1(1) = 0.5*(((2.*dxa(0,j)+dxa(-1,j))*q(0,j)-dxa(0,j)*q(-1,j))/(dxa(-1,j)+dxa(0,j)) &
663  + ((2.*dxa(1,j)+dxa( 2,j))*q(1,j)-dxa(1,j)*q( 2,j))/(dxa(1, j)+dxa(2,j)))
664  a1(2) = c3*q(1,j) + c2*q(2,j) +c1*q(3,j)
665  endif
666  if ( (ie+1)==npx ) then
667  a1(npx-1) = c1*q(npx-3,j) + c2*q(npx-2,j) + c3*q(npx-1,j)
668  a1(npx) = 0.5*(((2.*dxa(npx-1,j)+dxa(npx-2,j))*q(npx-1,j)-dxa(npx-1,j)*q(npx-2,j))/(dxa(npx-2,j)+dxa(npx-1,j)) &
669  + ((2.*dxa(npx, j)+dxa(npx+1,j))*q(npx, j)-dxa(npx, j)*q(npx+1,j))/(dxa(npx, j)+dxa(npx+1,j)))
670  a1(npx+1) = c3*q(npx,j) + c2*q(npx+1,j) + c1*q(npx+2,j)
671  endif
672  endif
673 
674  if ( filter_type == 0 ) then
675  do i=is-1, ie+1
676  if( abs(3.*(a1(i)+a1(i+1)-2.*q(i,j))) > abs(a1(i)-a1(i+1)) ) then
677  extm(i) = .true.
678  else
679  extm(i) = .false.
680  endif
681  enddo
682  else
683  do i=is-1, ie+1
684  if ( (a1(i)-q(i,j))*(a1(i+1)-q(i,j)) > 0. ) then
685  extm(i) = .true.
686  else
687  extm(i) = .false.
688  endif
689  enddo
690  endif
691 
692  do i=is,ie+1
693  ddx(i,j) = (q(i-1,j)-q(i,j))/dxc(i,j)
694  if ( extm(i-1).and.extm(i) ) then
695  ddx(i,j) = 0.5*(sin_sg(3,i-1,j)+sin_sg(1,i,j))*dy(i,j)*ddx(i,j)
696  elseif ( abs(ddx(i,j)) > m_slope ) then
697  fac = min(1., max(0.1,(abs(ddx(i,j))-m_slope)/m_slope ) )
698  ddx(i,j) = fac*0.5*(sin_sg(3,i-1,j)+sin_sg(1,i,j))*dy(i,j)*ddx(i,j)
699  else
700  ddx(i,j) = 0.
701  endif
702  enddo
703 333 continue
704 
705 ! y-diffusive flux:
706  do j=js1,je2
707  do i=is,ie
708  a2(i,j) = p1*(q(i,j-1)+q(i,j)) + p2*(q(i,j-2)+q(i,j+1))
709  enddo
710  enddo
711  if ( .not. nested ) then
712  if( js==1 ) then
713  do i=is,ie
714  a2(i,0) = c1*q(i,-2) + c2*q(i,-1) + c3*q(i,0)
715  a2(i,1) = 0.5*(((2.*dya(i,0)+dya(i,-1))*q(i,0)-dya(i,0)*q(i,-1))/(dya(i,-1)+dya(i,0)) &
716  + ((2.*dya(i,1)+dya(i, 2))*q(i,1)-dya(i,1)*q(i, 2))/(dya(i, 1)+dya(i,2)))
717  a2(i,2) = c3*q(i,1) + c2*q(i,2) + c1*q(i,3)
718  enddo
719  endif
720  if( (je+1)==npy ) then
721  do i=is,ie
722  a2(i,npy-1) = c1*q(i,npy-3) + c2*q(i,npy-2) + c3*q(i,npy-1)
723  a2(i,npy) = 0.5*(((2.*dya(i,npy-1)+dya(i,npy-2))*q(i,npy-1)-dya(i,npy-1)*q(i,npy-2))/(dya(i,npy-2)+dya(i,npy-1)) &
724  + ((2.*dya(i,npy)+dya(i,npy+1))*q(i,npy)-dya(i,npy)*q(i,npy+1))/(dya(i,npy)+dya(i,npy+1)))
725  a2(i,npy+1) = c3*q(i,npy) + c2*q(i,npy+1) + c1*q(i,npy+2)
726  enddo
727  endif
728  endif
729 
730  if ( filter_type == 0 ) then
731  do j=js-1,je+1
732  do i=is,ie
733  if( abs(3.*(a2(i,j)+a2(i,j+1)-2.*q(i,j))) > abs(a2(i,j)-a2(i,j+1)) ) then
734  ext2(i,j) = .true.
735  else
736  ext2(i,j) = .false.
737  endif
738  enddo
739  enddo
740  else
741  do j=js-1,je+1
742  do i=is,ie
743  if ( (a2(i,j)-q(i,j))*(a2(i,j+1)-q(i,j)) > 0. ) then
744  ext2(i,j) = .true.
745  else
746  ext2(i,j) = .false.
747  endif
748  enddo
749  enddo
750  endif
751 
752  do j=js,je+1
753  do i=is,ie
754  ddy(i,j) = (q(i,j-1)-q(i,j))/dyc(i,j)
755  if ( ext2(i,j-1) .and. ext2(i,j) ) then
756  ddy(i,j) = 0.5*(sin_sg(4,i,j-1)+sin_sg(2,i,j))*dx(i,j)*ddy(i,j)
757  elseif ( abs(ddy(i,j))>m_slope ) then
758  fac = min(1., max(0.1,(abs(ddy(i,j))-m_slope)/m_slope))
759  ddy(i,j) = fac*0.5*(sin_sg(4,i,j-1)+sin_sg(2,i,j))*dx(i,j)*ddy(i,j)
760  else
761  ddy(i,j) = 0.
762  endif
763  enddo
764  enddo
765 
766  if ( zero_ocean ) then
767 ! Limit diffusive flux over water cells:
768  do j=js,je
769  do i=is,ie+1
770  ddx(i,j) = max(0., min(oro(i-1,j), oro(i,j))) * ddx(i,j)
771  enddo
772  enddo
773  do j=js,je+1
774  do i=is,ie
775  ddy(i,j) = max(0., min(oro(i,j-1), oro(i,j))) * ddy(i,j)
776  enddo
777  enddo
778  endif
779 
780  do j=js,je
781  do i=is,ie
782  q(i,j) = q(i,j) + cd/area(i,j)*(ddx(i,j)-ddx(i+1,j)+ddy(i,j)-ddy(i,j+1))
783  enddo
784  enddo
785 777 continue
786 
787 ! Check slope
788  if ( check_slope ) then
789  call mpp_update_domains(q, domain)
790  do j=js,je
791  do i=is,ie+1
792  ddx(i,j) = (q(i,j) - q(i-1,j))/dxc(i,j)
793  ddx(i,j) = abs(ddx(i,j))
794  enddo
795  enddo
796  do j=js,je+1
797  do i=is,ie
798  ddy(i,j) = (q(i,j) - q(i,j-1))/dyc(i,j)
799  ddy(i,j) = abs(ddy(i,j))
800  enddo
801  enddo
802  do j=js,je
803  do i=is,ie
804  a3(i,j) = max( ddx(i,j), ddx(i+1,j), ddy(i,j), ddy(i,j+1) )
805  enddo
806  enddo
807  call global_mx(a3, 0, smin, smax, bd)
808  if ( is_master() ) write(*,*) 'After filter: Max_slope=', smax
809  endif
810 
811  end subroutine two_delta_filter
812 
813 
814 
815  subroutine del2_cubed_sphere(npx, npy, q, area, dx, dy, dxc, dyc, sin_sg, nmax, cd, zero_ocean, oro, nested, domain, bd)
816  type(fv_grid_bounds_type), intent(IN) :: bd
817  integer, intent(in):: npx, npy
818  integer, intent(in):: nmax
819  real(kind=R_GRID), intent(in):: cd
820  logical, intent(in):: zero_ocean
821  ! INPUT arrays
822  real(kind=R_GRID), intent(in)::area(bd%isd:bd%ied, bd%jsd:bd%jed)
823  real, intent(in):: dx(bd%isd:bd%ied, bd%jsd:bd%jed+1)
824  real, intent(in):: dy(bd%isd:bd%ied+1,bd%jsd:bd%jed)
825  real, intent(in):: dxc(bd%isd:bd%ied+1,bd%jsd:bd%jed)
826  real, intent(in):: dyc(bd%isd:bd%ied, bd%jsd:bd%jed+1)
827  real, intent(IN):: sin_sg(bd%isd:bd%ied,bd%jsd:bd%jed,9)
828  real, intent(in):: oro(bd%isd:bd%ied, bd%jsd:bd%jed) ! 0==water, 1==land
829  logical, intent(IN) :: nested
830  type(domain2d), intent(INOUT) :: domain
831  ! OUTPUT arrays
832  real, intent(inout):: q(bd%is-ng:bd%ie+ng, bd%js-ng:bd%je+ng)
833 ! Local:
834  real ddx(bd%is:bd%ie+1,bd%js:bd%je), ddy(bd%is:bd%ie,bd%js:bd%je+1)
835  integer i,j,n
836 
837  integer :: is, ie, js, je
838  integer :: isd, ied, jsd, jed
839 
840  is = bd%is
841  ie = bd%ie
842  js = bd%js
843  je = bd%je
844  isd = bd%isd
845  ied = bd%ied
846  jsd = bd%jsd
847  jed = bd%jed
848 
849 
850  call mpp_update_domains(q,domain,whalo=ng,ehalo=ng,shalo=ng,nhalo=ng)
851 
852 ! First step: average the corners:
853  if ( is==1 .and. js==1 .and. .not. nested) then
854  q(1,1) = (q(1,1)*area(1,1)+q(0,1)*area(0,1)+q(1,0)*area(1,0)) &
855  / ( area(1,1)+ area(0,1)+ area(1,0) )
856  q(0,1) = q(1,1)
857  q(1,0) = q(1,1)
858  endif
859  if ( (ie+1)==npx .and. js==1 .and. .not. nested) then
860  q(ie, 1) = (q(ie,1)*area(ie,1)+q(npx,1)*area(npx,1)+q(ie,0)*area(ie,0)) &
861  / ( area(ie,1)+ area(npx,1)+ area(ie,0))
862  q(npx,1) = q(ie,1)
863  q(ie, 0) = q(ie,1)
864  endif
865  if ( (ie+1)==npx .and. (je+1)==npy .and. .not. nested ) then
866  q(ie, je) = (q(ie,je)*area(ie,je)+q(npx,je)*area(npx,je)+q(ie,npy)*area(ie,npy)) &
867  / ( area(ie,je)+ area(npx,je)+ area(ie,npy))
868  q(npx,je) = q(ie,je)
869  q(ie,npy) = q(ie,je)
870  endif
871  if ( is==1 .and. (je+1)==npy .and. .not. nested) then
872  q(1, je) = (q(1,je)*area(1,je)+q(0,je)*area(0,je)+q(1,npy)*area(1,npy)) &
873  / ( area(1,je)+ area(0,je)+ area(1,npy))
874  q(0, je) = q(1,je)
875  q(1,npy) = q(1,je)
876  endif
877 
878 
879  do n=1,nmax
880  if( n>1 ) call mpp_update_domains(q,domain,whalo=ng,ehalo=ng,shalo=ng,nhalo=ng)
881  do j=js,je
882  do i=is,ie+1
883  ddx(i,j) = 0.5*(sin_sg(i-1,j,3)+sin_sg(i,j,1))*dy(i,j)*(q(i-1,j)-q(i,j))/dxc(i,j)
884  enddo
885  enddo
886  do j=js,je+1
887  do i=is,ie
888  ddy(i,j) = dx(i,j)*(q(i,j-1)-q(i,j))/dyc(i,j) &
889  *0.5*(sin_sg(i,j-1,4)+sin_sg(i,j,2))
890  enddo
891  enddo
892 
893  if ( zero_ocean ) then
894 ! Limit diffusive flux over ater cells:
895  do j=js,je
896  do i=is,ie+1
897  ddx(i,j) = max(0., min(oro(i-1,j), oro(i,j))) * ddx(i,j)
898  enddo
899  enddo
900  do j=js,je+1
901  do i=is,ie
902  ddy(i,j) = max(0., min(oro(i,j-1), oro(i,j))) * ddy(i,j)
903  enddo
904  enddo
905  endif
906 
907  do j=js,je
908  do i=is,ie
909  q(i,j) = q(i,j) + cd/area(i,j)*(ddx(i,j)-ddx(i+1,j)+ddy(i,j)-ddy(i,j+1))
910  enddo
911  enddo
912  enddo
913 
914  end subroutine del2_cubed_sphere
915 
916 
917  subroutine del4_cubed_sphere(npx, npy, q, area, dx, dy, dxc, dyc, sin_sg, nmax, zero_ocean, oro, nested, domain, bd)
918  type(fv_grid_bounds_type), intent(IN) :: bd
919  integer, intent(in):: npx, npy, nmax
920  logical, intent(in):: zero_ocean
921  real, intent(in):: oro(bd%isd:bd%ied, bd%jsd:bd%jed) ! 0==water, 1==land
922  real(kind=R_GRID), intent(in)::area(bd%isd:bd%ied, bd%jsd:bd%jed)
923  real, intent(in):: dx(bd%isd:bd%ied, bd%jsd:bd%jed+1)
924  real, intent(in):: dy(bd%isd:bd%ied+1,bd%jsd:bd%jed)
925  real, intent(in):: dxc(bd%isd:bd%ied+1,bd%jsd:bd%jed)
926  real, intent(in):: dyc(bd%isd:bd%ied, bd%jsd:bd%jed+1)
927  real, intent(IN):: sin_sg(bd%isd:bd%ied,bd%jsd:bd%jed,9)
928  real, intent(inout):: q(bd%is-ng:bd%ie+ng, bd%js-ng:bd%je+ng)
929  logical, intent(IN) :: nested
930  type(domain2d), intent(INOUT) :: domain
931 ! diffusivity
932  real :: diff(bd%is-3:bd%ie+2,bd%js-3:bd%je+2)
933 ! diffusive fluxes:
934  real :: fx1(bd%is:bd%ie+1,bd%js:bd%je), fy1(bd%is:bd%ie,bd%js:bd%je+1)
935  real :: fx2(bd%is:bd%ie+1,bd%js:bd%je), fy2(bd%is:bd%ie,bd%js:bd%je+1)
936  real :: fx4(bd%is:bd%ie+1,bd%js:bd%je), fy4(bd%is:bd%ie,bd%js:bd%je+1)
937  real, dimension(bd%isd:bd%ied,bd%jsd:bd%jed):: d2, win, wou
938  real, dimension(bd%is:bd%ie,bd%js:bd%je):: qlow, qmin, qmax, q0
939  real, parameter:: esl = 1.e-20
940  integer i,j, n
941 
942  integer :: is, ie, js, je
943  integer :: isd, ied, jsd, jed
944 
945  is = bd%is
946  ie = bd%ie
947  js = bd%js
948  je = bd%je
949  isd = bd%isd
950  ied = bd%ied
951  jsd = bd%jsd
952  jed = bd%jed
953 
954  !On a nested grid the haloes are not filled. Set to zero.
955  d2 = 0.
956  win = 0.
957  wou = 0.
958 
959  do j=js-1,je+1
960  do i=is-1,ie+1
961  diff(i,j) = cd4*area(i,j) ! area dependency is needed for stretched grid
962  enddo
963  enddo
964 
965  do j=js,je
966  do i=is,ie
967  q0(i,j) = q(i,j)
968  enddo
969  enddo
970 
971  do n=1,nmax
972  call mpp_update_domains(q,domain)
973 
974 ! First step: average the corners:
975  if ( is==1 .and. js==1 .and. .not. nested) then
976  q(1,1) = (q(1,1)*area(1,1)+q(0,1)*area(0,1)+q(1,0)*area(1,0)) &
977  / ( area(1,1)+ area(0,1)+ area(1,0) )
978  q(0,1) = q(1,1)
979  q(1,0) = q(1,1)
980  q(0,0) = q(1,1)
981  endif
982  if ( (ie+1)==npx .and. js==1 .and. .not. nested) then
983  q(ie, 1) = (q(ie,1)*area(ie,1)+q(npx,1)*area(npx,1)+q(ie,0)*area(ie,0)) &
984  / ( area(ie,1)+ area(npx,1)+ area(ie,0))
985  q(npx,1) = q(ie,1)
986  q(ie, 0) = q(ie,1)
987  q(npx,0) = q(ie,1)
988  endif
989  if ( (ie+1)==npx .and. (je+1)==npy .and. .not. nested) then
990  q(ie, je) = (q(ie,je)*area(ie,je)+q(npx,je)*area(npx,je)+q(ie,npy)*area(ie,npy)) &
991  / ( area(ie,je)+ area(npx,je)+ area(ie,npy))
992  q(npx, je) = q(ie,je)
993  q(ie, npy) = q(ie,je)
994  q(npx,npy) = q(ie,je)
995  endif
996  if ( is==1 .and. (je+1)==npy .and. .not. nested) then
997  q(1, je) = (q(1,je)*area(1,je)+q(0,je)*area(0,je)+q(1,npy)*area(1,npy)) &
998  / ( area(1,je)+ area(0,je)+ area(1,npy))
999  q(0, je) = q(1,je)
1000  q(1,npy) = q(1,je)
1001  q(0,npy) = q(1,je)
1002  endif
1003 
1004  do j=js,je
1005  do i=is,ie
1006  qmin(i,j) = min(q0(i,j), q(i-1,j-1), q(i,j-1), q(i+1,j-1), &
1007  q(i-1,j ), q(i,j ), q(i+1,j ), &
1008  q(i-1,j+1), q(i,j+1), q(i+1,j+1) )
1009  qmax(i,j) = max(peak_fac*q0(i,j), q(i-1,j-1), q(i,j-1), q(i+1,j-1), &
1010  q(i-1,j ), q(i,j ), q(i+1,j ), &
1011  q(i-1,j+1), q(i,j+1), q(i+1,j+1) )
1012  enddo
1013  enddo
1014 
1015 !--------------
1016 ! Compute del-2
1017 !--------------
1018 ! call copy_corners(q, npx, npy, 1)
1019  do j=js,je
1020  do i=is,ie+1
1021  fx2(i,j) = 0.25*(diff(i-1,j)+diff(i,j))*dy(i,j)*(q(i-1,j)-q(i,j))/dxc(i,j) &
1022  *(sin_sg(i,j,1)+sin_sg(i-1,j,3))
1023  enddo
1024  enddo
1025 
1026 ! call copy_corners(q, npx, npy, 2)
1027  do j=js,je+1
1028  do i=is,ie
1029  fy2(i,j) = 0.25*(diff(i,j-1)+diff(i,j))*dx(i,j)*(q(i,j-1)-q(i,j))/dyc(i,j) &
1030  *(sin_sg(i,j,2)+sin_sg(i,j-1,4))
1031  enddo
1032  enddo
1033 
1034  do j=js,je
1035  do i=is,ie
1036  d2(i,j) = (fx2(i,j)-fx2(i+1,j)+fy2(i,j)-fy2(i,j+1)) / area(i,j)
1037  enddo
1038  enddo
1039 
1040 ! qlow == low order monotonic solution
1041  if ( zero_ocean ) then
1042 ! Limit diffusive flux over ater cells:
1043  do j=js,je
1044  do i=is,ie+1
1045  fx1(i,j) = max(0., min(oro(i-1,j), oro(i,j))) * fx2(i,j)
1046  enddo
1047  enddo
1048  do j=js,je+1
1049  do i=is,ie
1050  fy1(i,j) = max(0., min(oro(i,j-1), oro(i,j))) * fy2(i,j)
1051  enddo
1052  enddo
1053  do j=js,je
1054  do i=is,ie
1055  qlow(i,j) = q(i,j) + (fx1(i,j)-fx1(i+1,j)+fy1(i,j)-fy1(i,j+1)) / area(i,j)
1056  d2(i,j) = diff(i,j) * d2(i,j)
1057  enddo
1058  enddo
1059  else
1060  do j=js,je
1061  do i=is,ie
1062  qlow(i,j) = q(i,j) + d2(i,j)
1063  d2(i,j) = diff(i,j) * d2(i,j)
1064  enddo
1065  enddo
1066  endif
1067 
1068  call mpp_update_domains(d2,domain)
1069 
1070 !---------------------
1071 ! Compute del4 fluxes:
1072 !---------------------
1073 ! call copy_corners(d2, npx, npy, 1)
1074  do j=js,je
1075  do i=is,ie+1
1076  fx4(i,j) = 0.5*(sin_sg(i-1,j,3)+sin_sg(i,j,1))*dy(i,j)*(d2(i,j)-d2(i-1,j))/dxc(i,j)-fx2(i,j)
1077  enddo
1078  enddo
1079 
1080 ! call copy_corners(d2, npx, npy, 2)
1081  do j=js,je+1
1082  do i=is,ie
1083  fy4(i,j) = dx(i,j)*(d2(i,j)-d2(i,j-1))/dyc(i,j) &
1084  *0.5*(sin_sg(i,j,2)+sin_sg(i,j-1,4))-fy2(i,j)
1085  enddo
1086  enddo
1087 
1088 !----------------
1089 ! Flux limitting:
1090 !----------------
1091  do j=js,je
1092  do i=is,ie
1093  win(i,j) = max(0.,fx4(i, j)) - min(0.,fx4(i+1,j)) + &
1094  max(0.,fy4(i, j)) - min(0.,fy4(i,j+1)) + esl
1095  wou(i,j) = max(0.,fx4(i+1,j)) - min(0.,fx4(i, j)) + &
1096  max(0.,fy4(i,j+1)) - min(0.,fy4(i, j)) + esl
1097  win(i,j) = max(0., qmax(i,j) - qlow(i,j)) / win(i,j)*area(i,j)
1098  wou(i,j) = max(0., qlow(i,j) - qmin(i,j)) / wou(i,j)*area(i,j)
1099  enddo
1100  enddo
1101 
1102  call mpp_update_domains(win,domain, complete=.false.)
1103  call mpp_update_domains(wou,domain, complete=.true.)
1104 
1105  do j=js,je
1106  do i=is,ie+1
1107  if ( fx4(i,j) > 0. ) then
1108  fx4(i,j) = min(1., wou(i-1,j), win(i,j)) * fx4(i,j)
1109  else
1110  fx4(i,j) = min(1., win(i-1,j), wou(i,j)) * fx4(i,j)
1111  endif
1112  enddo
1113  enddo
1114  do j=js,je+1
1115  do i=is,ie
1116  if ( fy4(i,j) > 0. ) then
1117  fy4(i,j) = min(1., wou(i,j-1), win(i,j)) * fy4(i,j)
1118  else
1119  fy4(i,j) = min(1., win(i,j-1), wou(i,j)) * fy4(i,j)
1120  endif
1121  enddo
1122  enddo
1123 
1124 
1125  if ( zero_ocean ) then
1126 ! Limit diffusive flux over ocean cells:
1127  do j=js,je
1128  do i=is,ie+1
1129  fx4(i,j) = max(0., min(oro(i-1,j), oro(i,j))) * fx4(i,j)
1130  enddo
1131  enddo
1132  do j=js,je+1
1133  do i=is,ie
1134  fy4(i,j) = max(0., min(oro(i,j-1), oro(i,j))) * fy4(i,j)
1135  enddo
1136  enddo
1137  endif
1138 
1139 ! Update:
1140  do j=js,je
1141  do i=is,ie
1142  q(i,j) = qlow(i,j) + (fx4(i,j)-fx4(i+1,j)+fy4(i,j)-fy4(i,j+1))/area(i,j)
1143  enddo
1144  enddo
1145 
1146  enddo ! end n-loop
1147 
1148  end subroutine del4_cubed_sphere
1149 
1150 
1151  subroutine map_to_cubed_raw(igh, im, jt, lat1, lon1, zs, ft, grid, agrid, &
1152  q2, f2, h2, npx, npy, jstart, jend, stretch_fac, &
1153  nested, npx_global, bd)
1155 ! Input
1156  type(fv_grid_bounds_type), intent(IN) :: bd
1157  integer, intent(in):: igh, im, jt
1158  integer, intent(in):: npx, npy, npx_global
1159  real, intent(in):: lat1(jt+1) ! original southern edge of the cell [-pi/2:pi/2]
1160  real, intent(in):: lon1(im+1) ! original western edge of the cell [0:2*pi]
1161  real(kind=4), intent(in), dimension(-igh:im+igh,jt):: zs, ft
1162  real(kind=R_GRID), intent(in):: grid(bd%isd:bd%ied+1, bd%jsd:bd%jed+1,2)
1163  real(kind=R_GRID), intent(in):: agrid(bd%isd:bd%ied, bd%jsd:bd%jed, 2)
1164  integer, intent(in):: jstart, jend
1165  real(kind=R_GRID), intent(IN) :: stretch_fac
1166  logical, intent(IN) :: nested
1167 ! Output
1168  real, intent(out):: q2(bd%isd:bd%ied,bd%jsd:bd%jed) ! Mapped data at the target resolution
1169  real, intent(out):: f2(bd%isd:bd%ied,bd%jsd:bd%jed) ! oro
1170  real, intent(out):: h2(bd%isd:bd%ied,bd%jsd:bd%jed) ! variances of terrain
1171 ! Local
1172  real :: lon_g(-igh:im+igh)
1173  real lat_g(jt), cos_g(jt)
1174  real(kind=R_GRID) e2(2)
1175  real(kind=R_GRID) grid3(3, bd%isd:bd%ied+1, bd%jsd:bd%jed+1)
1176  real(kind=R_GRID), dimension(3):: p1, p2, p3, p4, pc, pp
1177  real(kind=R_GRID), dimension(3):: vp_12, vp_23, vp_34, vp_14
1178  integer i,j, np, k
1179  integer ii, jj, i1, i2, j1, j2, min_pts
1180  real th1, aa, asum, qsum, fsum, hsum, lon_w, lon_e, lat_s, lat_n, r2d
1181  real qsp, fsp, hsp
1182  real qnp, fnp, hnp
1183  real delg, th0, tmp1, prod1, prod2, prod3, prod4
1184  integer ig_lon, jp
1185  integer:: lat_crit
1186  real pi5, pi2
1187 
1188  integer :: is, ie, js, je
1189  integer :: isd, ied, jsd, jed
1190 
1191  is = bd%is
1192  ie = bd%ie
1193  js = bd%js
1194  je = bd%je
1195  isd = bd%isd
1196  ied = bd%ied
1197  jsd = bd%jsd
1198  jed = bd%jed
1199 
1200  pi2 = pi + pi
1201  pi5 = 0.5 * pi
1202  r2d = 180./pi
1203 
1204  do i=1,im
1205  lon_g(i) = 0.5*(lon1(i)+lon1(i+1))
1206  enddo
1207  do i=-igh,0
1208  lon_g(i) = lon_g(i+im)
1209  enddo
1210  do i=im+1,im+igh
1211  lon_g(i) = lon_g(i-im)
1212  enddo
1213 
1214  do j=1,jt
1215  lat_g(j) = 0.5*(lat1(j)+lat1(j+1))
1216  cos_g(j) = cos( lat_g(j) )
1217  enddo
1218 
1219  do j=jsd,jed+1
1220  do i=isd,ied+1
1221  call latlon2xyz(real(grid(i,j,1:2),kind=R_GRID), grid3(1,i,j))
1222  enddo
1223  enddo
1224 
1225  if(is_master()) write(*,*) 'surf_map: Search started ....'
1226 
1227 ! stretch_fac * pi5/(npx-1) / (pi/nlat)
1228  lat_crit = nint( stretch_fac*real(nlat)/real(npx_global-1) )
1229  lat_crit = min( jt, max( 4, lat_crit ) )
1230 
1231  if ( jstart==1 ) then
1232  write(*,*) mpp_pe(), 'lat_crit=', r2d*lat_g(lat_crit)
1233  elseif ( jend==nlat ) then
1234 ! write(*,*) mpp_pe(), 'lat_crit=', r2d*lat_g(jt-lat_crit+1)
1235  endif
1236 
1237 !----
1238 ! SP:
1239 !----
1240  iF ( jstart == 1 ) then
1241  asum = 0.
1242  qsum = 0.
1243  fsum = 0.
1244  hsum = 0.
1245  do j=1,lat_crit
1246  aa = cos_g(j)
1247  do i=1,im
1248  asum = asum + aa
1249  qsum = qsum + zs(i,j)*aa
1250  fsum = fsum + ft(i,j)*aa
1251  enddo
1252  enddo
1253  qsp = qsum / asum
1254  fsp = fsum / asum
1255  hsum = 0.
1256  np = 0
1257  do j=1,lat_crit
1258  do i=1,im
1259  np = np + 1
1260  hsum = hsum + (qsp-zs(i,j))**2
1261  enddo
1262  enddo
1263  hsp = hsum / real(np)
1264 ! write(*,*) 'SP GID, zs_sp, f_sp, hsp=', mpp_pe(), qsp, fsp, hsp
1265  endif
1266 !----
1267 ! NP:
1268 !----
1269  iF ( jend == nlat ) then
1270  asum = 0.
1271  qsum = 0.
1272  fsum = 0.
1273  hsum = 0.
1274  do jp=jend-lat_crit+1, jend
1275  j = jp - jstart + 1
1276  aa = cos_g(j)
1277  do i=1,im
1278  asum = asum + aa
1279  qsum = qsum + zs(i,j)*aa
1280  fsum = fsum + ft(i,j)*aa
1281  enddo
1282  enddo
1283  qnp = qsum / asum
1284  fnp = fsum / asum
1285  hsum = 0.
1286  np = 0
1287  do jp=jend-lat_crit+1, jend
1288  j = jp - jstart + 1
1289  do i=1,im
1290  np = np + 1
1291  hsum = hsum + (qnp-zs(i,j))**2
1292  enddo
1293  enddo
1294  hnp = hsum / real(np)
1295 ! write(*,*) 'NP GID, zs_np, f_np, hnp=', mpp_pe(), qnp, fnp, hnp
1296  endif
1297 
1298  min_pts = 999999
1299  do j=jsd,jed
1300  do i=isd,ied
1301 
1302  !Do not go into corners
1303  if (((i < is .and. j < js) .or. &
1304  (i < is .and. j > je) .or. &
1305  (i > ie .and. j < js) .or. &
1306  (i > ie .and. j > je)) .and. .not. nested) then
1307  q2(i,j) = 1.e25
1308  f2(i,j) = 1.e25
1309  h2(i,j) = 1.e25
1310  goto 4444
1311  end if
1312 
1313  if ( agrid(i,j,2) < -pi5+stretch_fac*pi5/real(npx_global-1) ) then
1314 ! SP:
1315  q2(i,j) = qsp
1316  f2(i,j) = fsp
1317  h2(i,j) = hsp
1318  goto 4444
1319  elseif ( agrid(i,j,2) > pi5-stretch_fac*pi5/real(npx_global-1) ) then
1320 ! NP:
1321  q2(i,j) = qnp
1322  f2(i,j) = fnp
1323  h2(i,j) = hnp
1324  goto 4444
1325  endif
1326 
1327  lat_s = min( grid(i,j,2), grid(i+1,j,2), grid(i,j+1,2), grid(i+1,j+1,2), agrid(i,j,2) )
1328  lat_n = max( grid(i,j,2), grid(i+1,j,2), grid(i,j+1,2), grid(i+1,j+1,2), agrid(i,j,2) )
1329 ! Enlarge the search zone:
1330  delg = max( 0.2*(lat_n-lat_s), pi/real(npx_global-1), pi2/real(nlat) )
1331  lat_s = max(-pi5, lat_s - delg)
1332  lat_n = min( pi5, lat_n + delg)
1333 
1334  j1 = nint( (pi5+lat_s)/(pi/real(nlat)) ) - 1
1335  if ( lat_s*r2d < (-90.+90./real(npx_global-1)) ) j1 = 1
1336  j1 = max(jstart, j1)
1337 
1338  j2 = nint( (pi5+lat_n)/(pi/real(nlat)) ) + 1
1339  if ( lat_n*r2d > (90.-90./real(npx_global-1)) ) j2 = nlat
1340  j2 = min(jend, j2)
1341 
1342  j1 = j1 - jstart + 1
1343  j2 = j2 - jstart + 1
1344 
1345  lon_w = min( grid(i,j,1), grid(i+1,j,1), grid(i,j+1,1), grid(i+1,j+1,1) )
1346  lon_e = max( grid(i,j,1), grid(i+1,j,1), grid(i,j+1,1), grid(i+1,j+1,1) )
1347 
1348  if ( (lon_e-lon_w) > pi ) then
1349  i1 = -nint( (pi2-lon_e)/pi2 * real(im) ) - 1
1350  i2 = nint( lon_w/pi2 * real(im) ) + 1
1351  else
1352  i1 = nint( lon_w/pi2 * real(im) ) - 1
1353  i2 = nint( lon_e/pi2 * real(im) ) + 1
1354  endif
1355 
1356 ! Enlarge the search zone:
1357  ig_lon = max(1, (i2-i1)/8)
1358  i1 = max( -igh, i1 - ig_lon)
1359  i2 = min(im+igh, i2 + ig_lon)
1360 
1361  np = 0
1362  qsum = 0.
1363  fsum = 0.
1364  hsum = 0.
1365  asum = 0.
1366 !
1367 ! 4----------3
1368 ! / /
1369 ! / pp /
1370 ! / /
1371 ! 1----------2
1372 !
1373  do k=1,3
1374  p1(k) = grid3(k,i, j)
1375  p2(k) = grid3(k,i+1,j)
1376  p3(k) = grid3(k,i+1,j+1)
1377  p4(k) = grid3(k,i,j+1)
1378  pc(k) = p1(k) + p2(k) + p3(k) + p4(k)
1379  enddo
1380  call normalize_vect( pc )
1381 
1382  th0 = min( v_prod(p1,p3), v_prod(p2, p4) )
1383  th1 = min( cos_grid, cos(0.25*acos(max(v_prod(p1,p3), v_prod(p2, p4)))))
1384 
1385  call vect_cross(vp_12, p1, p2)
1386  call vect_cross(vp_23, p2, p3)
1387  call vect_cross(vp_34, p3, p4)
1388  call vect_cross(vp_14, p1, p4)
1389 
1390  prod1 = v_prod(p3, vp_12)
1391  prod2 = v_prod(p1, vp_23)
1392  prod3 = v_prod(p1, vp_34)
1393  prod4 = v_prod(p2, vp_14)
1394 
1395  do jj=j1,j2
1396  aa = cos_g(jj)
1397  e2(2) = lat_g(jj)
1398  do ii=i1,i2
1399  e2(1) = lon_g(ii)
1400  call latlon2xyz(e2, pp)
1401 ! Check two extrems:
1402  tmp1 = v_prod(pp, pc)
1403 ! Points that are close to center:
1404  if ( tmp1 > th1 ) goto 1111 ! inside
1405 ! check to exclude points too far away:
1406  if ( tmp1 < th0 ) goto 2222 ! outside
1407 ! Check if inside the polygon
1408  if ( v_prod(pp, vp_12)*prod1 < 0. ) goto 2222
1409  if ( v_prod(pp, vp_23)*prod2 < 0. ) goto 2222
1410  if ( v_prod(pp, vp_34)*prod3 < 0. ) goto 2222
1411  if ( v_prod(pp, vp_14)*prod4 < 0. ) goto 2222
1412 1111 np = np + 1
1413  qsum = qsum + zs(ii,jj)*aa
1414  fsum = fsum + ft(ii,jj)*aa
1415  hsum = hsum + zs(ii,jj)**2
1416  asum = asum + aa
1417 2222 continue
1418  enddo
1419  enddo
1420 
1421  if ( np > 0 ) then
1422  q2(i,j) = qsum / asum
1423  f2(i,j) = fsum / asum
1424  h2(i,j) = hsum / real(np) - q2(i,j)**2
1425  min_pts = min(min_pts, np)
1426  else
1427  write(*,*) 'min and max lat_g is ', r2d*minval(lat_g), r2d*maxval(lat_g), mpp_pe()
1428  write(*,*) 'Surf_map failed for GID=', mpp_pe(), i,j, '(lon,lat)=', r2d*agrid(i,j,1),r2d*agrid(i,j,2)
1429  write(*,*) '[jstart, jend]', jstart, jend
1430  call mpp_error(fatal,'Surf_map failed')
1431  endif
1432 4444 continue
1433  enddo
1434  enddo
1435 
1436  if(is_master()) write(*,*) 'surf_map: minimum pts per cell (master PE)=', min_pts
1437  if ( min_pts <3 ) then
1438  if(is_master()) write(*,*) 'Warning: too few points used in creating the cell mean terrain !!!'
1439  endif
1440 
1441  end subroutine map_to_cubed_raw
1442 
1443 
1444 #ifdef JUNK
1445  logical function inside_p4(p1, p2, p3, p4, pp, th0)
1447  real, intent(in):: p1(3), p2(3), p3(3), p4(3)
1448  real, intent(in):: pp(3)
1449  real, intent(in):: th0
1450 ! A * B = |A| |B| cos(angle)
1451 ! Local:
1452  real(kind=R_GRID) v1(3), v2(3), vp(3)
1453  real(kind=R_GRID) tmp1
1454  integer k
1455 
1456  inside_p4 = .false.
1457 
1458 ! 1st check: to exclude points too far away:
1459  do k=1,3
1460  vp(k) = p1(k) + p2(k) + p3(k) + p4(k)
1461  enddo
1462  call normalize_vect( vp )
1463 
1464  tmp1 = v_prod(pp, vp)
1465  if ( tmp1 < th0 ) then ! outside
1466  return
1467  endif
1468 
1469 ! 1st segment: 1---2
1470  call vect_cross(vp, p1, p2)
1471  if ( v_prod(pp, vp)*v_prod(p3, vp) < 0. ) return
1472 ! 2nd segment: 2---3
1473  call vect_cross(vp, p2, p3)
1474  if ( v_prod(pp, vp)*v_prod(p1, vp) < 0. ) return
1475 ! 3rd segment: 3---4
1476  call vect_cross(vp, p3, p4)
1477  if ( v_prod(pp, vp)*v_prod(p1, vp) < 0. ) return
1478 ! 4th segment: 1---4
1479  call vect_cross(vp, p1, p4)
1480  if ( v_prod(pp, vp)*v_prod(p2, vp) < 0. ) return
1481 
1482  inside_p4 = .true.
1483 
1484  end function inside_p4
1485 #endif
1486 
1487  subroutine handle_err(status)
1488 #include <netcdf.inc>
1489  integer status
1490 
1491  if (status .ne. nf_noerr) then
1492  print *, nf_strerror(status)
1493  stop 'Stopped'
1494  endif
1495 
1496  end subroutine handle_err
1497 
1498  subroutine remove_ice_sheets (lon, lat, lfrac, bd )
1499 !---------------------------------
1500 ! Bruce Wyman's fix for Antarctic
1501 !---------------------------------
1502  type(fv_grid_bounds_type), intent(IN) :: bd
1503  real(kind=R_GRID), intent(in) :: lon(bd%isd:bd%ied,bd%jsd:bd%jed), lat(bd%isd:bd%ied,bd%jsd:bd%jed)
1504  real, intent(inout) :: lfrac(bd%isd:bd%ied,bd%jsd:bd%jed)
1505 
1506  integer :: is, ie, js, je
1507  integer :: isd, ied, jsd, jed
1508 
1509 ! lon = longitude in radians
1510 ! lat = latitude in radians
1511 ! lfrac = land-sea mask (land=1, sea=0)
1512 
1513  integer :: i, j
1514  real :: dtr, phs, phn
1515 
1516  is = bd%is
1517  ie = bd%ie
1518  js = bd%js
1519  je = bd%je
1520  isd = bd%isd
1521  ied = bd%ied
1522  jsd = bd%jsd
1523  jed = bd%jed
1524 
1525  dtr = acos(0.)/90.
1526  phs = -83.9999*dtr
1527 ! phn = -78.9999*dtr
1528  phn = -76.4*dtr
1529 
1530  do j = jsd, jed
1531  do i = isd, ied
1532  if ( lat(i,j) < phn ) then
1533  ! replace all below this latitude
1534  if ( lat(i,j) < phs ) then
1535  lfrac(i,j) = 1.0
1536  cycle
1537  endif
1538  ! replace between 270 and 360 deg
1539  if ( sin(lon(i,j)) < 0. .and. cos(lon(i,j)) > 0.) then
1540  lfrac(i,j) = 1.0
1541  cycle
1542  endif
1543  endif
1544  enddo
1545  enddo
1546  end subroutine remove_ice_sheets
1547 
1548 
1549 !#######################################################################
1550 ! reads the namelist file, write namelist to log file,
1551 ! and initializes constants
1552 
1553 subroutine read_namelist
1555  integer :: unit, ierr, io
1556 ! real :: dtr, ght
1557 
1558 ! read namelist
1559 
1560  if (namelist_read) return
1561 
1562 #ifdef INTERNAL_FILE_NML
1563  read (input_nml_file, nml=surf_map_nml, iostat=io)
1564  ierr = check_nml_error(io,'surf_map_nml')
1565 #else
1566  unit = open_namelist_file( )
1567  ierr=1
1568  do while (ierr /= 0)
1569  read (unit, nml=surf_map_nml, iostat=io, end=10)
1570  ierr = check_nml_error(io,'surf_map_nml')
1571  enddo
1572  10 call close_file (unit)
1573 #endif
1574 
1575 ! write version and namelist to log file
1576 
1577  if (mpp_pe() == mpp_root_pe()) then
1578  unit = stdlog()
1579  write (unit, nml=surf_map_nml)
1580  endif
1581 
1582  namelist_read = .true.
1583 
1584 end subroutine read_namelist
1585 
1586 subroutine zonal_mean(im, p, zmean)
1587 ! replace p with its zonal mean
1588  integer, intent(in):: im
1589  real(kind=4), intent(inout):: p(im)
1590  real, intent(out):: zmean
1591  integer i
1592 
1593  zmean = 0.
1594  do i=1,im
1595  zmean = zmean + p(i)
1596  enddo
1597  zmean = zmean / real(im)
1598  do i=1,im
1599  p(i) = zmean
1600  enddo
1601 end subroutine zonal_mean
1602 
1603 
1604  end module fv_surf_map_nlm_mod
Definition: fms.F90:20
real, parameter, public radius
Radius of the Earth [m].
Definition: constants.F90:72
subroutine, public del4_cubed_sphere(npx, npy, q, area, dx, dy, dxc, dyc, sin_sg, nmax, zero_ocean, oro, nested, domain, bd)
real, dimension(:,:), allocatable, public zs_g
character(len=128) surf_file
real(kind=8), parameter, public pi_8
Ratio of circle circumference to diameter [N/A].
Definition: constants.F90:73
real(kind=r_grid) da_min
character(len=6) surf_format
real(kind=r_grid) function, public v_prod(v1, v2)
subroutine, public global_mx(q, n_g, qmin, qmax, bd)
subroutine two_delta_filter(npx, npy, q, area, dx, dy, dxa, dya, dxc, dyc, sin_sg, cd, zero_ocean, check_slope, filter_type, oro, nested, domain, bd, ntmax)
void vect_cross(const double *p1, const double *p2, double *e)
Definition: mosaic_util.c:648
real, dimension(:,:), allocatable, public oro_g
Definition: mpp.F90:39
integer function, public check_nml_error(IOSTAT, NML_NAME)
Definition: fms.F90:658
character(len=input_str_length), dimension(:), allocatable, target, public input_nml_file
Definition: mpp.F90:1378
subroutine handle_err(status)
integer, parameter, public ng
subroutine timing_on(blk_name)
integer, parameter, public r_grid
subroutine, public surfdrv(npx, npy, grid, agrid, area, dx, dy, dxa, dya, dxc, dyc, sin_sg, phis, stretch_fac, nested, npx_global, domain, grid_number, bd)
real function, public great_circle_dist(q1, q2, radius)
subroutine remove_ice_sheets(lon, lat, lfrac, bd)
logical function inside_p4(p1, p2, p3, p4, pp, th0)
subroutine zonal_mean(im, p, zmean)
real, parameter, public grav
Acceleration due to gravity [m/s^2].
Definition: constants.F90:76
real function, public g_sum(domain, p, ifirst, ilast, jfirst, jlast, ngc, area, mode, reproduce)
void normalize_vect(double *e)
Definition: mosaic_util.c:679
subroutine, public fv3_zs_filter(bd, isd, ied, jsd, jed, npx, npy, npx_global, stretch_fac, nested, domain, area, dxa, dya, dx, dy, dxc, dyc, grid, agrid, sin_sg, phis, oro)
#define max(a, b)
Definition: mosaic_util.h:33
real, dimension(:,:), allocatable, public sgh_g
subroutine, public del2_cubed_sphere(npx, npy, q, area, dx, dy, dxc, dyc, sin_sg, nmax, cd, zero_ocean, oro, nested, domain, bd)
#define min(a, b)
Definition: mosaic_util.h:32
void latlon2xyz(int size, const double *lon, const double *lat, double *x, double *y, double *z)
Definition: mosaic_util.c:211
subroutine, public error_mesg(routine, message, level)
Definition: fms.F90:529
subroutine map_to_cubed_raw(igh, im, jt, lat1, lon1, zs, ft, grid, agrid, q2, f2, h2, npx, npy, jstart, jend, stretch_fac, nested, npx_global, bd)
real(fp), parameter, public pi
subroutine timing_off(blk_name)
character(len=3) grid_string