cpld_gridgen  1.13.0
All Data Structures Files Functions Variables Pages
postwgts.F90
Go to the documentation of this file.
1 
7 module postwgts
8 
9  use esmf
10 
11  use gengrid_kinds, only : cl,cm,cs
12  use grdvars, only : nv
13  use charstrings, only : dirout, res, logmsg
14  use netcdf
15 
16  implicit none
17 
18 contains
22 
23  subroutine make_postwgts
24 
25  ! local variables
26  character(len=CL) :: fsrc, fdst, fwgt
27 
28  character(len=CM), dimension(2) :: methodname = (/'conserve', 'bilinear'/)
29 
30  type(ESMF_RouteHandle) :: rh
31  type(ESMF_RegridMethod_Flag) :: method
32  ! the number of possible destination grids depends on the source grid resolution
33  integer :: k,rc,nd,ndest
34  character(len=CS), allocatable, dimension(:) :: destgrds
35 
36  !---------------------------------------------------------------------
37  ! set the destination grids
38  !---------------------------------------------------------------------
39 
40  if(trim(res) .eq. '500')then
41  ndest = 1
42  allocate(destgrds(ndest))
43  destgrds = (/'5p00'/)
44  end if
45  if(trim(res) .eq. '100')then
46  ndest = 2
47  allocate(destgrds(ndest))
48  destgrds = (/'5p00', '1p00'/)
49  end if
50  if(trim(res) .eq. '050')then
51  ndest = 3
52  allocate(destgrds(ndest))
53  destgrds = (/'5p00', '1p00', '0p50'/)
54  end if
55  if(trim(res) .eq. '025')then
56  ndest = 4
57  allocate(destgrds(ndest))
58  destgrds = (/'5p00', '1p00', '0p50', '0p25'/)
59  end if
60 
61  !---------------------------------------------------------------------
62  ! use ESMF to create the weights from the Ct tripole to the rectilinear
63  ! grids with conservative and bilinear methods for post; the source
64  ! file is always Ct
65  !---------------------------------------------------------------------
66 
67  do nd = 1,ndest
68  fsrc = trim(dirout)//'/'//'Ct.mx'//trim(res)//'_SCRIP.nc'
69  fdst = trim(dirout)//'/rect.'//trim(destgrds(nd))//'_SCRIP.nc'
70 
71  do k = 1,size(methodname)
72  if(trim(methodname(k)) .eq. 'bilinear')method=esmf_regridmethod_bilinear
73  if(trim(methodname(k)) .eq. 'conserve')method=esmf_regridmethod_conserve
74 
75  fwgt = trim(dirout)//'/'//'tripole.mx'//trim(res)//'.Ct.to.rect.'//trim(destgrds(nd)) &
76  //'.'//trim(methodname(k))//'.nc'
77  logmsg = 'creating weight file '//trim(fwgt)
78  print '(a)',trim(logmsg)
79 
80  call esmf_regridweightgen(srcfile=trim(fsrc),dstfile=trim(fdst), &
81  weightfile=trim(fwgt), regridmethod=method, &
82  ignoredegenerate=.true., regridroutehandle=rh, &
83  unmappedaction=esmf_unmappedaction_ignore, rc=rc)
84  if (esmf_logfounderror(rctocheck=rc, msg=esmf_logerr_passthru, &
85  line=__line__, file=__file__)) call esmf_finalize(endflag=esmf_end_abort)
86  end do
87  end do
88 
89  deallocate(destgrds)
90 
91  end subroutine make_postwgts
92 end module postwgts