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, staggerlocs, logmsg
14  use netcdf
15 
16  implicit none
17 
18 contains
23 
24  subroutine make_postwgts
25 
26  ! local variables
27  character(len=CL) :: fsrc, fdst, fwgt
28  character(len= 2) :: cstagger
29 
30  character(len=CM), dimension(2) :: methodname = (/'conserve', 'bilinear'/)
31 
32  type(esmf_routehandle) :: rh
33  type(esmf_regridmethod_flag) :: method
34  ! the number of possible destination grids depends on the source grid resolution
35  integer :: k,rc,nd,ndest
36  character(len=CS), allocatable, dimension(:) :: destgrds
37 
38  !---------------------------------------------------------------------
39  ! set the destination grids
40  !---------------------------------------------------------------------
41 
42  if(trim(res) .eq. '500')then
43  ndest = 1
44  allocate(destgrds(ndest))
45  destgrds = (/'5p00'/)
46  end if
47  if(trim(res) .eq. '100')then
48  ndest = 2
49  allocate(destgrds(ndest))
50  destgrds = (/'5p00', '1p00'/)
51  end if
52  if(trim(res) .eq. '050')then
53  ndest = 3
54  allocate(destgrds(ndest))
55  destgrds = (/'5p00', '1p00', '0p50'/)
56  end if
57  if(trim(res) .eq. '025')then
58  ndest = 4
59  allocate(destgrds(ndest))
60  destgrds = (/'5p00', '1p00', '0p50', '0p25'/)
61  end if
62 
63  !---------------------------------------------------------------------
64  ! use ESMF to create the weights for unstaggering the points onto
65  ! the Ct staggers for post; the destination is always Ct
66  !---------------------------------------------------------------------
67 
68  method=esmf_regridmethod_bilinear
69  fdst = trim(dirout)//'/'//'Ct.mx'//trim(res)//'_SCRIP.nc'
70  do k = 2,nv
71  cstagger = trim(staggerlocs(k))
72  fsrc = trim(dirout)//'/'//trim(cstagger)//'.mx'//trim(res)//'_SCRIP.nc'
73  fwgt = trim(dirout)//'/'//'tripole.mx'//trim(res)//'.'//trim(cstagger)//'.to.Ct.bilinear.nc'
74  logmsg = 'creating weight file '//trim(fwgt)
75  print '(a)',trim(logmsg)
76 
77  call esmf_regridweightgen(srcfile=trim(fsrc),dstfile=trim(fdst), &
78  weightfile=trim(fwgt), regridmethod=method, &
79  ignoredegenerate=.true., &
80  unmappedaction=esmf_unmappedaction_ignore, rc=rc)
81  if (esmf_logfounderror(rctocheck=rc, msg=esmf_logerr_passthru, &
82  line=__line__, file=__file__)) call esmf_finalize(endflag=esmf_end_abort)
83  end do
84 
85  !---------------------------------------------------------------------
86  ! use ESMF to create the weights from the Ct tripole to the rectilinear
87  ! grids with conservative and bilinear methods for post; the source
88  ! file is always Ct
89  !---------------------------------------------------------------------
90 
91  do nd = 1,ndest
92  fsrc = trim(dirout)//'/'//'Ct.mx'//trim(res)//'_SCRIP.nc'
93  fdst = trim(dirout)//'/rect.'//trim(destgrds(nd))//'_SCRIP.nc'
94 
95  do k = 1,size(methodname)
96  if(trim(methodname(k)) .eq. 'bilinear')method=esmf_regridmethod_bilinear
97  if(trim(methodname(k)) .eq. 'conserve')method=esmf_regridmethod_conserve
98 
99  fwgt = trim(dirout)//'/'//'tripole.mx'//trim(res)//'.Ct.to.rect.'//trim(destgrds(nd)) &
100  //'.'//trim(methodname(k))//'.nc'
101  logmsg = 'creating weight file '//trim(fwgt)
102  print '(a)',trim(logmsg)
103 
104  call esmf_regridweightgen(srcfile=trim(fsrc),dstfile=trim(fdst), &
105  weightfile=trim(fwgt), regridmethod=method, &
106  ignoredegenerate=.true., regridroutehandle=rh, &
107  unmappedaction=esmf_unmappedaction_ignore, rc=rc)
108  if (esmf_logfounderror(rctocheck=rc, msg=esmf_logerr_passthru, &
109  line=__line__, file=__file__)) call esmf_finalize(endflag=esmf_end_abort)
110  end do
111  end do
112 
113  deallocate(destgrds)
114 
115  end subroutine make_postwgts
116 end module postwgts
subroutine make_postwgts
Create the ESMF weights files to remap velocity points from their native stagger location to the cent...
Definition: postwgts.F90:24