cpld_gridgen  1.13.0
 All Data Structures Files Functions Variables Pages
tripolegrid.F90
Go to the documentation of this file.
1 
7 
8 module 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
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 
26 contains
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
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  !bathymetry
71  rc = nf90_def_var(ncid, 'depth', nf90_float, dim2, id)
72  rc = nf90_put_att(ncid, id, 'units', 'm')
73 
74  dim2(:) = (/idimid, jdimid/)
75  do ii = 1,ncoord
76  rc = nf90_def_var(ncid, trim(fixvars(ii)%var_name), nf90_double, dim2, id)
77  rc = nf90_put_att(ncid, id, 'units', trim(fixvars(ii)%unit_name))
78  rc = nf90_put_att(ncid, id, 'long_name', trim(fixvars(ii)%long_name))
79  if(trim(fixvars(ii)%var_name(1:3)) .eq. "lon")then
80  rc = nf90_put_att(ncid, id, 'lon_bounds', trim(fixvars(ii)%vertices))
81  else
82  rc = nf90_put_att(ncid, id, 'lat_bounds', trim(fixvars(ii)%vertices))
83  endif
84  enddo
85 
86  dim3(:) = (/idimid, jdimid, kdimid/)
87  do ii = ncoord+1,ncoord+nverts
88  rc = nf90_def_var(ncid, trim(fixvars(ii)%var_name), nf90_double, dim3, id)
89  rc = nf90_put_att(ncid, id, 'units', trim(fixvars(ii)%unit_name))
90  rc = nf90_put_att(ncid, id, 'long_name', trim(fixvars(ii)%long_name))
91  enddo
92 
93  rc = nf90_put_att(ncid, nf90_global, 'history', trim(history))
94  rc = nf90_enddef(ncid)
95 
96  rc = nf90_inq_varid(ncid, 'wet', id)
97  rc = nf90_put_var(ncid, id, int(wet4))
98 
99  rc = nf90_inq_varid(ncid, 'area', id)
100  rc = nf90_put_var(ncid, id, areact)
101 
102  rc = nf90_inq_varid(ncid,'anglet', id)
103  rc = nf90_put_var(ncid, id, anglet)
104 
105  rc = nf90_inq_varid(ncid, 'depth', id)
106  rc = nf90_put_var(ncid, id, dp4)
107 
108  rc = nf90_inq_varid(ncid, 'lonCt', id)
109  rc = nf90_put_var(ncid, id, lonct)
110 
111  rc = nf90_inq_varid(ncid, 'latCt', id)
112  rc = nf90_put_var(ncid, id, latct)
113 
114  rc = nf90_inq_varid(ncid, 'lonCv', id)
115  rc = nf90_put_var(ncid, id, loncv)
116 
117  rc = nf90_inq_varid(ncid, 'latCv', id)
118  rc = nf90_put_var(ncid, id, latcv)
119 
120  rc = nf90_inq_varid(ncid, 'lonCu', id)
121  rc = nf90_put_var(ncid, id, loncu)
122 
123  rc = nf90_inq_varid(ncid, 'latCu', id)
124  rc = nf90_put_var(ncid, id, latcu)
125 
126  rc = nf90_inq_varid(ncid, 'lonBu', id)
127  rc = nf90_put_var(ncid, id, lonbu)
128 
129  rc = nf90_inq_varid(ncid, 'latBu', id)
130  rc = nf90_put_var(ncid, id, latbu)
131 
132  ! vertices
133  rc = nf90_inq_varid(ncid, 'lonCt_vert', id)
134  rc = nf90_put_var(ncid, id, lonct_vert)
135 
136  rc = nf90_inq_varid(ncid, 'latCt_vert', id)
137  rc = nf90_put_var(ncid, id, latct_vert)
138 
139  rc = nf90_inq_varid(ncid, 'lonCv_vert', id)
140  rc = nf90_put_var(ncid, id, loncv_vert)
141 
142  rc = nf90_inq_varid(ncid, 'latCv_vert', id)
143  rc = nf90_put_var(ncid, id, latcv_vert)
144 
145  rc = nf90_inq_varid(ncid, 'lonCu_vert', id)
146  rc = nf90_put_var(ncid, id, loncu_vert)
147 
148  rc = nf90_inq_varid(ncid, 'latCu_vert', id)
149  rc = nf90_put_var(ncid, id, latcu_vert)
150 
151  rc = nf90_inq_varid(ncid, 'lonBu_vert', id)
152  rc = nf90_put_var(ncid, id, lonbu_vert)
153 
154  rc = nf90_inq_varid(ncid, 'latBu_vert', id)
155  rc = nf90_put_var(ncid, id, latbu_vert)
156 
157  rc = nf90_close(ncid)
158 
159  end subroutine write_tripolegrid
160 end module tripolegrid
subroutine, public write_tripolegrid(fname)
Write the tripole grid file.
Definition: tripolegrid.F90:33
subroutine fixvars_typedefine
Define the variables written to the tripole grid file.
Definition: vartypedefs.F90:34