cpld_gridgen 1.14.0
Loading...
Searching...
No Matches
tripolegrid.F90
Go to the documentation of this file.
1
7
8module tripolegrid
9
10 use gengrid_kinds, only: dbl_kind,int_kind,cm
11 use grdvars, only: ni,nj,nv,nverts,ncoord
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
16 use grdvars, only: wet4,areact,anglet,dp4,angle,angchk
17 use charstrings, only: logmsg,history
18 use vartypedefs, only: maxvars, fixvars, fixvars_typedefine
19 use netcdf
20
21 implicit none
22 private
23
24 public write_tripolegrid
25
26contains
32
33 subroutine write_tripolegrid(fname)
34
35 character(len=*), intent(in) :: fname
36
37 ! local variables
38 integer :: ii,id,rc, ncid, dim2(2),dim3(3)
39 integer :: idimid,jdimid,kdimid
40
41 !---------------------------------------------------------------------
42 ! create the netcdf file
43 !---------------------------------------------------------------------
44
45 ! define the output variables and file name
46 call fixvars_typedefine
47
48 ! create the file
49 ! 64_bit offset reqd for 008 grid
50 ! produces b4b results for smaller grids
51 rc = nf90_create(trim(fname), nf90_64bit_offset, ncid)
52 logmsg = '==> writing tripole grid to '//trim(fname)
53 print '(a)', trim(logmsg)
54 if(rc .ne. 0)print '(a)', 'nf90_create = '//trim(nf90_strerror(rc))
55
56 rc = nf90_def_dim(ncid, 'ni', ni, idimid)
57 rc = nf90_def_dim(ncid, 'nj', nj, jdimid)
58 rc = nf90_def_dim(ncid, 'nv', nv, kdimid)
59
60 !mask
61 dim2(:) = (/idimid, jdimid/)
62 rc = nf90_def_var(ncid, 'wet', nf90_int, dim2, id)
63 rc = nf90_put_att(ncid, id, 'units', 'nd')
64 !area
65 rc = nf90_def_var(ncid, 'area', nf90_double, dim2, id)
66 rc = nf90_put_att(ncid, id, 'units', 'm2')
67 !angleT
68 rc = nf90_def_var(ncid, 'anglet', nf90_double, dim2, id)
69 rc = nf90_put_att(ncid, id, 'units', 'radians')
70 !angle (angBu)
71 rc = nf90_def_var(ncid, 'angle', nf90_double, dim2, id)
72 rc = nf90_put_att(ncid, id, 'units', 'radians')
73 !angchk
74 rc = nf90_def_var(ncid, 'angchk', nf90_double, dim2, id)
75 rc = nf90_put_att(ncid, id, 'units', 'radians')
76 !bathymetry
77 rc = nf90_def_var(ncid, 'depth', nf90_float, dim2, id)
78 rc = nf90_put_att(ncid, id, 'units', 'm')
79
80 dim2(:) = (/idimid, jdimid/)
81 do ii = 1,ncoord
82 rc = nf90_def_var(ncid, trim(fixvars(ii)%var_name), nf90_double, dim2, id)
83 rc = nf90_put_att(ncid, id, 'units', trim(fixvars(ii)%unit_name))
84 rc = nf90_put_att(ncid, id, 'long_name', trim(fixvars(ii)%long_name))
85 if(trim(fixvars(ii)%var_name(1:3)) .eq. "lon")then
86 rc = nf90_put_att(ncid, id, 'lon_bounds', trim(fixvars(ii)%vertices))
87 else
88 rc = nf90_put_att(ncid, id, 'lat_bounds', trim(fixvars(ii)%vertices))
89 endif
90 enddo
91
92 dim3(:) = (/idimid, jdimid, kdimid/)
93 do ii = ncoord+1,ncoord+nverts
94 rc = nf90_def_var(ncid, trim(fixvars(ii)%var_name), nf90_double, dim3, id)
95 rc = nf90_put_att(ncid, id, 'units', trim(fixvars(ii)%unit_name))
96 rc = nf90_put_att(ncid, id, 'long_name', trim(fixvars(ii)%long_name))
97 enddo
98
99 rc = nf90_put_att(ncid, nf90_global, 'history', trim(history))
100 rc = nf90_enddef(ncid)
101
102 rc = nf90_inq_varid(ncid, 'wet', id)
103 rc = nf90_put_var(ncid, id, int(wet4))
104
105 rc = nf90_inq_varid(ncid, 'area', id)
106 rc = nf90_put_var(ncid, id, areact)
107
108 rc = nf90_inq_varid(ncid,'anglet', id)
109 rc = nf90_put_var(ncid, id, anglet)
110
111 rc = nf90_inq_varid(ncid, 'angle', id)
112 rc = nf90_put_var(ncid, id, angle)
113
114 rc = nf90_inq_varid(ncid,'angchk', id)
115 rc = nf90_put_var(ncid, id, angchk)
116
117 rc = nf90_inq_varid(ncid, 'depth', id)
118 rc = nf90_put_var(ncid, id, dp4)
119
120 rc = nf90_inq_varid(ncid, 'lonCt', id)
121 rc = nf90_put_var(ncid, id, lonct)
122
123 rc = nf90_inq_varid(ncid, 'latCt', id)
124 rc = nf90_put_var(ncid, id, latct)
125
126 rc = nf90_inq_varid(ncid, 'lonCv', id)
127 rc = nf90_put_var(ncid, id, loncv)
128
129 rc = nf90_inq_varid(ncid, 'latCv', id)
130 rc = nf90_put_var(ncid, id, latcv)
131
132 rc = nf90_inq_varid(ncid, 'lonCu', id)
133 rc = nf90_put_var(ncid, id, loncu)
134
135 rc = nf90_inq_varid(ncid, 'latCu', id)
136 rc = nf90_put_var(ncid, id, latcu)
137
138 rc = nf90_inq_varid(ncid, 'lonBu', id)
139 rc = nf90_put_var(ncid, id, lonbu)
140
141 rc = nf90_inq_varid(ncid, 'latBu', id)
142 rc = nf90_put_var(ncid, id, latbu)
143
144 ! vertices
145 rc = nf90_inq_varid(ncid, 'lonCt_vert', id)
146 rc = nf90_put_var(ncid, id, lonct_vert)
147
148 rc = nf90_inq_varid(ncid, 'latCt_vert', id)
149 rc = nf90_put_var(ncid, id, latct_vert)
150
151 rc = nf90_inq_varid(ncid, 'lonCv_vert', id)
152 rc = nf90_put_var(ncid, id, loncv_vert)
153
154 rc = nf90_inq_varid(ncid, 'latCv_vert', id)
155 rc = nf90_put_var(ncid, id, latcv_vert)
156
157 rc = nf90_inq_varid(ncid, 'lonCu_vert', id)
158 rc = nf90_put_var(ncid, id, loncu_vert)
159
160 rc = nf90_inq_varid(ncid, 'latCu_vert', id)
161 rc = nf90_put_var(ncid, id, latcu_vert)
162
163 rc = nf90_inq_varid(ncid, 'lonBu_vert', id)
164 rc = nf90_put_var(ncid, id, lonbu_vert)
165
166 rc = nf90_inq_varid(ncid, 'latBu_vert', id)
167 rc = nf90_put_var(ncid, id, latbu_vert)
168
169 rc = nf90_close(ncid)
170
171 end subroutine write_tripolegrid
172end module tripolegrid