11 use init_mod ,
only : nlevs, nxr, nyr
12 use init_mod ,
only : debug, logunit
13 use arrays_mod ,
only : b2d, c2d, b3d, rgb2d, rgc2d, rgb3d
14 use arrays_mod ,
only : nbilin2d, nconsd2d, nbilin3d, bilin2d, consd2d, bilin3d
15 use utils_mod ,
only : nf90_err
29 subroutine setup_icerestart(fin, fout)
31 character(len=*),
intent(in) :: fin, fout
33 integer :: istep1, myear, mmonth, mday, msec
34 integer :: ncid,varid,n
35 integer :: idimid,jdimid,kdimid
36 character(len=20) :: vname
37 character(len=40) :: subname =
'setup_icerestart' 39 if (debug)
write(logunit,
'(a)')
'enter '//trim(subname)
41 call nf90_err(nf90_open(trim(fin), nf90_nowrite, ncid),
'open: '//trim(fin))
43 call nf90_err(nf90_get_att(ncid, nf90_global,
'istep1', istep1), &
44 'get global attribute istep1 '//trim(fin))
45 call nf90_err(nf90_get_att(ncid, nf90_global,
'myear', myear), &
46 'get global attribute myear '//trim(fin))
47 call nf90_err(nf90_get_att(ncid, nf90_global,
'mmonth', mmonth), &
48 'get global attribute mmonth '//trim(fin))
49 call nf90_err(nf90_get_att(ncid, nf90_global,
'mday', mday), &
50 'get global attribute mday '//trim(fin))
51 call nf90_err(nf90_get_att(ncid, nf90_global,
'msec', msec), &
52 'get global attribute msec '//trim(fin))
53 call nf90_err(nf90_close(ncid),
'close: '//trim(fin))
56 call nf90_err(nf90_create(trim(fout), nf90_clobber, ncid),
'create: '//trim(fout))
57 call nf90_err(nf90_def_dim(ncid,
'ni', nxr, idimid),
'define dimension: ni')
58 call nf90_err(nf90_def_dim(ncid,
'nj', nyr, jdimid),
'define dimension: nj')
59 call nf90_err(nf90_def_dim(ncid,
'ncat', nlevs, kdimid),
'define dimension: ncat')
61 if (
allocated(b2d))
then 63 vname = trim(b2d(n)%var_name)
64 call nf90_err(nf90_def_var(ncid, vname, nf90_double, (/idimid,jdimid/), varid), &
65 'define variable: '// vname)
68 if (
allocated(c2d))
then 70 vname = trim(c2d(n)%var_name)
71 call nf90_err(nf90_def_var(ncid, vname, nf90_double, (/idimid,jdimid/), varid), &
72 'define variable: '// vname)
75 if (
allocated(b3d))
then 77 vname = trim(b3d(n)%var_name)
78 call nf90_err(nf90_def_var(ncid, vname, nf90_double, (/idimid,jdimid,kdimid/), varid), &
79 'define variable: '// vname)
83 call nf90_err(nf90_put_att(ncid, nf90_global,
'istep1', istep1), &
84 'put global attribute istep1')
85 call nf90_err(nf90_put_att(ncid, nf90_global,
'myear', myear), &
86 'put global attribute myear')
87 call nf90_err(nf90_put_att(ncid, nf90_global,
'mmonth', mmonth), &
88 'put global attribute mmonth')
89 call nf90_err(nf90_put_att(ncid, nf90_global,
'mday', mday), &
90 'put global attribute mday')
91 call nf90_err(nf90_put_att(ncid, nf90_global,
'msec', msec), &
92 'put global attribute msec')
93 call nf90_err(nf90_enddef(ncid),
'enddef: '// trim(fout))
94 call nf90_err(nf90_close(ncid),
'close: '// trim(fout))
96 if (debug)
write(logunit,
'(a)')
'exit '//trim(subname)
97 end subroutine setup_icerestart
106 subroutine setup_ocnrestart(fin, fout, bathy)
108 character(len=*),
intent(in) :: fin, fout
109 real(kind=8),
intent(in) :: bathy(:)
111 real(kind=8) :: timestamp
112 character(len= 40) :: timeunit
113 character(len= 20) :: vname, vunit
114 character(len=120) :: vlong
115 real(kind=8),
allocatable :: Layer(:)
116 real(kind=8),
allocatable :: out3d(:,:,:)
118 integer :: k,n,ncid,varid
119 integer :: idimid,jdimid,kdimid,edimid,timid
121 allocate(out3d(nxr,nyr,nlevs+1)); out3d = 0.0
123 out3d(:,:,k) = reshape(-bathy(:), (/nxr,nyr/))
126 call nf90_err(nf90_open(trim(fin), nf90_nowrite, ncid),
'open: '//trim(fin))
128 call nf90_err(nf90_inq_varid(ncid,
'Time', varid), &
129 'get variable Id: Time '//trim(fin))
130 call nf90_err(nf90_get_var(ncid, varid, timestamp), &
131 'get variable: timestamp '//trim(fin))
132 call nf90_err(nf90_get_att(ncid, varid,
'units', timeunit), &
133 'get variable attribute : units '//trim(fin))
135 allocate(layer(nlevs)) ; layer = 0.0
136 call nf90_err(nf90_inq_varid(ncid,
'Layer', varid), &
137 'get variable Id: Layer '//trim(fin))
138 call nf90_err(nf90_get_var(ncid, varid, layer), &
139 'get variable: Layer '//trim(fin))
140 call nf90_err(nf90_close(ncid),
'close: '//trim(fin))
142 call nf90_err(nf90_create(trim(fout), nf90_clobber, ncid),
'create: '//trim(fout))
143 call nf90_err(nf90_def_dim(ncid,
'nx', nxr, idimid), &
144 'define dimension: nx')
145 call nf90_err(nf90_def_dim(ncid,
'ny', nyr, jdimid), &
146 'define dimension: ny')
147 call nf90_err(nf90_def_dim(ncid,
'Layer', nlevs, kdimid), &
148 'define dimension: Layer')
149 call nf90_err(nf90_def_dim(ncid,
'Interface', nlevs+1, edimid), &
150 'define dimension: Interface')
151 call nf90_err(nf90_def_dim(ncid,
'Time', nf90_unlimited, timid), &
152 'define dimension: Time')
154 call nf90_err(nf90_def_var(ncid,
'Time', nf90_double, (/timid/), varid), &
155 'define variable: Time')
156 call nf90_err(nf90_put_att(ncid, varid,
'units', trim(timeunit)), &
157 'put variable attribute: units')
159 call nf90_err(nf90_def_var(ncid,
'Layer', nf90_double, (/kdimid/), varid), &
160 'define variable: Layer')
161 call nf90_err(nf90_put_att(ncid, varid,
'units',
'm'), &
162 'put variable attribute: units')
164 if (
allocated(b2d))
then 166 vname = trim(b2d(n)%var_name)
167 vunit = trim(b2d(n)%units)
168 vlong = trim(b2d(n)%long_name)
169 call nf90_err(nf90_def_var(ncid, vname, nf90_double, (/idimid,jdimid,timid/), varid), &
170 'define variable: '// vname)
171 call nf90_err(nf90_put_att(ncid, varid,
'units', vunit), &
172 'put variable attribute: units')
173 call nf90_err(nf90_put_att(ncid, varid,
'long_name', vlong), &
174 'put variable attribute: long_name')
177 if (
allocated(c2d))
then 179 vname = trim(c2d(n)%var_name)
180 vunit = trim(c2d(n)%units)
181 vlong = trim(c2d(n)%long_name)
182 call nf90_err(nf90_def_var(ncid, vname, nf90_double, (/idimid,jdimid,timid/), varid), &
183 'define variable: '// vname)
184 call nf90_err(nf90_put_att(ncid, varid,
'units', vunit), &
185 'put variable attribute: units' )
186 call nf90_err(nf90_put_att(ncid, varid,
'long_name', vlong), &
187 'put variable attribute: long_name' )
191 if (
allocated(b3d))
then 193 vname = trim(b3d(n)%var_name)
194 vunit = trim(b3d(n)%units)
195 vlong = trim(b3d(n)%long_name)
196 if (vname .eq.
'eta')
then 197 call nf90_err(nf90_def_var(ncid, vname, nf90_double, (/idimid,jdimid,edimid,timid/), varid), &
198 'define variable: '// vname)
200 call nf90_err(nf90_def_var(ncid, vname, nf90_double, (/idimid,jdimid,kdimid,timid/), varid), &
201 'define variable: '// vname)
203 call nf90_err(nf90_put_att(ncid, varid,
'units', vunit), &
204 'put variable attribute: units' )
205 call nf90_err(nf90_put_att(ncid, varid,
'long_name', vlong), &
206 'put variable attribute: long_name' )
209 call nf90_err(nf90_enddef(ncid),
'enddef: '// trim(fout))
212 call nf90_err(nf90_inq_varid(ncid,
'Time', varid), &
213 'get variable Id: Time')
214 call nf90_err(nf90_put_var(ncid, varid, timestamp), &
215 'put variable: Time')
217 call nf90_err(nf90_inq_varid(ncid,
'Layer', varid), &
218 'get variable Id: Layer')
219 call nf90_err(nf90_put_var(ncid, varid, layer) , &
220 'put variable: Layer')
222 call nf90_err(nf90_inq_varid(ncid,
'eta', varid), &
223 'get variable Id: eta')
224 call nf90_err(nf90_put_var(ncid, varid, out3d) , &
226 call nf90_err(nf90_close(ncid),
'close: '//trim(fout))
228 end subroutine setup_ocnrestart
229 end module restarts_mod