11 use grdvars, only: ni,nj,nv,mastertask
12 use grdvars, only: lonct,latct,lonct_vert,latct_vert
13 use grdvars, only: loncu,latcu,loncu_vert,latcu_vert
14 use grdvars, only: loncv,latcv,loncv_vert,latcv_vert
15 use grdvars, only: lonbu,latbu,lonbu_vert,latbu_vert
36 character(len=*) ,
intent(in) :: fname
37 character(len=*) ,
intent(in) :: cstagger
38 integer(int_kind),
optional,
intent(in) :: imask(:,:)
41 integer,
parameter :: grid_rank = 2
43 integer :: ii,n,id,rc, ncid, dim2(2),dim1(1)
44 integer :: idimid,jdimid,kdimid
46 integer,
dimension(grid_rank) :: gdims
47 integer(int_kind),
dimension(ni*nj) :: cnmask
48 real(dbl_kind),
dimension(ni*nj) :: cnlons, cnlats
49 real(dbl_kind),
dimension(nv,ni*nj) :: crlons, crlats
51 real(dbl_kind),
dimension(ni,nj) :: tmp
53 character(len=2) :: vtype
54 character(len=CM) :: vname
55 character(len=CM) :: vunit
62 if(trim(cstagger) .eq.
'Ct')
then
63 cnlons = reshape(lonct, (/ni*nj/))
64 cnlats = reshape(latct, (/ni*nj/))
66 tmp(:,:) = lonct_vert(:,:,n)
67 crlons(n,:) = reshape(tmp, (/ni*nj/))
68 tmp(:,:) = latct_vert(:,:,n)
69 crlats(n,:) = reshape(tmp, (/ni*nj/))
73 if(trim(cstagger) .eq.
'Cu')
then
74 cnlons = reshape(loncu, (/ni*nj/))
75 cnlats = reshape(latcu, (/ni*nj/))
77 tmp(:,:) = loncu_vert(:,:,n)
78 crlons(n,:) = reshape(tmp, (/ni*nj/))
79 tmp(:,:) = latcu_vert(:,:,n)
80 crlats(n,:) = reshape(tmp, (/ni*nj/))
84 if(trim(cstagger) .eq.
'Cv')
then
85 cnlons = reshape(loncv, (/ni*nj/))
86 cnlats = reshape(latcv, (/ni*nj/))
88 tmp(:,:) = loncv_vert(:,:,n)
89 crlons(n,:) = reshape(tmp, (/ni*nj/))
90 tmp(:,:) = latcv_vert(:,:,n)
91 crlats(n,:) = reshape(tmp, (/ni*nj/))
95 if(trim(cstagger) .eq.
'Bu')
then
96 cnlons = reshape(lonbu, (/ni*nj/))
97 cnlats = reshape(latbu, (/ni*nj/))
99 tmp(:,:) = lonbu_vert(:,:,n)
100 crlons(n,:) = reshape(tmp, (/ni*nj/))
101 tmp(:,:) = latbu_vert(:,:,n)
102 crlats(n,:) = reshape(tmp, (/ni*nj/))
106 if(present(imask))
then
107 cnmask = reshape(imask, (/ni*nj/))
121 rc = nf90_create(trim(fname), nf90_64bit_offset, ncid)
123 logmsg =
'==> writing SCRIP grid to '//trim(fname)
124 print
'(a)',trim(logmsg)
125 if(rc .ne. 0)print
'(a)',
'nf90_create = '//trim(nf90_strerror(rc))
128 rc = nf90_def_dim(ncid,
'grid_size', ni*nj, idimid)
129 rc = nf90_def_dim(ncid,
'grid_corners', nv, jdimid)
130 rc = nf90_def_dim(ncid,
'grid_rank', grid_rank, kdimid)
134 rc = nf90_def_var(ncid,
'grid_dims', nf90_int, dim1, id)
137 rc = nf90_def_var(ncid,
'grid_imask', nf90_int, dim1, id)
138 rc = nf90_put_att(ncid, id,
'units',
'unitless')
142 vname = trim(scripvars(ii)%var_name)
143 vunit = trim(scripvars(ii)%unit_name)
144 vtype = trim(scripvars(ii)%var_type)
146 if(vtype .eq.
'r8')rc = nf90_def_var(ncid, vname, nf90_double, dim1, id)
147 if(vtype .eq.
'r4')rc = nf90_def_var(ncid, vname, nf90_float, dim1, id)
148 if(vtype .eq.
'i4')rc = nf90_def_var(ncid, vname, nf90_int, dim1, id)
149 rc = nf90_put_att(ncid, id,
'units', vunit)
154 vname = trim(scripvars(ii)%var_name)
155 vunit = trim(scripvars(ii)%unit_name)
156 vtype = trim(scripvars(ii)%var_type)
157 dim2(:) = (/jdimid,idimid/)
158 if(vtype .eq.
'r8')rc = nf90_def_var(ncid, vname, nf90_double, dim2, id)
159 if(vtype .eq.
'r4')rc = nf90_def_var(ncid, vname, nf90_float, dim2, id)
160 if(vtype .eq.
'i4')rc = nf90_def_var(ncid, vname, nf90_int, dim2, id)
161 rc = nf90_put_att(ncid, id,
'units', vunit)
163 rc = nf90_enddef(ncid)
165 rc = nf90_inq_varid(ncid,
'grid_dims', id)
166 rc = nf90_put_var(ncid, id, gdims)
167 rc = nf90_inq_varid(ncid,
'grid_imask', id)
168 rc = nf90_put_var(ncid, id, cnmask)
170 rc = nf90_inq_varid(ncid,
'grid_center_lon', id)
171 rc = nf90_put_var(ncid, id, cnlons)
172 rc = nf90_inq_varid(ncid,
'grid_center_lat', id)
173 rc = nf90_put_var(ncid, id, cnlats)
175 rc = nf90_inq_varid(ncid,
'grid_corner_lon', id)
176 rc = nf90_put_var(ncid, id, crlons)
177 rc = nf90_inq_varid(ncid,
'grid_corner_lat', id)
178 rc = nf90_put_var(ncid, id, crlats)
180 rc = nf90_close(ncid)
subroutine, public write_scripgrid(fname, cstagger, imask)
Write a SCRIP grid file.
subroutine scripvars_typedefine
Define the variables written to any SCRIP grid file.