-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathio_code.f90
295 lines (240 loc) · 11 KB
/
io_code.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
!>@author
!>Paul Connolly, The University of Manchester
!>@brief
!>io routines
module io_module
use numerics_type
private
public :: output_2d
contains
!>@author
!>Paul J. Connolly, The University of Manchester
!>@brief
!>output 1 time-step of model
!>@param[in] time in seconds
!>@param[in] nq number of q fields
!>@param[in] nprec number of precipitation fields
!>@param[in] ip number of horizontal levels
!>@param[in] kp number of vertical levels
!>@param[in] q_name: name of q-variables
!>@param[in] q, precip, theta, pressure, x,xn,z,zn, temperature,u,w
!>@param[inout] new_file
subroutine output_2d(time,nq,nprec,ip,kp,q_name, &
q,precip,theta,p,x,xn,z,zn,t,u,w,new_file)
use numerics_type
use netcdf
use variables, only : io1, outputfile
implicit none
real(wp), intent(in) :: time
integer(i4b), intent(in) :: nq,nprec,kp,ip
character(len=20), dimension(nq) :: q_name
real(wp), dimension(kp,ip,nq), intent(in) :: q
real(wp), dimension(kp,ip,nprec), intent(in) :: precip
real(wp), dimension(kp,ip), intent(in) :: theta, p, t, u, w
real(wp), dimension(ip), intent(in) :: x,xn
real(wp), dimension(kp), intent(in) :: z,zn
logical, intent(inout) :: new_file
! output to netcdf file
if(new_file) then
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! open / create the netcdf file !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
call check( nf90_create(outputfile, NF90_CLOBBER, io1%ncid) )
! define dimensions (netcdf hands back a handle)
call check( nf90_def_dim(io1%ncid, "times", NF90_UNLIMITED, io1%x_dimid) )
call check( nf90_def_dim(io1%ncid, "nq", nq, io1%nq_dimid) )
call check( nf90_def_dim(io1%ncid, "nprec", nprec, io1%nprec_dimid) )
call check( nf90_def_dim(io1%ncid, "ip", ip, io1%i_dimid) )
call check( nf90_def_dim(io1%ncid, "kp", kp, io1%k_dimid) )
call check( nf90_def_dim(io1%ncid, "l_q_names", 20, io1%lq_dimid) )
! close the file, freeing up any internal netCDF resources
! associated with the file, and flush any buffers
call check( nf90_close(io1%ncid) )
! now define some variables, units, etc
call check( nf90_open(outputfile, NF90_WRITE, io1%ncid) )
! define mode
call check( nf90_redef(io1%ncid) )
! define variable: q_names
call check( nf90_def_var(io1%ncid, "q_names", NF90_CHAR, &
(/io1%lq_dimid,io1%nq_dimid/), io1%varid) )
! define variable: time
call check( nf90_def_var(io1%ncid, "time", NF90_DOUBLE, &
(/io1%x_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "time", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "seconds") )
! define variable: q
call check( nf90_def_var(io1%ncid, "q", NF90_DOUBLE, &
(/io1%nq_dimid, io1%k_dimid, io1%i_dimid, io1%x_dimid/), io1%varid) )
! call check( nf90_def_var(io1%ncid, "q", NF90_DOUBLE, &
! (/io1%k_dimid, io1%i_dimid,io1%nq_dimid, io1%x_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "q", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "kg or number per kg") )
! define variable: precip
call check( nf90_def_var(io1%ncid, "precip", NF90_DOUBLE, &
(/io1%nprec_dimid, io1%k_dimid, io1%i_dimid, io1%x_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "precip", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "mm hr-1") )
! define variable: theta
call check( nf90_def_var(io1%ncid, "theta", NF90_DOUBLE, &
(/io1%k_dimid, io1%i_dimid, io1%x_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "theta", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "K") )
! define variable: p
call check( nf90_def_var(io1%ncid, "p", NF90_DOUBLE, &
(/io1%k_dimid, io1%i_dimid, io1%x_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "p", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "Pa") )
! define variable: x
call check( nf90_def_var(io1%ncid, "x", NF90_DOUBLE, &
(/io1%i_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "x", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "m") )
! define variable: xn
call check( nf90_def_var(io1%ncid, "xn", NF90_DOUBLE, &
(/io1%i_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "xn", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "m") )
! define variable: z
call check( nf90_def_var(io1%ncid, "z", NF90_DOUBLE, &
(/io1%k_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "z", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "m") )
! define variable: zn
call check( nf90_def_var(io1%ncid, "zn", NF90_DOUBLE, &
(/io1%k_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "zn", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "m") )
! define variable: t
call check( nf90_def_var(io1%ncid, "t", NF90_DOUBLE, &
(/io1%k_dimid, io1%i_dimid, io1%x_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "t", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "K") )
! define variable: u
call check( nf90_def_var(io1%ncid, "u", NF90_DOUBLE, &
(/io1%k_dimid, io1%i_dimid, io1%x_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "u", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "ms-1") )
! define variable: w
call check( nf90_def_var(io1%ncid, "w", NF90_DOUBLE, &
(/io1%k_dimid, io1%i_dimid, io1%x_dimid/), io1%varid) )
! get id to a_dimid
call check( nf90_inq_varid(io1%ncid, "w", io1%a_dimid) )
! units
call check( nf90_put_att(io1%ncid, io1%a_dimid, &
"units", "ms-1") )
call check( nf90_enddef(io1%ncid) )
call check( nf90_close(io1%ncid) )
new_file=.false.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! write x,xn,z,zn data to file !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
call check( nf90_open(outputfile, NF90_WRITE, io1%ncid) )
call check( nf90_inq_varid(io1%ncid, "q_names", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, q_name, start = (/1,1/)))
call check( nf90_inq_varid(io1%ncid, "x", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, x, start = (/1/)))
call check( nf90_inq_varid(io1%ncid, "xn", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, xn, start = (/1/)))
call check( nf90_inq_varid(io1%ncid, "z", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, z, start = (/1/)))
call check( nf90_inq_varid(io1%ncid, "zn", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, zn, start = (/1/)))
call check( nf90_close(io1%ncid) )
endif
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! write data to file !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
call check( nf90_open(outputfile, NF90_WRITE, io1%ncid) )
! write variable: time
call check( nf90_inq_varid(io1%ncid, "time", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, time, &
start = (/io1%icur/)))
!
! write variable: q
call check( nf90_inq_varid(io1%ncid, "q", io1%varid ) )
! had to reshape twice
call check( nf90_put_var(io1%ncid, io1%varid, &
reshape(reshape(q,[nq,ip,kp],order=[3,2,1]),[nq,kp,ip],order=[1,3,2] ), &
start = (/1,1,1,io1%icur/)))
! call check( nf90_put_var(io1%ncid, io1%varid, q, &
! start = (/1,1,1,io1%icur/)))
! write variable: precip
call check( nf90_inq_varid(io1%ncid, "precip", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, &
reshape(reshape(precip,[nprec,ip,kp],order=[3,2,1]),[nprec,kp,ip],order=[1,3,2] ), &
start = (/1,1,1,io1%icur/)))
! call check( nf90_put_var(io1%ncid, io1%varid, precip, &
! start = (/1,1,1,io1%icur/)))
! write variable: theta
call check( nf90_inq_varid(io1%ncid, "theta", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, theta, &
start = (/1,1,io1%icur/)))
! write variable: p
call check( nf90_inq_varid(io1%ncid, "p", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, p, &
start = (/1,1,io1%icur/)))
! write variable: u
call check( nf90_inq_varid(io1%ncid, "u", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, u, &
start = (/1,1,io1%icur/)))
! write variable: w
call check( nf90_inq_varid(io1%ncid, "w", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, w, &
start = (/1,1,io1%icur/)))
! write variable: t
call check( nf90_inq_varid(io1%ncid, "t", io1%varid ) )
call check( nf90_put_var(io1%ncid, io1%varid, t, &
start = (/1,1,io1%icur/)))
call check( nf90_close(io1%ncid) )
!
!
io1%icur=io1%icur+1
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
end subroutine output_2d
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! HELPER ROUTINE !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine check(status)
use netcdf
use numerics_type
integer(I4B), intent ( in) :: status
if(status /= nf90_noerr) then
print *, trim(nf90_strerror(status))
stop "Stopped"
end if
end subroutine check
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
end module io_module