1 | MODULE fpmetbinary_mod |
---|
2 | |
---|
3 | !***************************************************************************** |
---|
4 | ! * |
---|
5 | ! Contains data and routines for dumping and loading processed met * |
---|
6 | ! fields. * |
---|
7 | ! Authors Don Morton (Don.Morton@borealscicomp.com) * |
---|
8 | ! Delia Arnold (deliona.arnold@gmail.com) * |
---|
9 | ! * |
---|
10 | ! 07 Oct 2016 * |
---|
11 | ! * |
---|
12 | ! Most of the data structures from com_mod.f90 that are dumped and * |
---|
13 | ! loaded have a final dimension of size two, so that they may hold data * |
---|
14 | ! from two met files. When we dump the contents into a .fp file, we * |
---|
15 | ! need to specify which of the two to dump. Likewise, when we load * |
---|
16 | ! from a .fp file, we need to specify which of the two possible indices * |
---|
17 | ! to load into. * |
---|
18 | ! * |
---|
19 | ! Note that these routines need more robustness. For example, what * |
---|
20 | ! what happens if the filename can't be read or written. Or, what * |
---|
21 | ! happens if a read or write fails in any way. Right now, it's crash * |
---|
22 | ! city. * |
---|
23 | ! * |
---|
24 | ! Recent enhancements (07 Oct 2016) DJM: * |
---|
25 | ! * |
---|
26 | ! - file format changed so that compiled dimensions are output, and * |
---|
27 | ! during input these same dimensions are compared with the dimensions * |
---|
28 | ! compiled into the flexpart that is reading it. A discrepancy * |
---|
29 | ! causes abort, so that time isn't wasted reading an incompatible * |
---|
30 | ! file. * |
---|
31 | ! * |
---|
32 | ! - file format changed so that first item is an 8-character string * |
---|
33 | ! depicting the version of the preprocessed file format. * |
---|
34 | ! An inconsistency between a detected and expected string results * |
---|
35 | ! in program abort. * |
---|
36 | ! * |
---|
37 | ! *** IMPORTANT *** - when the format of the preprocessed output is * |
---|
38 | ! modified in any way, be sure to change the version string below, * |
---|
39 | ! PREPROC_FORMAT_VERSION_STR, so that attempts to read the output * |
---|
40 | ! with a different format version will cause an abort. * |
---|
41 | ! * |
---|
42 | !***************************************************************************** |
---|
43 | |
---|
44 | USE com_mod |
---|
45 | USE conv_mod |
---|
46 | USE par_mod, ONLY : nxmax, nymax, nzmax, nuvzmax, nwzmax, numclass |
---|
47 | |
---|
48 | USE netcdf |
---|
49 | |
---|
50 | IMPLICIT NONE |
---|
51 | |
---|
52 | ! Users may want to change these IO Unit values if they conflict with other parts |
---|
53 | ! of code |
---|
54 | INTEGER, PARAMETER :: IOUNIT_DUMP = 33, IOUNIT_LOAD = 34, & |
---|
55 | IOUNIT_TEXTOUT = 35 |
---|
56 | |
---|
57 | ! When a change is made to the format of the preprocessed file, such that |
---|
58 | ! this routine will not be able to read a previous version, this version |
---|
59 | ! string should be modified |
---|
60 | CHARACTER(LEN=12), PARAMETER :: PREPROC_FORMAT_VERSION_STR = 'FP_p-9.3.1'//char(0) |
---|
61 | |
---|
62 | PRIVATE IOUNIT_DUMP, IOUNIT_LOAD, IOUNIT_TEXTOUT, fpio, & |
---|
63 | & PREPROC_FORMAT_VERSION_STR |
---|
64 | |
---|
65 | |
---|
66 | CONTAINS |
---|
67 | |
---|
68 | !***************************************************************************** |
---|
69 | ! * |
---|
70 | ! Subroutines fpmetbinary_dump() and fpmetbinary_load() provide the * |
---|
71 | ! public interface to * |
---|
72 | ! this module functionality. I created the PRIVATE fpio() because I * |
---|
73 | ! wanted all interactions with variables to be in one place. The read * |
---|
74 | ! and write operations need to be done in exactly the same sequence, so * |
---|
75 | ! I felt like keeping them in the same routine would at least allow for * |
---|
76 | ! coders to more easily compare the two sequences than if they were * |
---|
77 | ! separate. * |
---|
78 | ! * |
---|
79 | ! As mentioned above, the dumps and loads will, for most variables, * |
---|
80 | ! need to refer to one of two index values for the last dimension of * |
---|
81 | ! the array. * |
---|
82 | ! * |
---|
83 | !***************************************************************************** |
---|
84 | |
---|
85 | |
---|
86 | SUBROUTINE fpmetbinary_dump(filename, cm_index) |
---|
87 | CHARACTER(LEN=*), INTENT(IN) :: filename ! Full path for file |
---|
88 | INTEGER, INTENT(IN) :: cm_index ! Index of last dimension in |
---|
89 | ! most com_mod variables. |
---|
90 | ! Should be 1 or 2 |
---|
91 | |
---|
92 | INTEGER millisecs_start, millisecs_stop, count_rate, count_max |
---|
93 | |
---|
94 | INTEGER :: ncretval, ncid ! NetCDF func return value, file id |
---|
95 | |
---|
96 | CALL SYSTEM_CLOCK(millisecs_start, count_rate, count_max) |
---|
97 | |
---|
98 | ! Create and open NC4 file for writing |
---|
99 | PRINT *, 'Opening NC4 file...' |
---|
100 | ncretval = nf90_create(filename // ".nc4", & |
---|
101 | & OR(NF90_CLOBBER, NF90_HDF5), & |
---|
102 | & ncid) |
---|
103 | |
---|
104 | OPEN(IOUNIT_DUMP, file=filename, action='WRITE', status='REPLACE', form="unformatted", access="stream") |
---|
105 | |
---|
106 | |
---|
107 | |
---|
108 | |
---|
109 | |
---|
110 | |
---|
111 | CALL fpio(IOUNIT_DUMP, ncid, 'DUMP', cm_index) |
---|
112 | CLOSE(IOUNIT_DUMP) |
---|
113 | |
---|
114 | PRINT *, 'Closing NC4 file...' |
---|
115 | ncretval = nf90_close(ncid) |
---|
116 | |
---|
117 | CALL SYSTEM_CLOCK(millisecs_stop, count_rate, count_max) |
---|
118 | |
---|
119 | !PRINT *, 'Dump walltime secs: ', (millisecs_stop-millisecs_start)/1000.0 |
---|
120 | END SUBROUTINE fpmetbinary_dump |
---|
121 | |
---|
122 | SUBROUTINE fpmetbinary_load(filename, cm_index) |
---|
123 | CHARACTER(LEN=*), INTENT(IN) :: filename ! Full path for file |
---|
124 | INTEGER, INTENT(IN) :: cm_index ! Index of last dimension in |
---|
125 | ! most com_mod variables. |
---|
126 | ! Should be 1 or 2 |
---|
127 | |
---|
128 | INTEGER :: ncretval, ncid ! NetCDF func return value, file id |
---|
129 | |
---|
130 | INTEGER millisecs_start, millisecs_stop, count_rate, count_max |
---|
131 | |
---|
132 | CALL SYSTEM_CLOCK(millisecs_start, count_rate, count_max) |
---|
133 | |
---|
134 | print *, "Opening nc file for reading" |
---|
135 | ncretval = nf90_open(filename // ".nc4", NF90_NOWRITE, ncid) |
---|
136 | |
---|
137 | |
---|
138 | |
---|
139 | OPEN(IOUNIT_LOAD, file=filename, action='READ', status='OLD', form="unformatted", access="stream") |
---|
140 | CALL fpio(IOUNIT_LOAD, ncid, 'LOAD', cm_index) |
---|
141 | CLOSE(IOUNIT_LOAD) |
---|
142 | CALL SYSTEM_CLOCK(millisecs_stop, count_rate, count_max) |
---|
143 | !PRINT *, 'Load walltime secs: ', (millisecs_stop-millisecs_start)/1000.0 |
---|
144 | END SUBROUTINE fpmetbinary_load |
---|
145 | |
---|
146 | SUBROUTINE fpmetbinary_zero(cm_index) |
---|
147 | INTEGER, INTENT(IN) :: cm_index ! Index of last dimension in |
---|
148 | ! most com_mod variables. |
---|
149 | ! Should be 1 or 2 |
---|
150 | |
---|
151 | |
---|
152 | ! Zeroes out, in local datastructures, the values dumped/loaded |
---|
153 | ! This was written primarily as a testing mechanism. |
---|
154 | ! The lines here correspond to READ and WRITE in the dump/load routines |
---|
155 | |
---|
156 | ! Scalar values |
---|
157 | nx=0.0; ny=0.0; nxmin1=0.0; nymin1=0.0; nxfield=0.0 |
---|
158 | nuvz=0.0; nwz=0.0; nz=0.0; nmixz=0.0; nlev_ec=0.0 |
---|
159 | dx=0.0; dy=0.0; xlon0=0.0; ylat0=0.0; dxconst=0.0; dyconst=0.0 |
---|
160 | |
---|
161 | ! Fixed fields, static in time |
---|
162 | oro=0.0; excessoro=0.0; lsm=0.0; xlanduse=0.0; height=0.0 |
---|
163 | |
---|
164 | ! 3d fields |
---|
165 | uu(:,:,:,cm_index) = 0.0 |
---|
166 | vv(:,:,:,cm_index) = 0.0 |
---|
167 | uupol(:,:,:,cm_index) = 0.0 |
---|
168 | vvpol(:,:,:,cm_index) = 0.0 |
---|
169 | ww(:,:,:,cm_index) = 0.0 |
---|
170 | tt(:,:,:,cm_index) = 0.0 |
---|
171 | qv(:,:,:,cm_index) = 0.0 |
---|
172 | pv(:,:,:,cm_index) = 0.0 |
---|
173 | rho(:,:,:,cm_index) = 0.0 |
---|
174 | drhodz(:,:,:,cm_index) = 0.0 |
---|
175 | tth(:,:,:,cm_index) = 0.0 |
---|
176 | qvh(:,:,:,cm_index) = 0.0 |
---|
177 | pplev(:,:,:,cm_index) = 0.0 |
---|
178 | clouds(:,:,:,cm_index) = 0.0 |
---|
179 | cloudsh(:,:,cm_index) = 0.0 |
---|
180 | |
---|
181 | ! 2d fields |
---|
182 | ps(:,:,:,cm_index) = 0.0 |
---|
183 | sd(:,:,:,cm_index) = 0.0 |
---|
184 | msl(:,:,:,cm_index) = 0.0 |
---|
185 | tcc(:,:,:,cm_index) = 0.0 |
---|
186 | u10(:,:,:,cm_index) = 0.0 |
---|
187 | v10(:,:,:,cm_index) = 0.0 |
---|
188 | tt2(:,:,:,cm_index) = 0.0 |
---|
189 | td2(:,:,:,cm_index) = 0.0 |
---|
190 | lsprec(:,:,:,cm_index) = 0.0 |
---|
191 | convprec(:,:,:,cm_index) = 0.0 |
---|
192 | sshf(:,:,:,cm_index) = 0.0 |
---|
193 | ssr(:,:,:,cm_index) = 0.0 |
---|
194 | surfstr(:,:,:,cm_index) = 0.0 |
---|
195 | ustar(:,:,:,cm_index) = 0.0 |
---|
196 | wstar(:,:,:,cm_index) = 0.0 |
---|
197 | hmix(:,:,:,cm_index) = 0.0 |
---|
198 | tropopause(:,:,:,cm_index) = 0.0 |
---|
199 | oli(:,:,:,cm_index) = 0.0 |
---|
200 | diffk(:,:,:,cm_index) = 0.0 |
---|
201 | vdep(:,:,:,cm_index) = 0.0 |
---|
202 | |
---|
203 | ! 1d fields |
---|
204 | z0(:) = 0.0 |
---|
205 | akm(:) = 0.0 |
---|
206 | bkm(:) = 0.0 |
---|
207 | akz(:) = 0.0 |
---|
208 | bkz(:) = 0.0 |
---|
209 | aknew(:) = 0.0 |
---|
210 | bknew(:) = 0.0 |
---|
211 | |
---|
212 | ! Nested, scalar values (for each nest) |
---|
213 | nxn(:) = 0.0 |
---|
214 | nyn(:) = 0.0 |
---|
215 | dxn(:) = 0.0 |
---|
216 | dyn(:) = 0.0 |
---|
217 | xlon0n(:) = 0.0 |
---|
218 | ylat0n(:) = 0.0 |
---|
219 | |
---|
220 | ! Nested fields, static in time |
---|
221 | oron=0.0; excessoron=0.0; lsmn=0.0; xlandusen=0.0 |
---|
222 | |
---|
223 | ! 3d nested fields |
---|
224 | uun(:,:,:,cm_index,:) = 0.0 |
---|
225 | wwn(:,:,:,cm_index,:) = 0.0 |
---|
226 | ttn(:,:,:,cm_index,:) = 0.0 |
---|
227 | qvn(:,:,:,cm_index,:) = 0.0 |
---|
228 | pvn(:,:,:,cm_index,:) = 0.0 |
---|
229 | cloudsn(:,:,:,cm_index,:) = 0.0 |
---|
230 | cloudsnh(:,:,cm_index,:) = 0.0 |
---|
231 | rhon(:,:,:,cm_index,:) = 0.0 |
---|
232 | drhodzn(:,:,:,cm_index,:) = 0.0 |
---|
233 | tthn(:,:,:,cm_index,:) = 0.0 |
---|
234 | qvhn(:,:,:,cm_index,:) = 0.0 |
---|
235 | |
---|
236 | ! 2d nested fields |
---|
237 | psn(:,:,:,cm_index,:) = 0.0 |
---|
238 | sdn(:,:,:,cm_index,:) = 0.0 |
---|
239 | msln(:,:,:,cm_index,:) = 0.0 |
---|
240 | tccn(:,:,:,cm_index,:) = 0.0 |
---|
241 | u10n(:,:,:,cm_index,:) = 0.0 |
---|
242 | v10n(:,:,:,cm_index,:) = 0.0 |
---|
243 | tt2n(:,:,:,cm_index,:) = 0.0 |
---|
244 | td2n(:,:,:,cm_index,:) = 0.0 |
---|
245 | lsprecn(:,:,:,cm_index,:) = 0.0 |
---|
246 | convprecn(:,:,:,cm_index,:) = 0.0 |
---|
247 | sshfn(:,:,:,cm_index,:) = 0.0 |
---|
248 | ssrn(:,:,:,cm_index,:) = 0.0 |
---|
249 | surfstrn(:,:,:,cm_index,:) = 0.0 |
---|
250 | ustarn(:,:,:,cm_index,:) = 0.0 |
---|
251 | wstarn(:,:,:,cm_index,:) = 0.0 |
---|
252 | hmixn(:,:,:,cm_index,:) = 0.0 |
---|
253 | tropopausen(:,:,:,cm_index,:) = 0.0 |
---|
254 | olin(:,:,:,cm_index,:) = 0.0 |
---|
255 | diffkn(:,:,:,cm_index,:) = 0.0 |
---|
256 | vdepn(:,:,:,cm_index,:) = 0.0 |
---|
257 | |
---|
258 | ! Auxiliary variables for nests |
---|
259 | xresoln(:) = 0.0 |
---|
260 | yresoln(:) = 0.0 |
---|
261 | xln(:) = 0.0 |
---|
262 | yln(:) = 0.0 |
---|
263 | xrn(:) = 0.0 |
---|
264 | yrn(:) = 0.0 |
---|
265 | |
---|
266 | ! Variables for polar stereographic projection |
---|
267 | xglobal=.FALSE.; sglobal=.FALSE.; nglobal=.FALSE. |
---|
268 | switchnorthg=0.0; switchsouthg=0.0 |
---|
269 | southpolemap(:) = 0.0 |
---|
270 | northpolemap(:) = 0.0 |
---|
271 | |
---|
272 | ! Variables declared in conv_mod (convection) |
---|
273 | pconv(:) = 0.0 |
---|
274 | phconv(:) = 0.0 |
---|
275 | dpr(:) = 0.0 |
---|
276 | pconv_hpa(:) = 0.0 |
---|
277 | phconv_hpa(:) = 0.0 |
---|
278 | ft(:) = 0.0 |
---|
279 | fq(:) = 0.0 |
---|
280 | fmass(:,:) = 0.0 |
---|
281 | sub(:) = 0.0 |
---|
282 | fmassfrac(:,:) = 0.0 |
---|
283 | cbaseflux(:,:) = 0.0 |
---|
284 | cbasefluxn(:,:,:) = 0.0 |
---|
285 | tconv(:) = 0.0 |
---|
286 | qconv(:) = 0.0 |
---|
287 | qsconv(:) = 0.0 |
---|
288 | psconv=0.0; tt2conv=0.0; td2conv=0.0 |
---|
289 | nconvlev=0.0; nconvtop=0.0 |
---|
290 | |
---|
291 | END SUBROUTINE fpmetbinary_zero |
---|
292 | |
---|
293 | SUBROUTINE fpio(iounit, ncid, op, cm_index) |
---|
294 | IMPLICIT NONE |
---|
295 | INTEGER, INTENT(IN) :: ncid ! NetCDF file id |
---|
296 | INTEGER, INTENT(IN) :: iounit |
---|
297 | CHARACTER(LEN=4), INTENT(IN) :: op ! Operation - DUMP or LOAD |
---|
298 | INTEGER, INTENT(IN) :: cm_index ! Index of last dimension in |
---|
299 | ! most com_mod variables. |
---|
300 | ! Should be 1 or 2 |
---|
301 | |
---|
302 | |
---|
303 | INTEGER :: ncret ! Return value from NetCDF calls |
---|
304 | INTEGER :: ncvarid ! NetCDF variable ID |
---|
305 | |
---|
306 | INTEGER :: nxmax_dimid, nymax_dimid, nzmax_dimid, nuvzmax_dimid, nwzmax_dimid |
---|
307 | |
---|
308 | INTEGER, DIMENSION(1) :: dim1dids ! Dimension IDs for 1D arrays |
---|
309 | INTEGER, DIMENSION(2) :: dim2dids ! Dimension IDs for 2D arrays |
---|
310 | INTEGER, DIMENSION(3) :: dim3dids ! Dimension IDs for 3D arrays |
---|
311 | |
---|
312 | ! These are used when loading in dimensions from NC file |
---|
313 | CHARACTER(LEN=NF90_MAX_NAME) :: nxmax_dimname, nymax_dimname, nzmax_dimname, & |
---|
314 | & nuvzmax_dimname, nwzmax_dimname |
---|
315 | |
---|
316 | ! These are temporary variables, used in the LOAD option, for |
---|
317 | ! comparing against the current values in FLEXPART of nxmax, nymax, ... |
---|
318 | INTEGER :: temp_nxmax, temp_nymax, temp_nzmax, & |
---|
319 | & temp_nuvzmax, temp_nwzmax |
---|
320 | |
---|
321 | CHARACTER(LEN=12) :: temp_preproc_format_version_str |
---|
322 | |
---|
323 | CHARACTER(LEN=128) :: errmesg |
---|
324 | |
---|
325 | INTEGER, PARAMETER :: DEF_LEVEL = 3 |
---|
326 | |
---|
327 | if (op == 'DUMP') THEN |
---|
328 | |
---|
329 | |
---|
330 | ! Write the preprocessing format version string |
---|
331 | WRITE (iounit) PREPROC_FORMAT_VERSION_STR |
---|
332 | |
---|
333 | ! Write the compiled max dimensions from par_mod - these are |
---|
334 | ! not meant to be reassigned during a LOAD, but used as "header" |
---|
335 | ! information to provide the structure of arrays |
---|
336 | WRITE (iounit) nxmax, nymax, nzmax, nuvzmax, nwzmax |
---|
337 | |
---|
338 | ncret = nf90_def_dim(ncid, 'nxmax', nxmax, nxmax_dimid) |
---|
339 | ncret = nf90_def_dim(ncid, 'nymax', nymax, nymax_dimid) |
---|
340 | ncret = nf90_def_dim(ncid, 'nzmax', nzmax, nzmax_dimid) |
---|
341 | ncret = nf90_def_dim(ncid, 'nuvzmax', nuvzmax, nuvzmax_dimid) |
---|
342 | ncret = nf90_def_dim(ncid, 'nwzmax', nwzmax, nwzmax_dimid) |
---|
343 | |
---|
344 | |
---|
345 | |
---|
346 | ! Scalar values |
---|
347 | WRITE(iounit) nx, ny, nxmin1, nymin1, nxfield |
---|
348 | WRITE(iounit) nuvz, nwz, nz, nmixz, nlev_ec |
---|
349 | WRITE(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst |
---|
350 | |
---|
351 | ncret = nf90_def_var(ncid, 'nx', NF90_INT, ncvarid) |
---|
352 | ncret = nf90_put_var(ncid, ncvarid, nx) |
---|
353 | |
---|
354 | ncret = nf90_def_var(ncid, 'ny', NF90_INT, ncvarid) |
---|
355 | ncret = nf90_put_var(ncid, ncvarid, ny) |
---|
356 | |
---|
357 | ncret = nf90_def_var(ncid, 'nxmin1', NF90_INT, ncvarid) |
---|
358 | ncret = nf90_put_var(ncid, ncvarid, nxmin1) |
---|
359 | |
---|
360 | ncret = nf90_def_var(ncid, 'nymin1', NF90_INT, ncvarid) |
---|
361 | ncret = nf90_put_var(ncid, ncvarid, nymin1) |
---|
362 | |
---|
363 | ncret = nf90_def_var(ncid, 'nxfield', NF90_INT, ncvarid) |
---|
364 | ncret = nf90_put_var(ncid, ncvarid, nxfield) |
---|
365 | |
---|
366 | ncret = nf90_def_var(ncid, 'nuvz', NF90_INT, ncvarid) |
---|
367 | ncret = nf90_put_var(ncid, ncvarid, nuvz) |
---|
368 | |
---|
369 | ncret = nf90_def_var(ncid, 'nwz', NF90_INT, ncvarid) |
---|
370 | ncret = nf90_put_var(ncid, ncvarid, nwz) |
---|
371 | |
---|
372 | ncret = nf90_def_var(ncid, 'nz', NF90_INT, ncvarid) |
---|
373 | ncret = nf90_put_var(ncid, ncvarid, nz) |
---|
374 | |
---|
375 | ncret = nf90_def_var(ncid, 'nmixz', NF90_INT, ncvarid) |
---|
376 | ncret = nf90_put_var(ncid, ncvarid, nmixz) |
---|
377 | |
---|
378 | ncret = nf90_def_var(ncid, 'nlev_', NF90_INT, ncvarid) |
---|
379 | ncret = nf90_put_var(ncid, ncvarid, nlev_ec) |
---|
380 | |
---|
381 | ncret = nf90_def_var(ncid, 'dx', NF90_FLOAT, ncvarid) |
---|
382 | ncret = nf90_put_var(ncid, ncvarid, dx) |
---|
383 | |
---|
384 | ncret = nf90_def_var(ncid, 'dy', NF90_FLOAT, ncvarid) |
---|
385 | ncret = nf90_put_var(ncid, ncvarid, dy) |
---|
386 | |
---|
387 | ncret = nf90_def_var(ncid, 'xlon0', NF90_FLOAT, ncvarid) |
---|
388 | ncret = nf90_put_var(ncid, ncvarid, xlon0) |
---|
389 | |
---|
390 | ncret = nf90_def_var(ncid, 'ylat0', NF90_FLOAT, ncvarid) |
---|
391 | ncret = nf90_put_var(ncid, ncvarid, ylat0) |
---|
392 | |
---|
393 | ncret = nf90_def_var(ncid, 'dxconst', NF90_FLOAT, ncvarid) |
---|
394 | ncret = nf90_put_var(ncid, ncvarid, dxconst) |
---|
395 | |
---|
396 | ncret = nf90_def_var(ncid, 'dyconst', NF90_FLOAT, ncvarid) |
---|
397 | ncret = nf90_put_var(ncid, ncvarid, dyconst) |
---|
398 | |
---|
399 | |
---|
400 | |
---|
401 | ! Fixed fields, static in time |
---|
402 | WRITE(iounit) oro, excessoro, lsm, xlanduse, height |
---|
403 | |
---|
404 | dim2dids = (/nxmax_dimid, nymax_dimid/) |
---|
405 | |
---|
406 | ncret = nf90_def_var(ncid, 'oro', NF90_FLOAT, & |
---|
407 | & dim2dids, ncvarid) |
---|
408 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
409 | & shuffle=0, & |
---|
410 | & deflate=1, & |
---|
411 | & deflate_level=DEF_LEVEL) |
---|
412 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
413 | & oro(0:nxmax-1, 0:nymax-1)) |
---|
414 | |
---|
415 | ncret = nf90_def_var(ncid, 'excessoro', NF90_FLOAT, & |
---|
416 | & dim2dids, ncvarid) |
---|
417 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
418 | & shuffle=0, & |
---|
419 | & deflate=1, & |
---|
420 | & deflate_level=DEF_LEVEL) |
---|
421 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
422 | & excessoro(0:nxmax-1, 0:nymax-1)) |
---|
423 | |
---|
424 | ncret = nf90_def_var(ncid, 'lsm', NF90_FLOAT, & |
---|
425 | & dim2dids, ncvarid) |
---|
426 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
427 | & shuffle=0, & |
---|
428 | & deflate=1, & |
---|
429 | & deflate_level=DEF_LEVEL) |
---|
430 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
431 | & lsm(0:nxmax-1, 0:nymax-1)) |
---|
432 | |
---|
433 | dim3dids = (/nxmax_dimid, nymax_dimid, numclass/) |
---|
434 | ! numclass comes from par_mod - number of land use classes |
---|
435 | ncret = nf90_def_var(ncid, 'xlanduse', NF90_FLOAT, & |
---|
436 | & dim3dids, ncvarid) |
---|
437 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
438 | & shuffle=0, & |
---|
439 | & deflate=1, & |
---|
440 | & deflate_level=DEF_LEVEL) |
---|
441 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
442 | & xlanduse(0:nxmax-1, 0:nymax-1, 1:numclass)) |
---|
443 | |
---|
444 | dim1dids = (/nzmax_dimid/) |
---|
445 | ncret = nf90_def_var(ncid, 'height', NF90_FLOAT, & |
---|
446 | & dim1dids, ncvarid) |
---|
447 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
448 | & shuffle=0, & |
---|
449 | & deflate=1, & |
---|
450 | & deflate_level=DEF_LEVEL) |
---|
451 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
452 | & height(1:nzmax)) |
---|
453 | |
---|
454 | |
---|
455 | |
---|
456 | |
---|
457 | ! 3d fields |
---|
458 | WRITE(iounit) uu(:,:,:,cm_index) |
---|
459 | WRITE(iounit) vv(:,:,:,cm_index) |
---|
460 | WRITE(iounit) uupol(:,:,:,cm_index) |
---|
461 | WRITE(iounit) vvpol(:,:,:,cm_index) |
---|
462 | WRITE(iounit) ww(:,:,:,cm_index) |
---|
463 | WRITE(iounit) tt(:,:,:,cm_index) |
---|
464 | WRITE(iounit) qv(:,:,:,cm_index) |
---|
465 | WRITE(iounit) pv(:,:,:,cm_index) |
---|
466 | WRITE(iounit) rho(:,:,:,cm_index) |
---|
467 | WRITE(iounit) drhodz(:,:,:,cm_index) |
---|
468 | WRITE(iounit) tth(:,:,:,cm_index) |
---|
469 | WRITE(iounit) qvh(:,:,:,cm_index) |
---|
470 | WRITE(iounit) pplev(:,:,:,cm_index) |
---|
471 | WRITE(iounit) clouds(:,:,:,cm_index) |
---|
472 | WRITE(iounit) cloudsh(:,:,cm_index) |
---|
473 | |
---|
474 | dim3dids = (/nxmax_dimid, nymax_dimid, nzmax_dimid/) |
---|
475 | |
---|
476 | ncret = nf90_def_var(ncid, 'uu', NF90_FLOAT, & |
---|
477 | & dim3dids, ncvarid) |
---|
478 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
479 | & shuffle=0, & |
---|
480 | & deflate=1, & |
---|
481 | & deflate_level=DEF_LEVEL) |
---|
482 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
483 | & uu(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)) |
---|
484 | |
---|
485 | ncret = nf90_def_var(ncid, 'vv', NF90_FLOAT, & |
---|
486 | & dim3dids, ncvarid) |
---|
487 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
488 | & shuffle=0, & |
---|
489 | & deflate=1, & |
---|
490 | & deflate_level=DEF_LEVEL) |
---|
491 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
492 | & vv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)) |
---|
493 | |
---|
494 | ncret = nf90_def_var(ncid, 'uupol', NF90_FLOAT, & |
---|
495 | & dim3dids, ncvarid) |
---|
496 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
497 | & shuffle=0, & |
---|
498 | & deflate=1, & |
---|
499 | & deflate_level=DEF_LEVEL) |
---|
500 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
501 | & uupol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)) |
---|
502 | |
---|
503 | ncret = nf90_def_var(ncid, 'vvpol', NF90_FLOAT, & |
---|
504 | & dim3dids, ncvarid) |
---|
505 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
506 | & shuffle=0, & |
---|
507 | & deflate=1, & |
---|
508 | & deflate_level=DEF_LEVEL) |
---|
509 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
510 | & vvpol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)) |
---|
511 | |
---|
512 | ncret = nf90_def_var(ncid, 'ww', NF90_FLOAT, & |
---|
513 | & dim3dids, ncvarid) |
---|
514 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
515 | & shuffle=0, & |
---|
516 | & deflate=1, & |
---|
517 | & deflate_level=DEF_LEVEL) |
---|
518 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
519 | & ww(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)) |
---|
520 | |
---|
521 | ncret = nf90_def_var(ncid, 'tt', NF90_FLOAT, & |
---|
522 | & dim3dids, ncvarid) |
---|
523 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
524 | & shuffle=0, & |
---|
525 | & deflate=1, & |
---|
526 | & deflate_level=DEF_LEVEL) |
---|
527 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
528 | & tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)) |
---|
529 | |
---|
530 | ncret = nf90_def_var(ncid, 'qv', NF90_FLOAT, & |
---|
531 | & dim3dids, ncvarid) |
---|
532 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
533 | & shuffle=0, & |
---|
534 | & deflate=1, & |
---|
535 | & deflate_level=DEF_LEVEL) |
---|
536 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
537 | & qv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)) |
---|
538 | |
---|
539 | ncret = nf90_def_var(ncid, 'pv', NF90_FLOAT, & |
---|
540 | & dim3dids, ncvarid) |
---|
541 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
542 | & shuffle=0, & |
---|
543 | & deflate=1, & |
---|
544 | & deflate_level=DEF_LEVEL) |
---|
545 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
546 | & pv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)) |
---|
547 | |
---|
548 | ncret = nf90_def_var(ncid, 'rho', NF90_FLOAT, & |
---|
549 | & dim3dids, ncvarid) |
---|
550 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
551 | & shuffle=0, & |
---|
552 | & deflate=1, & |
---|
553 | & deflate_level=DEF_LEVEL) |
---|
554 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
555 | & rho(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)) |
---|
556 | |
---|
557 | ncret = nf90_def_var(ncid, 'drhodz', NF90_FLOAT, & |
---|
558 | & dim3dids, ncvarid) |
---|
559 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
560 | & shuffle=0, & |
---|
561 | & deflate=1, & |
---|
562 | & deflate_level=DEF_LEVEL) |
---|
563 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
564 | & drhodz(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)) |
---|
565 | |
---|
566 | ncret = nf90_def_var(ncid, 'clouds', NF90_BYTE, & |
---|
567 | & dim3dids, ncvarid) |
---|
568 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
569 | & shuffle=0, & |
---|
570 | & deflate=1, & |
---|
571 | & deflate_level=DEF_LEVEL) |
---|
572 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
573 | & clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)) |
---|
574 | |
---|
575 | |
---|
576 | |
---|
577 | ! Note the change in z dimension for the following |
---|
578 | dim3dids = (/nxmax_dimid, nymax_dimid, nuvzmax_dimid/) |
---|
579 | |
---|
580 | ncret = nf90_def_var(ncid, 'tth', NF90_FLOAT, & |
---|
581 | & dim3dids, ncvarid) |
---|
582 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
583 | & shuffle=0, & |
---|
584 | & deflate=1, & |
---|
585 | & deflate_level=DEF_LEVEL) |
---|
586 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
587 | & tth(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index)) |
---|
588 | |
---|
589 | ncret = nf90_def_var(ncid, 'qvh', NF90_FLOAT, & |
---|
590 | & dim3dids, ncvarid) |
---|
591 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
592 | & shuffle=0, & |
---|
593 | & deflate=1, & |
---|
594 | & deflate_level=DEF_LEVEL) |
---|
595 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
596 | & qvh(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index)) |
---|
597 | |
---|
598 | ncret = nf90_def_var(ncid, 'pplev', NF90_FLOAT, & |
---|
599 | & dim3dids, ncvarid) |
---|
600 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
601 | & shuffle=0, & |
---|
602 | & deflate=1, & |
---|
603 | & deflate_level=DEF_LEVEL) |
---|
604 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
605 | & pplev(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index)) |
---|
606 | |
---|
607 | |
---|
608 | dim2dids = (/nxmax_dimid, nymax_dimid/) |
---|
609 | ncret = nf90_def_var(ncid, 'cloudsh', NF90_INT, & |
---|
610 | & dim2dids, ncvarid) |
---|
611 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
612 | & shuffle=0, & |
---|
613 | & deflate=1, & |
---|
614 | & deflate_level=DEF_LEVEL) |
---|
615 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
616 | & cloudsh(0:nxmax-1, 0:nymax-1, cm_index)) |
---|
617 | |
---|
618 | |
---|
619 | |
---|
620 | PRINT *, 'SUM(tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', & |
---|
621 | & SUM(tt(0:nxmax-1,0:nymax-1,1:nzmax, cm_index)) |
---|
622 | |
---|
623 | PRINT *, 'SUM(clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', & |
---|
624 | & SUM(clouds(0:nxmax-1,0:nymax-1,1:nzmax, cm_index)) |
---|
625 | |
---|
626 | |
---|
627 | |
---|
628 | ! 2d fields |
---|
629 | WRITE(iounit) ps(:,:,:,cm_index) |
---|
630 | WRITE(iounit) sd(:,:,:,cm_index) |
---|
631 | WRITE(iounit) msl(:,:,:,cm_index) |
---|
632 | WRITE(iounit) tcc(:,:,:,cm_index) |
---|
633 | WRITE(iounit) u10(:,:,:,cm_index) |
---|
634 | WRITE(iounit) v10(:,:,:,cm_index) |
---|
635 | WRITE(iounit) tt2(:,:,:,cm_index) |
---|
636 | WRITE(iounit) td2(:,:,:,cm_index) |
---|
637 | WRITE(iounit) lsprec(:,:,:,cm_index) |
---|
638 | WRITE(iounit) convprec(:,:,:,cm_index) |
---|
639 | WRITE(iounit) sshf(:,:,:,cm_index) |
---|
640 | WRITE(iounit) ssr(:,:,:,cm_index) |
---|
641 | WRITE(iounit) surfstr(:,:,:,cm_index) |
---|
642 | WRITE(iounit) ustar(:,:,:,cm_index) |
---|
643 | WRITE(iounit) wstar(:,:,:,cm_index) |
---|
644 | WRITE(iounit) hmix(:,:,:,cm_index) |
---|
645 | WRITE(iounit) tropopause(:,:,:,cm_index) |
---|
646 | WRITE(iounit) oli(:,:,:,cm_index) |
---|
647 | WRITE(iounit) diffk(:,:,:,cm_index) |
---|
648 | WRITE(iounit) vdep(:,:,:,cm_index) |
---|
649 | |
---|
650 | dim2dids = (/nxmax_dimid, nymax_dimid/) |
---|
651 | |
---|
652 | ncret = nf90_def_var(ncid, 'ps', NF90_FLOAT, & |
---|
653 | & dim2dids, ncvarid) |
---|
654 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
655 | & shuffle=0, & |
---|
656 | & deflate=1, & |
---|
657 | & deflate_level=DEF_LEVEL) |
---|
658 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
659 | & ps(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
660 | |
---|
661 | ncret = nf90_def_var(ncid, 'sd', NF90_FLOAT, & |
---|
662 | & dim2dids, ncvarid) |
---|
663 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
664 | & shuffle=0, & |
---|
665 | & deflate=1, & |
---|
666 | & deflate_level=DEF_LEVEL) |
---|
667 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
668 | & sd(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
669 | |
---|
670 | ncret = nf90_def_var(ncid, 'msl', NF90_FLOAT, & |
---|
671 | & dim2dids, ncvarid) |
---|
672 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
673 | & shuffle=0, & |
---|
674 | & deflate=1, & |
---|
675 | & deflate_level=DEF_LEVEL) |
---|
676 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
677 | & msl(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
678 | |
---|
679 | ncret = nf90_def_var(ncid, 'tcc', NF90_FLOAT, & |
---|
680 | & dim2dids, ncvarid) |
---|
681 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
682 | & shuffle=0, & |
---|
683 | & deflate=1, & |
---|
684 | & deflate_level=DEF_LEVEL) |
---|
685 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
686 | & tcc(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
687 | |
---|
688 | ncret = nf90_def_var(ncid, 'u10', NF90_FLOAT, & |
---|
689 | & dim2dids, ncvarid) |
---|
690 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
691 | & shuffle=0, & |
---|
692 | & deflate=1, & |
---|
693 | & deflate_level=DEF_LEVEL) |
---|
694 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
695 | & u10(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
696 | |
---|
697 | ncret = nf90_def_var(ncid, 'v10', NF90_FLOAT, & |
---|
698 | & dim2dids, ncvarid) |
---|
699 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
700 | & shuffle=0, & |
---|
701 | & deflate=1, & |
---|
702 | & deflate_level=DEF_LEVEL) |
---|
703 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
704 | & v10(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
705 | |
---|
706 | ncret = nf90_def_var(ncid, 'tt2', NF90_FLOAT, & |
---|
707 | & dim2dids, ncvarid) |
---|
708 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
709 | & shuffle=0, & |
---|
710 | & deflate=1, & |
---|
711 | & deflate_level=DEF_LEVEL) |
---|
712 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
713 | & tt2(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
714 | |
---|
715 | ncret = nf90_def_var(ncid, 'td2', NF90_FLOAT, & |
---|
716 | & dim2dids, ncvarid) |
---|
717 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
718 | & shuffle=0, & |
---|
719 | & deflate=1, & |
---|
720 | & deflate_level=DEF_LEVEL) |
---|
721 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
722 | & td2(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
723 | |
---|
724 | ncret = nf90_def_var(ncid, 'lsprec', NF90_FLOAT, & |
---|
725 | & dim2dids, ncvarid) |
---|
726 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
727 | & shuffle=0, & |
---|
728 | & deflate=1, & |
---|
729 | & deflate_level=DEF_LEVEL) |
---|
730 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
731 | & lsprec(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
732 | |
---|
733 | ncret = nf90_def_var(ncid, 'convprec', NF90_FLOAT, & |
---|
734 | & dim2dids, ncvarid) |
---|
735 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
736 | & shuffle=0, & |
---|
737 | & deflate=1, & |
---|
738 | & deflate_level=DEF_LEVEL) |
---|
739 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
740 | & convprec(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
741 | |
---|
742 | ncret = nf90_def_var(ncid, 'sshf', NF90_FLOAT, & |
---|
743 | & dim2dids, ncvarid) |
---|
744 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
745 | & shuffle=0, & |
---|
746 | & deflate=1, & |
---|
747 | & deflate_level=DEF_LEVEL) |
---|
748 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
749 | & sshf(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
750 | |
---|
751 | ncret = nf90_def_var(ncid, 'ssr', NF90_FLOAT, & |
---|
752 | & dim2dids, ncvarid) |
---|
753 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
754 | & shuffle=0, & |
---|
755 | & deflate=1, & |
---|
756 | & deflate_level=DEF_LEVEL) |
---|
757 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
758 | & ssr(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
759 | |
---|
760 | ncret = nf90_def_var(ncid, 'surfstr', NF90_FLOAT, & |
---|
761 | & dim2dids, ncvarid) |
---|
762 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
763 | & shuffle=0, & |
---|
764 | & deflate=1, & |
---|
765 | & deflate_level=DEF_LEVEL) |
---|
766 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
767 | & surfstr(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
768 | |
---|
769 | ncret = nf90_def_var(ncid, 'ustar', NF90_FLOAT, & |
---|
770 | & dim2dids, ncvarid) |
---|
771 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
772 | & shuffle=0, & |
---|
773 | & deflate=1, & |
---|
774 | & deflate_level=DEF_LEVEL) |
---|
775 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
776 | & ustar(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
777 | |
---|
778 | ncret = nf90_def_var(ncid, 'wstar', NF90_FLOAT, & |
---|
779 | & dim2dids, ncvarid) |
---|
780 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
781 | & shuffle=0, & |
---|
782 | & deflate=1, & |
---|
783 | & deflate_level=DEF_LEVEL) |
---|
784 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
785 | & wstar(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
786 | |
---|
787 | ncret = nf90_def_var(ncid, 'hmix', NF90_FLOAT, & |
---|
788 | & dim2dids, ncvarid) |
---|
789 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
790 | & shuffle=0, & |
---|
791 | & deflate=1, & |
---|
792 | & deflate_level=DEF_LEVEL) |
---|
793 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
794 | & hmix(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
795 | |
---|
796 | ncret = nf90_def_var(ncid, 'tropopause', NF90_FLOAT, & |
---|
797 | & dim2dids, ncvarid) |
---|
798 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
799 | & shuffle=0, & |
---|
800 | & deflate=1, & |
---|
801 | & deflate_level=DEF_LEVEL) |
---|
802 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
803 | & tropopause(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
804 | |
---|
805 | ncret = nf90_def_var(ncid, 'oli', NF90_FLOAT, & |
---|
806 | & dim2dids, ncvarid) |
---|
807 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
808 | & shuffle=0, & |
---|
809 | & deflate=1, & |
---|
810 | & deflate_level=DEF_LEVEL) |
---|
811 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
812 | & oli(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
813 | |
---|
814 | ncret = nf90_def_var(ncid, 'diffk', NF90_FLOAT, & |
---|
815 | & dim2dids, ncvarid) |
---|
816 | ncret = nf90_def_var_deflate(ncid, ncvarid, & |
---|
817 | & shuffle=0, & |
---|
818 | & deflate=1, & |
---|
819 | & deflate_level=DEF_LEVEL) |
---|
820 | ncret = nf90_put_var(ncid, ncvarid, & |
---|
821 | & diffk(0:nxmax-1, 0:nymax-1, 1, cm_index)) |
---|
822 | |
---|
823 | |
---|
824 | |
---|
825 | PRINT *, 'SUM(ps(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', & |
---|
826 | & SUM(ps(0:nxmax-1,0:nymax-1,1, cm_index)) |
---|
827 | |
---|
828 | PRINT *, 'SUM(wstar(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', & |
---|
829 | & SUM(wstar(0:nxmax-1,0:nymax-1,1, cm_index)) |
---|
830 | |
---|
831 | |
---|
832 | |
---|
833 | |
---|
834 | |
---|
835 | ! 1d fields |
---|
836 | WRITE(iounit) z0(:) |
---|
837 | WRITE(iounit) akm(:) |
---|
838 | WRITE(iounit) bkm(:) |
---|
839 | WRITE(iounit) akz(:) |
---|
840 | WRITE(iounit) bkz(:) |
---|
841 | WRITE(iounit) aknew(:) |
---|
842 | WRITE(iounit) bknew(:) |
---|
843 | |
---|
844 | ! Nested, scalar values (for each nest) |
---|
845 | WRITE(iounit) nxn(:) |
---|
846 | WRITE(iounit) nyn(:) |
---|
847 | WRITE(iounit) dxn(:) |
---|
848 | WRITE(iounit) dyn(:) |
---|
849 | WRITE(iounit) xlon0n(:) |
---|
850 | WRITE(iounit) ylat0n(:) |
---|
851 | |
---|
852 | ! Nested fields, static over time |
---|
853 | WRITE(iounit) oron, excessoron, lsmn, xlandusen |
---|
854 | |
---|
855 | ! 3d nested fields |
---|
856 | WRITE(iounit) uun(:,:,:,cm_index,:) |
---|
857 | WRITE(iounit) vvn(:,:,:,cm_index,:) |
---|
858 | WRITE(iounit) wwn(:,:,:,cm_index,:) |
---|
859 | WRITE(iounit) ttn(:,:,:,cm_index,:) |
---|
860 | WRITE(iounit) qvn(:,:,:,cm_index,:) |
---|
861 | WRITE(iounit) pvn(:,:,:,cm_index,:) |
---|
862 | WRITE(iounit) cloudsn(:,:,:,cm_index,:) |
---|
863 | WRITE(iounit) cloudsnh(:,:,cm_index,:) |
---|
864 | WRITE(iounit) rhon(:,:,:,cm_index,:) |
---|
865 | WRITE(iounit) drhodzn(:,:,:,cm_index,:) |
---|
866 | WRITE(iounit) tthn(:,:,:,cm_index,:) |
---|
867 | WRITE(iounit) qvhn(:,:,:,cm_index,:) |
---|
868 | |
---|
869 | ! 2d nested fields |
---|
870 | WRITE(iounit) psn(:,:,:,cm_index,:) |
---|
871 | WRITE(iounit) sdn(:,:,:,cm_index,:) |
---|
872 | WRITE(iounit) msln(:,:,:,cm_index,:) |
---|
873 | WRITE(iounit) tccn(:,:,:,cm_index,:) |
---|
874 | WRITE(iounit) u10n(:,:,:,cm_index,:) |
---|
875 | WRITE(iounit) v10n(:,:,:,cm_index,:) |
---|
876 | WRITE(iounit) tt2n(:,:,:,cm_index,:) |
---|
877 | WRITE(iounit) td2n(:,:,:,cm_index,:) |
---|
878 | WRITE(iounit) lsprecn(:,:,:,cm_index,:) |
---|
879 | WRITE(iounit) convprecn(:,:,:,cm_index,:) |
---|
880 | WRITE(iounit) sshfn(:,:,:,cm_index,:) |
---|
881 | WRITE(iounit) ssrn(:,:,:,cm_index,:) |
---|
882 | WRITE(iounit) surfstrn(:,:,:,cm_index,:) |
---|
883 | WRITE(iounit) ustarn(:,:,:,cm_index,:) |
---|
884 | WRITE(iounit) wstarn(:,:,:,cm_index,:) |
---|
885 | WRITE(iounit) hmixn(:,:,:,cm_index,:) |
---|
886 | WRITE(iounit) tropopausen(:,:,:,cm_index,:) |
---|
887 | WRITE(iounit) olin(:,:,:,cm_index,:) |
---|
888 | WRITE(iounit) diffkn(:,:,:,cm_index,:) |
---|
889 | WRITE(iounit) vdepn(:,:,:,cm_index,:) |
---|
890 | |
---|
891 | ! Auxiliary variables for nests |
---|
892 | WRITE(iounit) xresoln(:) |
---|
893 | WRITE(iounit) yresoln(:) |
---|
894 | WRITE(iounit) xln(:) |
---|
895 | WRITE(iounit) yln(:) |
---|
896 | WRITE(iounit) xrn(:) |
---|
897 | WRITE(iounit) yrn(:) |
---|
898 | |
---|
899 | ! Variables for polar stereographic projection |
---|
900 | WRITE(iounit) xglobal, sglobal, nglobal |
---|
901 | WRITE(iounit) switchnorthg, switchsouthg |
---|
902 | WRITE(iounit) southpolemap(:) |
---|
903 | WRITE(iounit) northpolemap(:) |
---|
904 | |
---|
905 | ! Variables declared in conv_mod (convection) |
---|
906 | WRITE(iounit) pconv(:) |
---|
907 | WRITE(iounit) phconv(:) |
---|
908 | WRITE(iounit) dpr(:) |
---|
909 | WRITE(iounit) pconv_hpa(:) |
---|
910 | WRITE(iounit) phconv_hpa(:) |
---|
911 | WRITE(iounit) ft(:) |
---|
912 | WRITE(iounit) fq(:) |
---|
913 | WRITE(iounit) fmass(:,:) |
---|
914 | WRITE(iounit) sub(:) |
---|
915 | WRITE(iounit) fmassfrac(:,:) |
---|
916 | WRITE(iounit) cbaseflux(:,:) |
---|
917 | WRITE(iounit) cbasefluxn(:,:,:) |
---|
918 | WRITE(iounit) tconv(:) |
---|
919 | WRITE(iounit) qconv(:) |
---|
920 | WRITE(iounit) qsconv(:) |
---|
921 | WRITE(iounit) psconv, tt2conv, td2conv |
---|
922 | WRITE(iounit) nconvlev, nconvtop |
---|
923 | |
---|
924 | ELSE IF (op == 'LOAD') THEN |
---|
925 | |
---|
926 | ! Read the preprocessed format version string and insure it |
---|
927 | ! matches this version |
---|
928 | READ (iounit) temp_preproc_format_version_str |
---|
929 | PRINT *, 'Reading preprocessed file format version: ', & |
---|
930 | & temp_preproc_format_version_str |
---|
931 | |
---|
932 | IF (TRIM(temp_preproc_format_version_str) == & |
---|
933 | & TRIM(PREPROC_FORMAT_VERSION_STR)) THEN |
---|
934 | CONTINUE |
---|
935 | ELSE |
---|
936 | ! PRINT *, '' GK: causes relocation truncated to fit: R_X86_64_32 |
---|
937 | PRINT *, 'Inconsistent preprocessing format version' |
---|
938 | PRINT *, 'Expected Version: ', PREPROC_FORMAT_VERSION_STR |
---|
939 | PRINT *, 'Detected Version: ', temp_preproc_format_version_str |
---|
940 | ! PRINT *, '' |
---|
941 | STOP |
---|
942 | END IF |
---|
943 | |
---|
944 | ! Read the compiled max dimensions that were dumped from par_mod |
---|
945 | ! when creating the fp file, so that we can compare against |
---|
946 | ! current FLEXPART dimensions - they need to be the same, or else |
---|
947 | ! we abort. |
---|
948 | READ (iounit) temp_nxmax, temp_nymax, temp_nzmax, & |
---|
949 | & temp_nuvzmax, temp_nwzmax |
---|
950 | |
---|
951 | ! Get dimensions |
---|
952 | ncret = nf90_inq_dimid(ncid, 'nxmax', nxmax_dimid) |
---|
953 | ncret = nf90_inquire_dimension(ncid, nxmax_dimid, nxmax_dimname, & |
---|
954 | & temp_nxmax) |
---|
955 | PRINT *, 'temp_nxmax: ', temp_nxmax |
---|
956 | |
---|
957 | ncret = nf90_inq_dimid(ncid, 'nymax', nymax_dimid) |
---|
958 | ncret = nf90_inquire_dimension(ncid, nymax_dimid, nymax_dimname, & |
---|
959 | & temp_nymax) |
---|
960 | PRINT *, 'temp_nymax: ', temp_nymax |
---|
961 | |
---|
962 | ncret = nf90_inq_dimid(ncid, 'nzmax', nzmax_dimid) |
---|
963 | ncret = nf90_inquire_dimension(ncid, nzmax_dimid, nzmax_dimname, & |
---|
964 | & temp_nzmax) |
---|
965 | PRINT *, 'temp_nzmax: ', temp_nzmax |
---|
966 | |
---|
967 | ncret = nf90_inq_dimid(ncid, 'nuvzmax', nuvzmax_dimid) |
---|
968 | ncret = nf90_inquire_dimension(ncid, nuvzmax_dimid, nuvzmax_dimname, & |
---|
969 | & temp_nuvzmax) |
---|
970 | PRINT *, 'temp_nuvzmax: ', temp_nuvzmax |
---|
971 | |
---|
972 | ncret = nf90_inq_dimid(ncid, 'nwzmax', nwzmax_dimid) |
---|
973 | ncret = nf90_inquire_dimension(ncid, nwzmax_dimid, nwzmax_dimname, & |
---|
974 | & temp_nwzmax) |
---|
975 | PRINT *, 'temp_nwzmax: ', temp_nwzmax |
---|
976 | |
---|
977 | |
---|
978 | |
---|
979 | IF ( (temp_nxmax == nxmax) .AND. (temp_nymax == nymax) .AND. & |
---|
980 | & (temp_nzmax == nzmax) .AND. & |
---|
981 | & (temp_nuvzmax == nuvzmax) .AND. & |
---|
982 | & (temp_nwzmax == nwzmax) ) THEN |
---|
983 | CONTINUE |
---|
984 | ELSE |
---|
985 | PRINT *, 'Incompatible dimensions between fp file and current FLEXPART!' |
---|
986 | ! PRINT *, '' |
---|
987 | PRINT *, ' FP file Compiled FP' |
---|
988 | PRINT *, 'nxmax: ', temp_nxmax, ' ', nxmax |
---|
989 | PRINT *, 'nymax: ', temp_nymax, ' ', nymax |
---|
990 | PRINT *, 'nzmax: ', temp_nzmax, ' ', nzmax |
---|
991 | PRINT *, 'nuvzmax: ', temp_nuvzmax, ' ', nuvzmax |
---|
992 | PRINT *, 'nwzmax: ', temp_nwzmax, ' ', nwzmax |
---|
993 | ! PRINT *, '' |
---|
994 | STOP |
---|
995 | END IF |
---|
996 | |
---|
997 | |
---|
998 | ! Scalar values |
---|
999 | READ(iounit) nx, ny, nxmin1, nymin1, nxfield |
---|
1000 | READ(iounit) nuvz, nwz, nz, nmixz, nlev_ec |
---|
1001 | READ(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst |
---|
1002 | |
---|
1003 | |
---|
1004 | |
---|
1005 | ! Get the varid , then read into scalar variable |
---|
1006 | ncret = nf90_inq_varid(ncid, 'nx', ncvarid) |
---|
1007 | ncret = nf90_get_var(ncid, ncvarid, nx) |
---|
1008 | |
---|
1009 | ncret = nf90_inq_varid(ncid, 'ny', ncvarid) |
---|
1010 | ncret = nf90_get_var(ncid, ncvarid, ny) |
---|
1011 | |
---|
1012 | ncret = nf90_inq_varid(ncid, 'nxmin1', ncvarid) |
---|
1013 | ncret = nf90_get_var(ncid, ncvarid, nxmin1) |
---|
1014 | |
---|
1015 | ncret = nf90_inq_varid(ncid, 'nymin1', ncvarid) |
---|
1016 | ncret = nf90_get_var(ncid, ncvarid, nymin1) |
---|
1017 | |
---|
1018 | ncret = nf90_inq_varid(ncid, 'nxfield', ncvarid) |
---|
1019 | ncret = nf90_get_var(ncid, ncvarid, nxfield) |
---|
1020 | |
---|
1021 | ncret = nf90_inq_varid(ncid, 'nuvz', ncvarid) |
---|
1022 | ncret = nf90_get_var(ncid, ncvarid, nuvz) |
---|
1023 | |
---|
1024 | ncret = nf90_inq_varid(ncid, 'nwz', ncvarid) |
---|
1025 | ncret = nf90_get_var(ncid, ncvarid, nwz) |
---|
1026 | |
---|
1027 | ncret = nf90_inq_varid(ncid, 'nz', ncvarid) |
---|
1028 | ncret = nf90_get_var(ncid, ncvarid, nz) |
---|
1029 | |
---|
1030 | ncret = nf90_inq_varid(ncid, 'nmixz', ncvarid) |
---|
1031 | ncret = nf90_get_var(ncid, ncvarid, nmixz) |
---|
1032 | |
---|
1033 | ncret = nf90_inq_varid(ncid, 'nlev_ec', ncvarid) |
---|
1034 | ncret = nf90_get_var(ncid, ncvarid, nlev_ec) |
---|
1035 | |
---|
1036 | ncret = nf90_inq_varid(ncid, 'dx', ncvarid) |
---|
1037 | ncret = nf90_get_var(ncid, ncvarid, dx) |
---|
1038 | |
---|
1039 | ncret = nf90_inq_varid(ncid, 'dy', ncvarid) |
---|
1040 | ncret = nf90_get_var(ncid, ncvarid, dy) |
---|
1041 | |
---|
1042 | ncret = nf90_inq_varid(ncid, 'xlon0', ncvarid) |
---|
1043 | ncret = nf90_get_var(ncid, ncvarid, xlon0) |
---|
1044 | |
---|
1045 | ncret = nf90_inq_varid(ncid, 'ylat0', ncvarid) |
---|
1046 | ncret = nf90_get_var(ncid, ncvarid, ylat0) |
---|
1047 | |
---|
1048 | ncret = nf90_inq_varid(ncid, 'dxconst', ncvarid) |
---|
1049 | ncret = nf90_get_var(ncid, ncvarid, dxconst) |
---|
1050 | |
---|
1051 | ncret = nf90_inq_varid(ncid, 'dyconst', ncvarid) |
---|
1052 | ncret = nf90_get_var(ncid, ncvarid, dyconst) |
---|
1053 | |
---|
1054 | |
---|
1055 | |
---|
1056 | |
---|
1057 | |
---|
1058 | |
---|
1059 | ! Fixed fields, static in time |
---|
1060 | READ(iounit) oro, excessoro, lsm, xlanduse, height |
---|
1061 | |
---|
1062 | ncret = nf90_inq_varid(ncid, 'oro', ncvarid) |
---|
1063 | ncret = nf90_get_var(ncid, ncvarid, oro(0:nxmax-1,0:nymax-1)) |
---|
1064 | |
---|
1065 | ncret = nf90_inq_varid(ncid, 'excessoro', ncvarid) |
---|
1066 | ncret = nf90_get_var(ncid, ncvarid, excessoro(0:nxmax-1,0:nymax-1)) |
---|
1067 | |
---|
1068 | ncret = nf90_inq_varid(ncid, 'lsm', ncvarid) |
---|
1069 | ncret = nf90_get_var(ncid, ncvarid, lsm(0:nxmax-1,0:nymax-1)) |
---|
1070 | |
---|
1071 | ncret = nf90_inq_varid(ncid, 'xlanduse', ncvarid) |
---|
1072 | ncret = nf90_get_var(ncid, ncvarid, xlanduse(0:nxmax-1,0:nymax-1, 1:numclass)) |
---|
1073 | |
---|
1074 | ncret = nf90_inq_varid(ncid, 'height', ncvarid) |
---|
1075 | ncret = nf90_get_var(ncid, ncvarid, height(1:nzmax)) |
---|
1076 | |
---|
1077 | |
---|
1078 | |
---|
1079 | |
---|
1080 | ! 3d fields |
---|
1081 | READ(iounit) uu(:,:,:,cm_index) |
---|
1082 | READ(iounit) vv(:,:,:,cm_index) |
---|
1083 | READ(iounit) uupol(:,:,:,cm_index) |
---|
1084 | READ(iounit) vvpol(:,:,:,cm_index) |
---|
1085 | READ(iounit) ww(:,:,:,cm_index) |
---|
1086 | READ(iounit) tt(:,:,:,cm_index) |
---|
1087 | READ(iounit) qv(:,:,:,cm_index) |
---|
1088 | READ(iounit) pv(:,:,:,cm_index) |
---|
1089 | READ(iounit) rho(:,:,:,cm_index) |
---|
1090 | READ(iounit) drhodz(:,:,:,cm_index) |
---|
1091 | READ(iounit) tth(:,:,:,cm_index) |
---|
1092 | READ(iounit) qvh(:,:,:,cm_index) |
---|
1093 | READ(iounit) pplev(:,:,:,cm_index) |
---|
1094 | READ(iounit) clouds(:,:,:,cm_index) |
---|
1095 | READ(iounit) cloudsh(:,:,cm_index) |
---|
1096 | |
---|
1097 | |
---|
1098 | |
---|
1099 | |
---|
1100 | ! Get the varid and read the variable into the array |
---|
1101 | ncret = nf90_inq_varid(ncid, 'uu', ncvarid) |
---|
1102 | ncret = nf90_get_var(ncid, ncvarid, uu(0:nxmax-1,0:nymax-1,1:nzmax,cm_index)) |
---|
1103 | |
---|
1104 | ncret = nf90_inq_varid(ncid, 'vv', ncvarid) |
---|
1105 | ncret = nf90_get_var(ncid, ncvarid, vv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index)) |
---|
1106 | |
---|
1107 | ncret = nf90_inq_varid(ncid, 'uupol', ncvarid) |
---|
1108 | ncret = nf90_get_var(ncid, ncvarid, uupol(0:nxmax-1,0:nymax-1,1:nzmax,cm_index)) |
---|
1109 | |
---|
1110 | ncret = nf90_inq_varid(ncid, 'vvpol', ncvarid) |
---|
1111 | ncret = nf90_get_var(ncid, ncvarid, vvpol(0:nxmax-1,0:nymax-1,1:nzmax,cm_index)) |
---|
1112 | |
---|
1113 | ncret = nf90_inq_varid(ncid, 'ww', ncvarid) |
---|
1114 | ncret = nf90_get_var(ncid, ncvarid, ww(0:nxmax-1,0:nymax-1,1:nzmax,cm_index)) |
---|
1115 | |
---|
1116 | ncret = nf90_inq_varid(ncid, 'tt', ncvarid) |
---|
1117 | ncret = nf90_get_var(ncid, ncvarid, tt(0:nxmax-1,0:nymax-1,1:nzmax,cm_index)) |
---|
1118 | |
---|
1119 | ncret = nf90_inq_varid(ncid, 'qv', ncvarid) |
---|
1120 | ncret = nf90_get_var(ncid, ncvarid, qv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index)) |
---|
1121 | |
---|
1122 | ncret = nf90_inq_varid(ncid, 'pv', ncvarid) |
---|
1123 | ncret = nf90_get_var(ncid, ncvarid, pv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index)) |
---|
1124 | |
---|
1125 | ncret = nf90_inq_varid(ncid, 'rho', ncvarid) |
---|
1126 | ncret = nf90_get_var(ncid, ncvarid, rho(0:nxmax-1,0:nymax-1,1:nzmax,cm_index)) |
---|
1127 | |
---|
1128 | ncret = nf90_inq_varid(ncid, 'drhodz', ncvarid) |
---|
1129 | ncret = nf90_get_var(ncid, ncvarid, drhodz(0:nxmax-1,0:nymax-1,1:nzmax,cm_index)) |
---|
1130 | |
---|
1131 | ncret = nf90_inq_varid(ncid, 'clouds', ncvarid) |
---|
1132 | ncret = nf90_get_var(ncid, ncvarid, clouds(0:nxmax-1,0:nzmax-1,1:nzmax,cm_index)) |
---|
1133 | |
---|
1134 | ncret = nf90_inq_varid(ncid, 'tth', ncvarid) |
---|
1135 | ncret = nf90_get_var(ncid, ncvarid, tth(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index)) |
---|
1136 | |
---|
1137 | ncret = nf90_inq_varid(ncid, 'qvh', ncvarid) |
---|
1138 | ncret = nf90_get_var(ncid, ncvarid, qvh(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index)) |
---|
1139 | |
---|
1140 | ncret = nf90_inq_varid(ncid, 'pplev', ncvarid) |
---|
1141 | ncret = nf90_get_var(ncid, ncvarid, pplev(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index)) |
---|
1142 | |
---|
1143 | ncret = nf90_inq_varid(ncid, 'cloudsh', ncvarid) |
---|
1144 | ncret = nf90_get_var(ncid, ncvarid, cloudsh(0:nxmax-1,0:nymax-1,cm_index)) |
---|
1145 | |
---|
1146 | |
---|
1147 | |
---|
1148 | |
---|
1149 | PRINT *, 'SUM(tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', & |
---|
1150 | & SUM(tt(0:nxmax-1,0:nymax-1,1:nzmax, cm_index)) |
---|
1151 | |
---|
1152 | |
---|
1153 | PRINT *, 'SUM(clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', & |
---|
1154 | & SUM(clouds(0:nxmax-1,0:nymax-1,1:nzmax, cm_index)) |
---|
1155 | |
---|
1156 | |
---|
1157 | |
---|
1158 | |
---|
1159 | ! 2d fields |
---|
1160 | READ(iounit) ps(:,:,:,cm_index) |
---|
1161 | READ(iounit) sd(:,:,:,cm_index) |
---|
1162 | READ(iounit) msl(:,:,:,cm_index) |
---|
1163 | READ(iounit) tcc(:,:,:,cm_index) |
---|
1164 | READ(iounit) u10(:,:,:,cm_index) |
---|
1165 | READ(iounit) v10(:,:,:,cm_index) |
---|
1166 | READ(iounit) tt2(:,:,:,cm_index) |
---|
1167 | READ(iounit) td2(:,:,:,cm_index) |
---|
1168 | READ(iounit) lsprec(:,:,:,cm_index) |
---|
1169 | READ(iounit) convprec(:,:,:,cm_index) |
---|
1170 | READ(iounit) sshf(:,:,:,cm_index) |
---|
1171 | READ(iounit) ssr(:,:,:,cm_index) |
---|
1172 | READ(iounit) surfstr(:,:,:,cm_index) |
---|
1173 | READ(iounit) ustar(:,:,:,cm_index) |
---|
1174 | READ(iounit) wstar(:,:,:,cm_index) |
---|
1175 | READ(iounit) hmix(:,:,:,cm_index) |
---|
1176 | READ(iounit) tropopause(:,:,:,cm_index) |
---|
1177 | READ(iounit) oli(:,:,:,cm_index) |
---|
1178 | READ(iounit) diffk(:,:,:,cm_index) |
---|
1179 | READ(iounit) vdep(:,:,:,cm_index) |
---|
1180 | |
---|
1181 | ! Get the varid and read the variable into the array |
---|
1182 | ncret = nf90_inq_varid(ncid, 'ps', ncvarid) |
---|
1183 | ncret = nf90_get_var(ncid, ncvarid, ps(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1184 | |
---|
1185 | ncret = nf90_inq_varid(ncid, 'sd', ncvarid) |
---|
1186 | ncret = nf90_get_var(ncid, ncvarid, sd(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1187 | |
---|
1188 | ncret = nf90_inq_varid(ncid, 'msl', ncvarid) |
---|
1189 | ncret = nf90_get_var(ncid, ncvarid, msl(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1190 | |
---|
1191 | ncret = nf90_inq_varid(ncid, 'tcc', ncvarid) |
---|
1192 | ncret = nf90_get_var(ncid, ncvarid, tcc(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1193 | |
---|
1194 | ncret = nf90_inq_varid(ncid, 'u10', ncvarid) |
---|
1195 | ncret = nf90_get_var(ncid, ncvarid, u10(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1196 | |
---|
1197 | ncret = nf90_inq_varid(ncid, 'v10', ncvarid) |
---|
1198 | ncret = nf90_get_var(ncid, ncvarid, v10(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1199 | |
---|
1200 | ncret = nf90_inq_varid(ncid, 'tt2', ncvarid) |
---|
1201 | ncret = nf90_get_var(ncid, ncvarid, tt2(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1202 | |
---|
1203 | ncret = nf90_inq_varid(ncid, 'td2', ncvarid) |
---|
1204 | ncret = nf90_get_var(ncid, ncvarid, td2(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1205 | |
---|
1206 | ncret = nf90_inq_varid(ncid, 'lsprec', ncvarid) |
---|
1207 | ncret = nf90_get_var(ncid, ncvarid, lsprec(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1208 | |
---|
1209 | ncret = nf90_inq_varid(ncid, 'convprec', ncvarid) |
---|
1210 | ncret = nf90_get_var(ncid, ncvarid, convprec(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1211 | |
---|
1212 | ncret = nf90_inq_varid(ncid, 'sshf', ncvarid) |
---|
1213 | ncret = nf90_get_var(ncid, ncvarid, sshf(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1214 | |
---|
1215 | ncret = nf90_inq_varid(ncid, 'ssr', ncvarid) |
---|
1216 | ncret = nf90_get_var(ncid, ncvarid, ssr(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1217 | |
---|
1218 | ncret = nf90_inq_varid(ncid, 'surfstr', ncvarid) |
---|
1219 | ncret = nf90_get_var(ncid, ncvarid, surfstr(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1220 | |
---|
1221 | ncret = nf90_inq_varid(ncid, 'ustar', ncvarid) |
---|
1222 | ncret = nf90_get_var(ncid, ncvarid, ustar(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1223 | |
---|
1224 | ncret = nf90_inq_varid(ncid, 'wstar', ncvarid) |
---|
1225 | ncret = nf90_get_var(ncid, ncvarid, wstar(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1226 | |
---|
1227 | ncret = nf90_inq_varid(ncid, 'hmix', ncvarid) |
---|
1228 | ncret = nf90_get_var(ncid, ncvarid, hmix(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1229 | |
---|
1230 | ncret = nf90_inq_varid(ncid, 'tropopause', ncvarid) |
---|
1231 | ncret = nf90_get_var(ncid, ncvarid, tropopause(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1232 | |
---|
1233 | ncret = nf90_inq_varid(ncid, 'oli', ncvarid) |
---|
1234 | ncret = nf90_get_var(ncid, ncvarid, oli(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1235 | |
---|
1236 | ncret = nf90_inq_varid(ncid, 'diffk', ncvarid) |
---|
1237 | ncret = nf90_get_var(ncid, ncvarid, diffk(0:nxmax-1,0:nymax-1,1,cm_index)) |
---|
1238 | |
---|
1239 | |
---|
1240 | |
---|
1241 | |
---|
1242 | PRINT *, 'SUM(ps(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', & |
---|
1243 | & SUM(ps(0:nxmax-1,0:nymax-1,1, cm_index)) |
---|
1244 | |
---|
1245 | PRINT *, 'SUM(wstar(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', & |
---|
1246 | & SUM(wstar(0:nxmax-1,0:nymax-1,1, cm_index)) |
---|
1247 | |
---|
1248 | |
---|
1249 | |
---|
1250 | |
---|
1251 | |
---|
1252 | |
---|
1253 | ! 1d fields |
---|
1254 | READ(iounit) z0(:) |
---|
1255 | READ(iounit) akm(:) |
---|
1256 | READ(iounit) bkm(:) |
---|
1257 | READ(iounit) akz(:) |
---|
1258 | READ(iounit) bkz(:) |
---|
1259 | READ(iounit) aknew(:) |
---|
1260 | READ(iounit) bknew(:) |
---|
1261 | |
---|
1262 | |
---|
1263 | ! Nested, scalar values (for each nest) |
---|
1264 | READ(iounit) nxn(:) |
---|
1265 | READ(iounit) nyn(:) |
---|
1266 | READ(iounit) dxn(:) |
---|
1267 | READ(iounit) dyn(:) |
---|
1268 | READ(iounit) xlon0n(:) |
---|
1269 | READ(iounit) ylat0n(:) |
---|
1270 | |
---|
1271 | |
---|
1272 | ! Nested fields, static over time |
---|
1273 | READ(iounit) oron, excessoron, lsmn, xlandusen |
---|
1274 | |
---|
1275 | ! 3d nested fields |
---|
1276 | READ(iounit) uun(:,:,:,cm_index,:) |
---|
1277 | READ(iounit) vvn(:,:,:,cm_index,:) |
---|
1278 | READ(iounit) wwn(:,:,:,cm_index,:) |
---|
1279 | READ(iounit) ttn(:,:,:,cm_index,:) |
---|
1280 | READ(iounit) qvn(:,:,:,cm_index,:) |
---|
1281 | READ(iounit) pvn(:,:,:,cm_index,:) |
---|
1282 | READ(iounit) cloudsn(:,:,:,cm_index,:) |
---|
1283 | READ(iounit) cloudsnh(:,:,cm_index,:) |
---|
1284 | READ(iounit) rhon(:,:,:,cm_index,:) |
---|
1285 | READ(iounit) drhodzn(:,:,:,cm_index,:) |
---|
1286 | READ(iounit) tthn(:,:,:,cm_index,:) |
---|
1287 | READ(iounit) qvhn(:,:,:,cm_index,:) |
---|
1288 | |
---|
1289 | ! 2d nested fields |
---|
1290 | READ(iounit) psn(:,:,:,cm_index,:) |
---|
1291 | READ(iounit) sdn(:,:,:,cm_index,:) |
---|
1292 | READ(iounit) msln(:,:,:,cm_index,:) |
---|
1293 | READ(iounit) tccn(:,:,:,cm_index,:) |
---|
1294 | READ(iounit) u10n(:,:,:,cm_index,:) |
---|
1295 | READ(iounit) v10n(:,:,:,cm_index,:) |
---|
1296 | READ(iounit) tt2n(:,:,:,cm_index,:) |
---|
1297 | READ(iounit) td2n(:,:,:,cm_index,:) |
---|
1298 | READ(iounit) lsprecn(:,:,:,cm_index,:) |
---|
1299 | READ(iounit) convprecn(:,:,:,cm_index,:) |
---|
1300 | READ(iounit) sshfn(:,:,:,cm_index,:) |
---|
1301 | READ(iounit) ssrn(:,:,:,cm_index,:) |
---|
1302 | READ(iounit) surfstrn(:,:,:,cm_index,:) |
---|
1303 | READ(iounit) ustarn(:,:,:,cm_index,:) |
---|
1304 | READ(iounit) wstarn(:,:,:,cm_index,:) |
---|
1305 | READ(iounit) hmixn(:,:,:,cm_index,:) |
---|
1306 | READ(iounit) tropopausen(:,:,:,cm_index,:) |
---|
1307 | READ(iounit) olin(:,:,:,cm_index,:) |
---|
1308 | READ(iounit) diffkn(:,:,:,cm_index,:) |
---|
1309 | READ(iounit) vdepn(:,:,:,cm_index,:) |
---|
1310 | |
---|
1311 | ! Auxiliary variables for nests |
---|
1312 | READ(iounit) xresoln(:) |
---|
1313 | READ(iounit) yresoln(:) |
---|
1314 | READ(iounit) xln(:) |
---|
1315 | READ(iounit) yln(:) |
---|
1316 | READ(iounit) xrn(:) |
---|
1317 | READ(iounit) yrn(:) |
---|
1318 | |
---|
1319 | ! Variables for polar stereographic projection |
---|
1320 | READ(iounit) xglobal, sglobal, nglobal |
---|
1321 | READ(iounit) switchnorthg, switchsouthg |
---|
1322 | READ(iounit) southpolemap(:) |
---|
1323 | READ(iounit) northpolemap(:) |
---|
1324 | |
---|
1325 | ! Variables declared in conv_mod (convection) |
---|
1326 | READ(iounit) pconv(:) |
---|
1327 | READ(iounit) phconv(:) |
---|
1328 | READ(iounit) dpr(:) |
---|
1329 | READ(iounit) pconv_hpa(:) |
---|
1330 | READ(iounit) phconv_hpa(:) |
---|
1331 | READ(iounit) ft(:) |
---|
1332 | READ(iounit) fq(:) |
---|
1333 | READ(iounit) fmass(:,:) |
---|
1334 | READ(iounit) sub(:) |
---|
1335 | READ(iounit) fmassfrac(:,:) |
---|
1336 | READ(iounit) cbaseflux(:,:) |
---|
1337 | READ(iounit) cbasefluxn(:,:,:) |
---|
1338 | READ(iounit) tconv(:) |
---|
1339 | READ(iounit) qconv(:) |
---|
1340 | READ(iounit) qsconv(:) |
---|
1341 | READ(iounit) psconv, tt2conv, td2conv |
---|
1342 | READ(iounit) nconvlev, nconvtop |
---|
1343 | |
---|
1344 | ELSE |
---|
1345 | STOP 'fpio(): Illegal operation' |
---|
1346 | |
---|
1347 | ENDIF |
---|
1348 | END SUBROUTINE fpio |
---|
1349 | |
---|
1350 | SUBROUTINE fpmetbinary_filetext(filename, cm_index) |
---|
1351 | |
---|
1352 | ! This is a utility subroutine meant to be used for testing purposes. |
---|
1353 | ! It facilitates the text output of variables read in from the |
---|
1354 | ! specified .fp file. This routine will easily cause the program |
---|
1355 | ! to crash due memory allocation issues, particularly when you are |
---|
1356 | ! trying to text print 3d arrays in a single formatted statetment. |
---|
1357 | |
---|
1358 | CHARACTER(LEN=*), INTENT(IN) :: filename |
---|
1359 | INTEGER, INTENT(IN) :: cm_index ! Index of last dimension in |
---|
1360 | ! most com_mod variables. |
---|
1361 | ! Should be 1 or 2 |
---|
1362 | |
---|
1363 | !OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', status='REPLACE', & |
---|
1364 | ! form="formatted", access="stream") |
---|
1365 | OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', & |
---|
1366 | form="formatted", access="APPEND") |
---|
1367 | |
---|
1368 | WRITE(IOUNIT_TEXTOUT, *) 'oro: ', oro |
---|
1369 | WRITE(IOUNIT_TEXTOUT, *) 'excessoro: ', excessoro |
---|
1370 | WRITE(IOUNIT_TEXTOUT, *) 'lsm: ', lsm |
---|
1371 | WRITE(IOUNIT_TEXTOUT, *) 'xlanduse: ', xlanduse |
---|
1372 | WRITE(IOUNIT_TEXTOUT, *) 'height: ', height |
---|
1373 | |
---|
1374 | WRITE(IOUNIT_TEXTOUT, *) 'uu: ', uu(:,:,:,cm_index) |
---|
1375 | WRITE(IOUNIT_TEXTOUT, *) 'vv: ', vv(:,:,:,cm_index) |
---|
1376 | WRITE(IOUNIT_TEXTOUT, *) 'uupol: ', uupol(:,:,:,cm_index) |
---|
1377 | WRITE(IOUNIT_TEXTOUT, *) 'vvpol: ', vvpol(:,:,:,cm_index) |
---|
1378 | WRITE(IOUNIT_TEXTOUT, *) 'ww: ', ww(:,:,:,cm_index) |
---|
1379 | WRITE(IOUNIT_TEXTOUT, *) 'tt: ', tt(:,:,:,cm_index) |
---|
1380 | WRITE(IOUNIT_TEXTOUT, *) 'qv: ', qv(:,:,:,cm_index) |
---|
1381 | WRITE(IOUNIT_TEXTOUT, *) 'pv: ', pv(:,:,:,cm_index) |
---|
1382 | WRITE(IOUNIT_TEXTOUT, *) 'rho: ', rho(:,:,:,cm_index) |
---|
1383 | WRITE(IOUNIT_TEXTOUT, *) 'drhodz: ', drhodz(:,:,:,cm_index) |
---|
1384 | WRITE(IOUNIT_TEXTOUT, *) 'tth: ', tth(:,:,:,cm_index) |
---|
1385 | WRITE(IOUNIT_TEXTOUT, *) 'qvh: ', qvh(:,:,:,cm_index) |
---|
1386 | WRITE(IOUNIT_TEXTOUT, *) 'pplev: ', pplev(:,:,:,cm_index) |
---|
1387 | WRITE(IOUNIT_TEXTOUT, *) 'clouds: ', clouds(:,:,:,cm_index) |
---|
1388 | WRITE(IOUNIT_TEXTOUT, *) 'cloudsh: ', cloudsh(:,:,cm_index) |
---|
1389 | |
---|
1390 | |
---|
1391 | |
---|
1392 | |
---|
1393 | CLOSE(IOUNIT_TEXTOUT) |
---|
1394 | END SUBROUTINE fpmetbinary_filetext |
---|
1395 | |
---|
1396 | |
---|
1397 | END MODULE fpmetbinary_mod |
---|