cpld_gridgen  1.7.0
 All Data Structures Files Functions Variables
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,mastertask,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  if(mastertask) then
53  logmsg = '==> writing tripole grid to '//trim(fname)
54  print '(a)', trim(logmsg)
55  if(rc .ne. 0)print '(a)', 'nf90_create = '//trim(nf90_strerror(rc))
56  end if
57 
58  rc = nf90_def_dim(ncid, 'ni', ni, idimid)
59  rc = nf90_def_dim(ncid, 'nj', nj, jdimid)
60  rc = nf90_def_dim(ncid, 'nv', nv, kdimid)
61 
62  !mask
63  dim2(:) = (/idimid, jdimid/)
64  rc = nf90_def_var(ncid, 'wet', nf90_int, dim2, id)
65  rc = nf90_put_att(ncid, id, 'units', 'nd')
66  !area
67  rc = nf90_def_var(ncid, 'area', nf90_double, dim2, id)
68  rc = nf90_put_att(ncid, id, 'units', 'm2')
69  !angleT
70  rc = nf90_def_var(ncid, 'anglet', nf90_double, dim2, id)
71  rc = nf90_put_att(ncid, id, 'units', 'radians')
72  !bathymetry
73  rc = nf90_def_var(ncid, 'depth', nf90_float, dim2, id)
74  rc = nf90_put_att(ncid, id, 'units', 'm')
75 
76  dim2(:) = (/idimid, jdimid/)
77  do ii = 1,ncoord
78  rc = nf90_def_var(ncid, trim(fixvars(ii)%var_name), nf90_double, dim2, id)
79  rc = nf90_put_att(ncid, id, 'units', trim(fixvars(ii)%unit_name))
80  rc = nf90_put_att(ncid, id, 'long_name', trim(fixvars(ii)%long_name))
81  if(trim(fixvars(ii)%var_name(1:3)) .eq. "lon")then
82  rc = nf90_put_att(ncid, id, 'lon_bounds', trim(fixvars(ii)%vertices))
83  else
84  rc = nf90_put_att(ncid, id, 'lat_bounds', trim(fixvars(ii)%vertices))
85  endif
86  enddo
87 
88  dim3(:) = (/idimid, jdimid, kdimid/)
89  do ii = ncoord+1,ncoord+nverts
90  rc = nf90_def_var(ncid, trim(fixvars(ii)%var_name), nf90_double, dim3, id)
91  rc = nf90_put_att(ncid, id, 'units', trim(fixvars(ii)%unit_name))
92  rc = nf90_put_att(ncid, id, 'long_name', trim(fixvars(ii)%long_name))
93  enddo
94 
95  rc = nf90_put_att(ncid, nf90_global, 'history', trim(history))
96  rc = nf90_enddef(ncid)
97 
98  rc = nf90_inq_varid(ncid, 'wet', id)
99  rc = nf90_put_var(ncid, id, int(wet4))
100 
101  rc = nf90_inq_varid(ncid, 'area', id)
102  rc = nf90_put_var(ncid, id, areact)
103 
104  rc = nf90_inq_varid(ncid,'anglet', id)
105  rc = nf90_put_var(ncid, id, anglet)
106 
107  rc = nf90_inq_varid(ncid, 'depth', id)
108  rc = nf90_put_var(ncid, id, dp4)
109 
110  rc = nf90_inq_varid(ncid, 'lonCt', id)
111  rc = nf90_put_var(ncid, id, lonct)
112 
113  rc = nf90_inq_varid(ncid, 'latCt', id)
114  rc = nf90_put_var(ncid, id, latct)
115 
116  rc = nf90_inq_varid(ncid, 'lonCv', id)
117  rc = nf90_put_var(ncid, id, loncv)
118 
119  rc = nf90_inq_varid(ncid, 'latCv', id)
120  rc = nf90_put_var(ncid, id, latcv)
121 
122  rc = nf90_inq_varid(ncid, 'lonCu', id)
123  rc = nf90_put_var(ncid, id, loncu)
124 
125  rc = nf90_inq_varid(ncid, 'latCu', id)
126  rc = nf90_put_var(ncid, id, latcu)
127 
128  rc = nf90_inq_varid(ncid, 'lonBu', id)
129  rc = nf90_put_var(ncid, id, lonbu)
130 
131  rc = nf90_inq_varid(ncid, 'latBu', id)
132  rc = nf90_put_var(ncid, id, latbu)
133 
134  ! vertices
135  rc = nf90_inq_varid(ncid, 'lonCt_vert', id)
136  rc = nf90_put_var(ncid, id, lonct_vert)
137 
138  rc = nf90_inq_varid(ncid, 'latCt_vert', id)
139  rc = nf90_put_var(ncid, id, latct_vert)
140 
141  rc = nf90_inq_varid(ncid, 'lonCv_vert', id)
142  rc = nf90_put_var(ncid, id, loncv_vert)
143 
144  rc = nf90_inq_varid(ncid, 'latCv_vert', id)
145  rc = nf90_put_var(ncid, id, latcv_vert)
146 
147  rc = nf90_inq_varid(ncid, 'lonCu_vert', id)
148  rc = nf90_put_var(ncid, id, loncu_vert)
149 
150  rc = nf90_inq_varid(ncid, 'latCu_vert', id)
151  rc = nf90_put_var(ncid, id, latcu_vert)
152 
153  rc = nf90_inq_varid(ncid, 'lonBu_vert', id)
154  rc = nf90_put_var(ncid, id, lonbu_vert)
155 
156  rc = nf90_inq_varid(ncid, 'latBu_vert', id)
157  rc = nf90_put_var(ncid, id, latbu_vert)
158 
159  rc = nf90_close(ncid)
160 
161  end subroutine write_tripolegrid
162 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