cpld_gridgen 1.14.0
Loading...
Searching...
No Matches
postwgts.F90
Go to the documentation of this file.
1
7module 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
18contains
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. '900')then
41 ndest = 1
42 allocate(destgrds(ndest))
43 destgrds = (/'9p00'/)
44 end if
45 if(trim(res) .eq. '500')then
46 ndest = 2
47 allocate(destgrds(ndest))
48 destgrds = (/'9p00', '5p00'/)
49 end if
50 if(trim(res) .eq. '100')then
51 ndest = 3
52 allocate(destgrds(ndest))
53 destgrds = (/'9p00', '5p00', '1p00'/)
54 end if
55 if(trim(res) .eq. '050')then
56 ndest = 4
57 allocate(destgrds(ndest))
58 destgrds = (/'9p00', '5p00', '1p00', '0p50'/)
59 end if
60 if(trim(res) .eq. '025')then
61 ndest = 5
62 allocate(destgrds(ndest))
63 destgrds = (/'9p00', '5p00', '1p00', '0p50', '0p25'/)
64 end if
65
66 !---------------------------------------------------------------------
67 ! use ESMF to create the weights from the Ct tripole to the rectilinear
68 ! grids with conservative and bilinear methods for post; the source
69 ! file is always Ct
70 !---------------------------------------------------------------------
71
72 do nd = 1,ndest
73 fsrc = trim(dirout)//'/'//'Ct.mx'//trim(res)//'_SCRIP.nc'
74 fdst = trim(dirout)//'/rect.'//trim(destgrds(nd))//'_SCRIP.nc'
75
76 do k = 1,size(methodname)
77 if(trim(methodname(k)) .eq. 'bilinear')method=esmf_regridmethod_bilinear
78 if(trim(methodname(k)) .eq. 'conserve')method=esmf_regridmethod_conserve
79
80 fwgt = trim(dirout)//'/'//'tripole.mx'//trim(res)//'.Ct.to.rect.'//trim(destgrds(nd)) &
81 //'.'//trim(methodname(k))//'.nc'
82 logmsg = 'creating weight file '//trim(fwgt)
83 print '(a)',trim(logmsg)
84
85 call esmf_regridweightgen(srcfile=trim(fsrc),dstfile=trim(fdst), &
86 weightfile=trim(fwgt), regridmethod=method, &
87 ignoredegenerate=.true., regridroutehandle=rh, &
88 unmappedaction=esmf_unmappedaction_ignore, rc=rc)
89 if (esmf_logfounderror(rctocheck=rc, msg=esmf_logerr_passthru, &
90 line=__line__, file=__file__)) call esmf_finalize(endflag=esmf_end_abort)
91 end do
92 end do
93
94 deallocate(destgrds)
95
96 end subroutine make_postwgts
97end module postwgts