cpld_gridgen  1.7.0
 All Data Structures Files Functions Variables
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, mastertask
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_regridmethod_flag) :: method
33  ! the number of possible destination grids depends on the source grid resolution
34  integer :: k,rc,nd,ndest
35  character(len=CS), allocatable, dimension(:) :: destgrds
36 
37 !---------------------------------------------------------------------
38 ! set the destination grids
39 !---------------------------------------------------------------------
40 
41  if(trim(res) .eq. '400')return
42 
43  if(trim(res) .eq. '100')then
44  ndest = 1
45  allocate(destgrds(ndest))
46  destgrds = (/'1p0 '/)
47  end if
48  if(trim(res) .eq. '050')then
49  ndest = 2
50  allocate(destgrds(ndest))
51  destgrds = (/'1p0 ', '0p5 '/)
52  end if
53  if(trim(res) .eq. '025')then
54  ndest = 3
55  allocate(destgrds(ndest))
56  destgrds = (/'1p0 ', '0p5 ', '0p25'/)
57  end if
58 
59 !---------------------------------------------------------------------
60 ! use ESMF to create the weights for unstaggering the points onto
61 ! the Ct staggers for post; the destination is always Ct
62 !---------------------------------------------------------------------
63 
64  method=esmf_regridmethod_bilinear
65  fdst = trim(dirout)//'/'//'Ct.mx'//trim(res)//'_SCRIP.nc'
66  do k = 2,nv
67  cstagger = trim(staggerlocs(k))
68  fsrc = trim(dirout)//'/'//trim(cstagger)//'.mx'//trim(res)//'_SCRIP.nc'
69  fwgt = trim(dirout)//'/'//'tripole.mx'//trim(res)//'.'//trim(cstagger)//'.to.Ct.bilinear.nc'
70  if(mastertask) then
71  logmsg = 'creating weight file '//trim(fwgt)
72  print '(a)',trim(logmsg)
73  end if
74 
75  call esmf_regridweightgen(srcfile=trim(fsrc),dstfile=trim(fdst), &
76  weightfile=trim(fwgt), regridmethod=method, &
77  ignoredegenerate=.true., &
78  unmappedaction=esmf_unmappedaction_ignore, rc=rc)
79  if (esmf_logfounderror(rctocheck=rc, msg=esmf_logerr_passthru, &
80  line=__line__, file=__file__)) call esmf_finalize(endflag=esmf_end_abort)
81  end do
82 
83 !---------------------------------------------------------------------
84 ! use ESMF to create the weights from the Ct tripole to the rectilinear
85 ! grids with conservative and bilinear methods for post; the source
86 ! file is always Ct
87 !---------------------------------------------------------------------
88 
89  do nd = 1,ndest
90  fsrc = trim(dirout)//'/'//'Ct.mx'//trim(res)//'_SCRIP.nc'
91  fdst = trim(dirout)//'/rect.'//trim(destgrds(nd))//'_SCRIP.nc'
92 
93  do k = 1,size(methodname)
94  if(trim(methodname(k)) .eq. 'bilinear')method=esmf_regridmethod_bilinear
95  if(trim(methodname(k)) .eq. 'conserve')method=esmf_regridmethod_conserve
96 
97  fwgt = trim(dirout)//'/'//'tripole.mx'//trim(res)//'.Ct.to.rect.'//trim(destgrds(nd)) &
98  //'.'//trim(methodname(k))//'.nc'
99  if(mastertask) then
100  logmsg = 'creating weight file '//trim(fwgt)
101  print '(a)',trim(logmsg)
102  end if
103 
104  call esmf_regridweightgen(srcfile=trim(fsrc),dstfile=trim(fdst), &
105  weightfile=trim(fwgt), regridmethod=method, &
106  ignoredegenerate=.true., &
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