source: flexpart.git/flexpart_code/fpmetbinary_mod.F90 @ 8651603

fp9.3.1-20161214-nc4
Last change on this file since 8651603 was 8651603, checked in by Don Morton <Don.Morton@…>, 7 years ago

Incremental save for branch fp9.3.1-20161214-nc4

  • Property mode set to 100644
File size: 75.7 KB
Line 
1MODULE 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, maxspec, &
47&                       maxnests, nxmaxn, nymaxn
48
49    USE netcdf
50
51    IMPLICIT NONE
52
53    ! Users may want to change these IO Unit values if they conflict with other parts
54    ! of code
55    INTEGER, PARAMETER :: IOUNIT_DUMP = 33, IOUNIT_LOAD = 34, &
56                          IOUNIT_TEXTOUT = 35
57
58    ! When a change is made to the format of the preprocessed file, such that
59    ! this routine will not be able to read a previous version, this version
60    ! string should be modified
61    CHARACTER(LEN=12), PARAMETER :: PREPROC_FORMAT_VERSION_STR = 'FP_p-9.3.1'//char(0)
62
63    PRIVATE IOUNIT_DUMP, IOUNIT_LOAD, IOUNIT_TEXTOUT, fpio,    &
64&           PREPROC_FORMAT_VERSION_STR
65
66
67CONTAINS
68
69  !*****************************************************************************
70  !                                                                            *
71  !    Subroutines fpmetbinary_dump() and fpmetbinary_load() provide the       *
72  !    public interface to                                                     *
73  !    this module functionality.  I created the PRIVATE fpio() because I      *
74  !    wanted all interactions with variables to be in one place.  The read    *
75  !    and write operations need to be done in exactly the same sequence, so   *
76  !    I felt like keeping them in the same routine would at least allow for   *
77  !    coders to more easily compare the two sequences than if they were       *
78  !    separate.                                                               *
79  !                                                                            *
80  !    As mentioned above, the dumps and loads will, for most variables,       *
81  !    need to refer to one of two index values for the last dimension of      *
82  !    the array.                                                              *
83  !                                                                            *
84  !*****************************************************************************
85
86
87    SUBROUTINE fpmetbinary_dump(filename, cm_index)
88        CHARACTER(LEN=*), INTENT(IN) :: filename  ! Full path for file
89        INTEGER, INTENT(IN) :: cm_index           ! Index of last dimension in
90                                                  ! most com_mod variables.
91                                                  ! Should be 1 or 2
92
93        INTEGER millisecs_start, millisecs_stop, count_rate, count_max
94
95        INTEGER :: ncretval, ncid          ! NetCDF func return value, file id
96
97        CALL SYSTEM_CLOCK(millisecs_start, count_rate, count_max)
98
99        ! Create and open NC4 file for writing
100        PRINT *, 'Opening NC4 file...'
101        ncretval = nf90_create(filename // ".nc4", &
102&                              OR(NF90_CLOBBER, NF90_HDF5), &
103&                              ncid)
104
105        OPEN(IOUNIT_DUMP, file=filename, action='WRITE', status='REPLACE', form="unformatted", access="stream")
106
107
108
109
110
111
112        CALL fpio(IOUNIT_DUMP, ncid, 'DUMP', cm_index)
113        CLOSE(IOUNIT_DUMP)
114
115        PRINT *, 'Closing NC4 file...'
116        ncretval = nf90_close(ncid)
117
118        CALL SYSTEM_CLOCK(millisecs_stop, count_rate, count_max)
119
120        !PRINT *, 'Dump walltime secs: ', (millisecs_stop-millisecs_start)/1000.0
121    END SUBROUTINE fpmetbinary_dump
122
123    SUBROUTINE fpmetbinary_load(filename, cm_index)
124        CHARACTER(LEN=*), INTENT(IN) :: filename  ! Full path for file
125        INTEGER, INTENT(IN) :: cm_index           ! Index of last dimension in
126                                                  ! most com_mod variables.
127                                                  ! Should be 1 or 2
128
129        INTEGER :: ncretval, ncid          ! NetCDF func return value, file id
130
131        INTEGER millisecs_start, millisecs_stop, count_rate, count_max
132
133        CALL SYSTEM_CLOCK(millisecs_start, count_rate, count_max)
134
135        print *, "Opening nc file for reading"
136        ncretval = nf90_open(filename // ".nc4", NF90_NOWRITE, ncid)
137
138
139
140        OPEN(IOUNIT_LOAD, file=filename, action='READ', status='OLD', form="unformatted", access="stream")
141        CALL fpio(IOUNIT_LOAD, ncid, 'LOAD', cm_index)
142        CLOSE(IOUNIT_LOAD)
143        CALL SYSTEM_CLOCK(millisecs_stop, count_rate, count_max)
144        !PRINT *, 'Load walltime secs: ', (millisecs_stop-millisecs_start)/1000.0
145    END SUBROUTINE fpmetbinary_load
146
147    SUBROUTINE fpmetbinary_zero(cm_index)
148        INTEGER, INTENT(IN) :: cm_index           ! Index of last dimension in
149                                                  ! most com_mod variables.
150                                                  ! Should be 1 or 2
151
152
153        ! Zeroes out, in local datastructures, the values dumped/loaded
154        ! This was written primarily as a testing mechanism.
155        ! DJM -- 17 February 2017 -- I don't think this routine has been used
156        ! for anything in recent past.  Might want to consider 86'ing it.
157        ! The lines here correspond to READ and WRITE in the dump/load routines
158
159        ! Scalar values
160        nx=0.0; ny=0.0; nxmin1=0.0; nymin1=0.0; nxfield=0.0
161        nuvz=0.0; nwz=0.0; nz=0.0; nmixz=0.0; nlev_ec=0.0
162        dx=0.0; dy=0.0; xlon0=0.0; ylat0=0.0; dxconst=0.0; dyconst=0.0
163
164        ! Fixed fields, static in time
165        oro=0.0; excessoro=0.0; lsm=0.0; xlanduse=0.0; height=0.0
166
167        ! 3d fields
168        uu(:,:,:,cm_index) = 0.0
169        vv(:,:,:,cm_index) = 0.0
170        uupol(:,:,:,cm_index) = 0.0
171        vvpol(:,:,:,cm_index) = 0.0
172        ww(:,:,:,cm_index) = 0.0
173        tt(:,:,:,cm_index) = 0.0
174        qv(:,:,:,cm_index) = 0.0
175        pv(:,:,:,cm_index) = 0.0
176        rho(:,:,:,cm_index) = 0.0
177        drhodz(:,:,:,cm_index) = 0.0
178        tth(:,:,:,cm_index) = 0.0
179        qvh(:,:,:,cm_index) = 0.0
180        pplev(:,:,:,cm_index) = 0.0
181        clouds(:,:,:,cm_index) = 0.0
182        cloudsh(:,:,cm_index) = 0.0
183
184        ! 2d fields
185        ps(:,:,:,cm_index) = 0.0
186        sd(:,:,:,cm_index) = 0.0
187        msl(:,:,:,cm_index) = 0.0
188        tcc(:,:,:,cm_index) = 0.0
189        u10(:,:,:,cm_index) = 0.0
190        v10(:,:,:,cm_index) = 0.0
191        tt2(:,:,:,cm_index) = 0.0
192        td2(:,:,:,cm_index) = 0.0
193        lsprec(:,:,:,cm_index) = 0.0
194        convprec(:,:,:,cm_index) = 0.0
195        sshf(:,:,:,cm_index) = 0.0
196        ssr(:,:,:,cm_index) = 0.0
197        surfstr(:,:,:,cm_index) = 0.0
198        ustar(:,:,:,cm_index) = 0.0
199        wstar(:,:,:,cm_index) = 0.0
200        hmix(:,:,:,cm_index) = 0.0
201        tropopause(:,:,:,cm_index) = 0.0
202        oli(:,:,:,cm_index) = 0.0
203        diffk(:,:,:,cm_index) = 0.0
204        vdep(:,:,:,cm_index) = 0.0
205
206        ! 1d fields
207        z0(:) = 0.0
208        akm(:) = 0.0
209        bkm(:) = 0.0
210        akz(:) = 0.0
211        bkz(:) = 0.0
212        aknew(:) = 0.0
213        bknew(:) = 0.0
214
215        ! Nested, scalar values (for each nest)
216        nxn(:) = 0.0
217        nyn(:) = 0.0
218        dxn(:) = 0.0
219        dyn(:) = 0.0
220        xlon0n(:) = 0.0
221        ylat0n(:) = 0.0
222
223        ! Nested fields, static in time
224        oron=0.0; excessoron=0.0; lsmn=0.0; xlandusen=0.0
225
226        ! 3d nested fields
227        uun(:,:,:,cm_index,:) = 0.0
228        wwn(:,:,:,cm_index,:) = 0.0
229        ttn(:,:,:,cm_index,:) = 0.0
230        qvn(:,:,:,cm_index,:) = 0.0
231        pvn(:,:,:,cm_index,:) = 0.0
232        cloudsn(:,:,:,cm_index,:) = 0.0
233        cloudsnh(:,:,cm_index,:) = 0.0
234        rhon(:,:,:,cm_index,:) = 0.0
235        drhodzn(:,:,:,cm_index,:) = 0.0
236        tthn(:,:,:,cm_index,:) = 0.0
237        qvhn(:,:,:,cm_index,:) = 0.0
238
239        ! 2d nested fields
240        psn(:,:,:,cm_index,:) = 0.0
241        sdn(:,:,:,cm_index,:) = 0.0
242        msln(:,:,:,cm_index,:) = 0.0
243        tccn(:,:,:,cm_index,:) = 0.0
244        u10n(:,:,:,cm_index,:) = 0.0
245        v10n(:,:,:,cm_index,:) = 0.0
246        tt2n(:,:,:,cm_index,:) = 0.0
247        td2n(:,:,:,cm_index,:) = 0.0
248        lsprecn(:,:,:,cm_index,:) = 0.0
249        convprecn(:,:,:,cm_index,:) = 0.0
250        sshfn(:,:,:,cm_index,:) = 0.0
251        ssrn(:,:,:,cm_index,:) = 0.0
252        surfstrn(:,:,:,cm_index,:) = 0.0
253        ustarn(:,:,:,cm_index,:) = 0.0
254        wstarn(:,:,:,cm_index,:) = 0.0
255        hmixn(:,:,:,cm_index,:) = 0.0
256        tropopausen(:,:,:,cm_index,:) = 0.0
257        olin(:,:,:,cm_index,:) = 0.0
258        diffkn(:,:,:,cm_index,:) = 0.0
259        vdepn(:,:,:,cm_index,:) = 0.0
260
261        ! Auxiliary variables for nests
262        xresoln(:) = 0.0
263        yresoln(:) = 0.0
264        xln(:) = 0.0
265        yln(:) = 0.0
266        xrn(:) = 0.0
267        yrn(:) = 0.0
268
269        ! Variables for polar stereographic projection
270        xglobal=.FALSE.; sglobal=.FALSE.; nglobal=.FALSE.
271        switchnorthg=0.0; switchsouthg=0.0
272        southpolemap(:) = 0.0
273        northpolemap(:) = 0.0
274
275        ! Variables declared in conv_mod (convection)
276        pconv(:) = 0.0
277        phconv(:) = 0.0
278        dpr(:) = 0.0
279        pconv_hpa(:) = 0.0
280        phconv_hpa(:) = 0.0
281        ft(:) = 0.0
282        fq(:) = 0.0
283        fmass(:,:) = 0.0
284        sub(:) = 0.0
285        fmassfrac(:,:) = 0.0
286        cbaseflux(:,:) = 0.0
287        cbasefluxn(:,:,:) = 0.0
288        tconv(:) = 0.0
289        qconv(:) = 0.0
290        qsconv(:) = 0.0
291        psconv=0.0; tt2conv=0.0; td2conv=0.0
292        nconvlev=0.0; nconvtop=0.0
293
294    END SUBROUTINE fpmetbinary_zero
295
296    SUBROUTINE fpio(iounit, ncid, op, cm_index)
297        IMPLICIT NONE
298        INTEGER, INTENT(IN) :: ncid               ! NetCDF file id
299        INTEGER, INTENT(IN) :: iounit
300        CHARACTER(LEN=4), INTENT(IN) :: op        ! Operation - DUMP or LOAD
301        INTEGER, INTENT(IN) :: cm_index           ! Index of last dimension in
302                                                  ! most com_mod variables.
303                                                  ! Should be 1 or 2
304
305
306        INTEGER :: ncret          ! Return value from NetCDF calls
307        INTEGER :: ncvarid          ! NetCDF variable ID
308
309        INTEGER :: nxmax_dimid, nymax_dimid, nzmax_dimid, nuvzmax_dimid, nwzmax_dimid, &
310&                  maxspec_dimid, numclass_dimid, maxnests_dimid, nxmaxn_dimid, nymaxn_dimid
311
312        INTEGER, DIMENSION(1) :: dim1dids    ! Dimension IDs for 1D arrays
313        INTEGER, DIMENSION(2) :: dim2dids    ! Dimension IDs for 2D arrays
314        INTEGER, DIMENSION(3) :: dim3dids    ! Dimension IDs for 3D arrays
315        INTEGER, DIMENSION(4) :: dim4dids    ! Dimension IDs for 4D arrays
316
317
318        ! These are used when loading in dimensions from NC file
319        CHARACTER(LEN=NF90_MAX_NAME) :: nxmax_dimname, nymax_dimname, nzmax_dimname, &
320&                                       nuvzmax_dimname, nwzmax_dimname,&
321&                                       maxspec_dimname, numclass_dimname,&
322&                                       maxnests_dimname, nxmaxn_dimname, nymaxn_dimname
323
324        ! These are temporary variables, used in the LOAD option, for
325        ! comparing against the current values in FLEXPART of nxmax, nymax, ...
326        INTEGER :: temp_nxmax, temp_nymax, temp_nzmax, &
327&                  temp_nuvzmax, temp_nwzmax, &
328&                  temp_maxspec, temp_numclass,&
329&                  temp_maxnests, temp_nxmaxn, temp_nymaxn
330
331        CHARACTER(LEN=12) :: temp_preproc_format_version_str
332
333        CHARACTER(LEN=128) :: errmesg
334
335        INTEGER, PARAMETER :: DEF_LEVEL = 3
336
337        if (op == 'DUMP') THEN
338
339
340            ! Write the preprocessing format version string
341            WRITE (iounit) PREPROC_FORMAT_VERSION_STR
342
343            ! Write the compiled max dimensions from par_mod - these are
344            ! not meant to be reassigned during a LOAD, but used as "header"
345            ! information to provide the structure of arrays
346            WRITE (iounit) nxmax, nymax, nzmax, nuvzmax, nwzmax
347
348            ncret = nf90_def_dim(ncid, 'nxmax', nxmax, nxmax_dimid)
349            ncret = nf90_def_dim(ncid, 'nymax', nymax, nymax_dimid)
350            ncret = nf90_def_dim(ncid, 'nzmax', nzmax, nzmax_dimid)
351            ncret = nf90_def_dim(ncid, 'nuvzmax', nuvzmax, nuvzmax_dimid)
352            ncret = nf90_def_dim(ncid, 'nwzmax', nwzmax, nwzmax_dimid)
353            ncret = nf90_def_dim(ncid, 'maxspec', maxspec, maxspec_dimid)
354            ncret = nf90_def_dim(ncid, 'numclass', numclass, numclass_dimid)
355
356            ! Scalar values
357            WRITE(iounit) nx, ny, nxmin1, nymin1, nxfield
358            WRITE(iounit) nuvz, nwz, nz, nmixz, nlev_ec
359            WRITE(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst
360
361            ncret = nf90_def_var(ncid, 'nx', NF90_INT, ncvarid)
362            ncret = nf90_put_var(ncid, ncvarid, nx)
363
364            ncret = nf90_def_var(ncid, 'ny', NF90_INT, ncvarid)
365            ncret = nf90_put_var(ncid, ncvarid, ny)
366
367            ncret = nf90_def_var(ncid, 'nxmin1', NF90_INT, ncvarid)
368            ncret = nf90_put_var(ncid, ncvarid, nxmin1)
369
370            ncret = nf90_def_var(ncid, 'nymin1', NF90_INT, ncvarid)
371            ncret = nf90_put_var(ncid, ncvarid, nymin1)
372
373            ncret = nf90_def_var(ncid, 'nxfield', NF90_INT, ncvarid)
374            ncret = nf90_put_var(ncid, ncvarid, nxfield)
375
376            ncret = nf90_def_var(ncid, 'nuvz', NF90_INT, ncvarid)
377            ncret = nf90_put_var(ncid, ncvarid, nuvz)
378
379            ncret = nf90_def_var(ncid, 'nwz', NF90_INT, ncvarid)
380            ncret = nf90_put_var(ncid, ncvarid, nwz)
381
382            ncret = nf90_def_var(ncid, 'nz', NF90_INT, ncvarid)
383            ncret = nf90_put_var(ncid, ncvarid, nz)
384
385            ncret = nf90_def_var(ncid, 'nmixz', NF90_INT, ncvarid)
386            ncret = nf90_put_var(ncid, ncvarid, nmixz)
387
388            ncret = nf90_def_var(ncid, 'nlev_ec', NF90_INT, ncvarid)
389            ncret = nf90_put_var(ncid, ncvarid, nlev_ec)
390
391            ncret = nf90_def_var(ncid, 'dx', NF90_FLOAT, ncvarid)
392            ncret = nf90_put_var(ncid, ncvarid, dx)
393
394            ncret = nf90_def_var(ncid, 'dy', NF90_FLOAT, ncvarid)
395            ncret = nf90_put_var(ncid, ncvarid, dy)
396
397            ncret = nf90_def_var(ncid, 'xlon0', NF90_FLOAT, ncvarid)
398            ncret = nf90_put_var(ncid, ncvarid, xlon0)
399
400            ncret = nf90_def_var(ncid, 'ylat0', NF90_FLOAT, ncvarid)
401            ncret = nf90_put_var(ncid, ncvarid, ylat0)
402
403            ncret = nf90_def_var(ncid, 'dxconst', NF90_FLOAT, ncvarid)
404            ncret = nf90_put_var(ncid, ncvarid, dxconst)
405
406            ncret = nf90_def_var(ncid, 'dyconst', NF90_FLOAT, ncvarid)
407            ncret = nf90_put_var(ncid, ncvarid, dyconst)
408
409
410
411            ! Fixed fields, static in time
412            WRITE(iounit) oro, excessoro, lsm, xlanduse, height
413
414            dim2dids = (/nxmax_dimid, nymax_dimid/)
415
416            ncret = nf90_def_var(ncid, 'oro', NF90_FLOAT, &
417&                                       dim2dids, ncvarid)
418            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
419&                                        shuffle=0,     &
420&                                        deflate=1,     &
421&                                        deflate_level=DEF_LEVEL)
422            ncret = nf90_put_var(ncid, ncvarid, &
423&                                oro(0:nxmax-1, 0:nymax-1))
424
425            ncret = nf90_def_var(ncid, 'excessoro', NF90_FLOAT, &
426&                                       dim2dids, ncvarid)
427            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
428&                                        shuffle=0,     &
429&                                        deflate=1,     &
430&                                        deflate_level=DEF_LEVEL)
431            ncret = nf90_put_var(ncid, ncvarid, &
432&                                excessoro(0:nxmax-1, 0:nymax-1))
433
434            ncret = nf90_def_var(ncid, 'lsm', NF90_FLOAT, &
435&                                       dim2dids, ncvarid)
436            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
437&                                        shuffle=0,     &
438&                                        deflate=1,     &
439&                                        deflate_level=DEF_LEVEL)
440            ncret = nf90_put_var(ncid, ncvarid, &
441&                                lsm(0:nxmax-1, 0:nymax-1))
442
443            dim3dids = (/nxmax_dimid, nymax_dimid, numclass_dimid/)
444            ! numclass comes from par_mod - number of land use classes
445            ncret = nf90_def_var(ncid, 'xlanduse', NF90_FLOAT, &
446&                                       dim3dids, 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&                                xlanduse(0:nxmax-1, 0:nymax-1, 1:numclass))
453
454            dim1dids = (/nzmax_dimid/)
455            ncret = nf90_def_var(ncid, 'height', NF90_FLOAT, &
456&                                       dim1dids, ncvarid)
457            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
458&                                        shuffle=0,     &
459&                                        deflate=1,     &
460&                                        deflate_level=DEF_LEVEL)
461            ncret = nf90_put_var(ncid, ncvarid, &
462&                                height(1:nzmax))
463
464
465
466
467            ! 3d fields
468            WRITE(iounit) uu(:,:,:,cm_index)
469            WRITE(iounit) vv(:,:,:,cm_index)
470            WRITE(iounit) uupol(:,:,:,cm_index)
471            WRITE(iounit) vvpol(:,:,:,cm_index)
472            WRITE(iounit) ww(:,:,:,cm_index)
473            WRITE(iounit) tt(:,:,:,cm_index)
474            WRITE(iounit) qv(:,:,:,cm_index)
475            WRITE(iounit) pv(:,:,:,cm_index)
476            WRITE(iounit) rho(:,:,:,cm_index)
477            WRITE(iounit) drhodz(:,:,:,cm_index)
478            WRITE(iounit) tth(:,:,:,cm_index)
479            WRITE(iounit) qvh(:,:,:,cm_index)
480            WRITE(iounit) pplev(:,:,:,cm_index)
481            WRITE(iounit) clouds(:,:,:,cm_index)
482            WRITE(iounit) cloudsh(:,:,cm_index)
483
484            dim3dids = (/nxmax_dimid, nymax_dimid, nzmax_dimid/)
485
486            ncret = nf90_def_var(ncid, 'uu', NF90_FLOAT, &
487&                                       dim3dids, ncvarid)
488            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
489&                                        shuffle=0,     &
490&                                        deflate=1,     &
491&                                        deflate_level=DEF_LEVEL)
492            ncret = nf90_put_var(ncid, ncvarid, &
493&                                uu(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
494
495            ncret = nf90_def_var(ncid, 'vv', NF90_FLOAT, &
496&                                       dim3dids, ncvarid)
497            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
498&                                        shuffle=0,     &
499&                                        deflate=1,     &
500&                                        deflate_level=DEF_LEVEL)
501            ncret = nf90_put_var(ncid, ncvarid, &
502&                                vv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
503
504            ncret = nf90_def_var(ncid, 'uupol', NF90_FLOAT, &
505&                                       dim3dids, ncvarid)
506            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
507&                                        shuffle=0,     &
508&                                        deflate=1,     &
509&                                        deflate_level=DEF_LEVEL)
510            ncret = nf90_put_var(ncid, ncvarid, &
511&                                uupol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
512
513            ncret = nf90_def_var(ncid, 'vvpol', NF90_FLOAT, &
514&                                       dim3dids, ncvarid)
515            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
516&                                        shuffle=0,     &
517&                                        deflate=1,     &
518&                                        deflate_level=DEF_LEVEL)
519            ncret = nf90_put_var(ncid, ncvarid, &
520&                                vvpol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
521
522            ncret = nf90_def_var(ncid, 'ww', NF90_FLOAT, &
523&                                       dim3dids, ncvarid)
524            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
525&                                        shuffle=0,     &
526&                                        deflate=1,     &
527&                                        deflate_level=DEF_LEVEL)
528            ncret = nf90_put_var(ncid, ncvarid, &
529&                                ww(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
530
531            ncret = nf90_def_var(ncid, 'tt', NF90_FLOAT, &
532&                                       dim3dids, ncvarid)
533            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
534&                                        shuffle=0,     &
535&                                        deflate=1,     &
536&                                        deflate_level=DEF_LEVEL)
537            ncret = nf90_put_var(ncid, ncvarid, &
538&                                tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
539
540            ncret = nf90_def_var(ncid, 'qv', NF90_FLOAT, &
541&                                       dim3dids, ncvarid)
542            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
543&                                        shuffle=0,     &
544&                                        deflate=1,     &
545&                                        deflate_level=DEF_LEVEL)
546            ncret = nf90_put_var(ncid, ncvarid, &
547&                                qv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
548
549            ncret = nf90_def_var(ncid, 'pv', NF90_FLOAT, &
550&                                       dim3dids, ncvarid)
551            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
552&                                        shuffle=0,     &
553&                                        deflate=1,     &
554&                                        deflate_level=DEF_LEVEL)
555            ncret = nf90_put_var(ncid, ncvarid, &
556&                                pv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
557
558            ncret = nf90_def_var(ncid, 'rho', NF90_FLOAT, &
559&                                       dim3dids, ncvarid)
560            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
561&                                        shuffle=0,     &
562&                                        deflate=1,     &
563&                                        deflate_level=DEF_LEVEL)
564            ncret = nf90_put_var(ncid, ncvarid, &
565&                                rho(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
566
567            ncret = nf90_def_var(ncid, 'drhodz', NF90_FLOAT, &
568&                                       dim3dids, ncvarid)
569            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
570&                                        shuffle=0,     &
571&                                        deflate=1,     &
572&                                        deflate_level=DEF_LEVEL)
573            ncret = nf90_put_var(ncid, ncvarid, &
574&                                drhodz(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
575
576            ncret = nf90_def_var(ncid, 'clouds', NF90_BYTE, &
577&                                       dim3dids, ncvarid)
578            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
579&                                        shuffle=0,     &
580&                                        deflate=1,     &
581&                                        deflate_level=DEF_LEVEL)
582            ncret = nf90_put_var(ncid, ncvarid, &
583&                                clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
584
585
586
587            ! Note the change in z dimension for the following
588            dim3dids = (/nxmax_dimid, nymax_dimid, nuvzmax_dimid/)
589
590            ncret = nf90_def_var(ncid, 'tth', NF90_FLOAT, &
591&                                       dim3dids, ncvarid)
592            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
593&                                        shuffle=0,     &
594&                                        deflate=1,     &
595&                                        deflate_level=DEF_LEVEL)
596            ncret = nf90_put_var(ncid, ncvarid, &
597&                                tth(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
598
599            ncret = nf90_def_var(ncid, 'qvh', NF90_FLOAT, &
600&                                       dim3dids, ncvarid)
601            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
602&                                        shuffle=0,     &
603&                                        deflate=1,     &
604&                                        deflate_level=DEF_LEVEL)
605            ncret = nf90_put_var(ncid, ncvarid, &
606&                                qvh(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
607
608            ncret = nf90_def_var(ncid, 'pplev', NF90_FLOAT, &
609&                                       dim3dids, ncvarid)
610            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
611&                                        shuffle=0,     &
612&                                        deflate=1,     &
613&                                        deflate_level=DEF_LEVEL)
614            ncret = nf90_put_var(ncid, ncvarid, &
615&                                pplev(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
616
617
618            dim2dids = (/nxmax_dimid, nymax_dimid/)
619            ncret = nf90_def_var(ncid, 'cloudsh', NF90_INT, &
620&                                       dim2dids, ncvarid)
621            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
622&                                        shuffle=0,     &
623&                                        deflate=1,     &
624&                                        deflate_level=DEF_LEVEL)
625            ncret = nf90_put_var(ncid, ncvarid, &
626&                                cloudsh(0:nxmax-1, 0:nymax-1, cm_index))
627
628
629
630            PRINT *, 'SUM(tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
631&                                        SUM(tt(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
632
633            PRINT *, 'SUM(clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
634&                                        SUM(clouds(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
635
636
637
638            ! 2d fields
639            WRITE(iounit) ps(:,:,:,cm_index)
640            WRITE(iounit) sd(:,:,:,cm_index)
641            WRITE(iounit) msl(:,:,:,cm_index)
642            WRITE(iounit) tcc(:,:,:,cm_index)
643            WRITE(iounit) u10(:,:,:,cm_index)
644            WRITE(iounit) v10(:,:,:,cm_index)
645            WRITE(iounit) tt2(:,:,:,cm_index)
646            WRITE(iounit) td2(:,:,:,cm_index)
647            WRITE(iounit) lsprec(:,:,:,cm_index)
648            WRITE(iounit) convprec(:,:,:,cm_index)
649            WRITE(iounit) sshf(:,:,:,cm_index)
650            WRITE(iounit) ssr(:,:,:,cm_index)
651            WRITE(iounit) surfstr(:,:,:,cm_index)
652            WRITE(iounit) ustar(:,:,:,cm_index)
653            WRITE(iounit) wstar(:,:,:,cm_index)
654            WRITE(iounit) hmix(:,:,:,cm_index)
655            WRITE(iounit) tropopause(:,:,:,cm_index)
656            WRITE(iounit) oli(:,:,:,cm_index)
657            WRITE(iounit) diffk(:,:,:,cm_index)
658            WRITE(iounit) vdep(:,:,:,cm_index)
659
660            dim2dids = (/nxmax_dimid, nymax_dimid/)
661
662            ncret = nf90_def_var(ncid, 'ps', NF90_FLOAT, &
663&                                       dim2dids, ncvarid)
664            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
665&                                        shuffle=0,     &
666&                                        deflate=1,     &
667&                                        deflate_level=DEF_LEVEL)
668            ncret = nf90_put_var(ncid, ncvarid, &
669&                                ps(0:nxmax-1, 0:nymax-1, 1, cm_index))
670
671            ncret = nf90_def_var(ncid, 'sd', NF90_FLOAT, &
672&                                       dim2dids, ncvarid)
673            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
674&                                        shuffle=0,     &
675&                                        deflate=1,     &
676&                                        deflate_level=DEF_LEVEL)
677            ncret = nf90_put_var(ncid, ncvarid, &
678&                                sd(0:nxmax-1, 0:nymax-1, 1, cm_index))
679
680            ncret = nf90_def_var(ncid, 'msl', NF90_FLOAT, &
681&                                       dim2dids, ncvarid)
682            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
683&                                        shuffle=0,     &
684&                                        deflate=1,     &
685&                                        deflate_level=DEF_LEVEL)
686            ncret = nf90_put_var(ncid, ncvarid, &
687&                                msl(0:nxmax-1, 0:nymax-1, 1, cm_index))
688
689            ncret = nf90_def_var(ncid, 'tcc', NF90_FLOAT, &
690&                                       dim2dids, ncvarid)
691            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
692&                                        shuffle=0,     &
693&                                        deflate=1,     &
694&                                        deflate_level=DEF_LEVEL)
695            ncret = nf90_put_var(ncid, ncvarid, &
696&                                tcc(0:nxmax-1, 0:nymax-1, 1, cm_index))
697
698            ncret = nf90_def_var(ncid, 'u10', NF90_FLOAT, &
699&                                       dim2dids, ncvarid)
700            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
701&                                        shuffle=0,     &
702&                                        deflate=1,     &
703&                                        deflate_level=DEF_LEVEL)
704            ncret = nf90_put_var(ncid, ncvarid, &
705&                                u10(0:nxmax-1, 0:nymax-1, 1, cm_index))
706
707            ncret = nf90_def_var(ncid, 'v10', NF90_FLOAT, &
708&                                       dim2dids, ncvarid)
709            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
710&                                        shuffle=0,     &
711&                                        deflate=1,     &
712&                                        deflate_level=DEF_LEVEL)
713            ncret = nf90_put_var(ncid, ncvarid, &
714&                                v10(0:nxmax-1, 0:nymax-1, 1, cm_index))
715
716            ncret = nf90_def_var(ncid, 'tt2', NF90_FLOAT, &
717&                                       dim2dids, ncvarid)
718            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
719&                                        shuffle=0,     &
720&                                        deflate=1,     &
721&                                        deflate_level=DEF_LEVEL)
722            ncret = nf90_put_var(ncid, ncvarid, &
723&                                tt2(0:nxmax-1, 0:nymax-1, 1, cm_index))
724
725            ncret = nf90_def_var(ncid, 'td2', NF90_FLOAT, &
726&                                       dim2dids, ncvarid)
727            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
728&                                        shuffle=0,     &
729&                                        deflate=1,     &
730&                                        deflate_level=DEF_LEVEL)
731            ncret = nf90_put_var(ncid, ncvarid, &
732&                                td2(0:nxmax-1, 0:nymax-1, 1, cm_index))
733
734            ncret = nf90_def_var(ncid, 'lsprec', NF90_FLOAT, &
735&                                       dim2dids, ncvarid)
736            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
737&                                        shuffle=0,     &
738&                                        deflate=1,     &
739&                                        deflate_level=DEF_LEVEL)
740            ncret = nf90_put_var(ncid, ncvarid, &
741&                                lsprec(0:nxmax-1, 0:nymax-1, 1, cm_index))
742
743            ncret = nf90_def_var(ncid, 'convprec', NF90_FLOAT, &
744&                                       dim2dids, ncvarid)
745            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
746&                                        shuffle=0,     &
747&                                        deflate=1,     &
748&                                        deflate_level=DEF_LEVEL)
749            ncret = nf90_put_var(ncid, ncvarid, &
750&                                convprec(0:nxmax-1, 0:nymax-1, 1, cm_index))
751
752            ncret = nf90_def_var(ncid, 'sshf', NF90_FLOAT, &
753&                                       dim2dids, ncvarid)
754            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
755&                                        shuffle=0,     &
756&                                        deflate=1,     &
757&                                        deflate_level=DEF_LEVEL)
758            ncret = nf90_put_var(ncid, ncvarid, &
759&                                sshf(0:nxmax-1, 0:nymax-1, 1, cm_index))
760
761            ncret = nf90_def_var(ncid, 'ssr', NF90_FLOAT, &
762&                                       dim2dids, ncvarid)
763            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
764&                                        shuffle=0,     &
765&                                        deflate=1,     &
766&                                        deflate_level=DEF_LEVEL)
767            ncret = nf90_put_var(ncid, ncvarid, &
768&                                ssr(0:nxmax-1, 0:nymax-1, 1, cm_index))
769
770            ncret = nf90_def_var(ncid, 'surfstr', NF90_FLOAT, &
771&                                       dim2dids, ncvarid)
772            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
773&                                        shuffle=0,     &
774&                                        deflate=1,     &
775&                                        deflate_level=DEF_LEVEL)
776            ncret = nf90_put_var(ncid, ncvarid, &
777&                                surfstr(0:nxmax-1, 0:nymax-1, 1, cm_index))
778
779            ncret = nf90_def_var(ncid, 'ustar', NF90_FLOAT, &
780&                                       dim2dids, ncvarid)
781            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
782&                                        shuffle=0,     &
783&                                        deflate=1,     &
784&                                        deflate_level=DEF_LEVEL)
785            ncret = nf90_put_var(ncid, ncvarid, &
786&                                ustar(0:nxmax-1, 0:nymax-1, 1, cm_index))
787
788            ncret = nf90_def_var(ncid, 'wstar', NF90_FLOAT, &
789&                                       dim2dids, ncvarid)
790            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
791&                                        shuffle=0,     &
792&                                        deflate=1,     &
793&                                        deflate_level=DEF_LEVEL)
794            ncret = nf90_put_var(ncid, ncvarid, &
795&                                wstar(0:nxmax-1, 0:nymax-1, 1, cm_index))
796
797            ncret = nf90_def_var(ncid, 'hmix', NF90_FLOAT, &
798&                                       dim2dids, ncvarid)
799            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
800&                                        shuffle=0,     &
801&                                        deflate=1,     &
802&                                        deflate_level=DEF_LEVEL)
803            ncret = nf90_put_var(ncid, ncvarid, &
804&                                hmix(0:nxmax-1, 0:nymax-1, 1, cm_index))
805
806            ncret = nf90_def_var(ncid, 'tropopause', NF90_FLOAT, &
807&                                       dim2dids, ncvarid)
808            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
809&                                        shuffle=0,     &
810&                                        deflate=1,     &
811&                                        deflate_level=DEF_LEVEL)
812            ncret = nf90_put_var(ncid, ncvarid, &
813&                                tropopause(0:nxmax-1, 0:nymax-1, 1, cm_index))
814
815            ncret = nf90_def_var(ncid, 'oli', NF90_FLOAT, &
816&                                       dim2dids, ncvarid)
817            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
818&                                        shuffle=0,     &
819&                                        deflate=1,     &
820&                                        deflate_level=DEF_LEVEL)
821            ncret = nf90_put_var(ncid, ncvarid, &
822&                                oli(0:nxmax-1, 0:nymax-1, 1, cm_index))
823
824            ncret = nf90_def_var(ncid, 'diffk', NF90_FLOAT, &
825&                                       dim2dids, ncvarid)
826            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
827&                                        shuffle=0,     &
828&                                        deflate=1,     &
829&                                        deflate_level=DEF_LEVEL)
830            ncret = nf90_put_var(ncid, ncvarid, &
831&                                diffk(0:nxmax-1, 0:nymax-1, 1, cm_index))
832
833
834
835            PRINT *, 'SUM(ps(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
836&                                        SUM(ps(0:nxmax-1,0:nymax-1,1, cm_index))
837
838            PRINT *, 'SUM(wstar(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
839&                                        SUM(wstar(0:nxmax-1,0:nymax-1,1, cm_index))
840
841
842            dim3dids = (/nxmax_dimid, nymax_dimid, maxspec_dimid/)
843
844            ncret = nf90_def_var(ncid, 'vdep', NF90_FLOAT, &
845&                                       dim3dids, ncvarid)
846            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
847&                                        shuffle=0,     &
848&                                        deflate=1,     &
849&                                        deflate_level=DEF_LEVEL)
850            ncret = nf90_put_var(ncid, ncvarid, &
851&                                vdep(0:nxmax-1, 0:nymax-1, 1:maxspec, cm_index))
852
853
854
855            ! 1d fields
856            WRITE(iounit) z0(:)
857            WRITE(iounit) akm(:)
858            WRITE(iounit) bkm(:)
859            WRITE(iounit) akz(:)
860            WRITE(iounit) bkz(:)
861            WRITE(iounit) aknew(:)
862            WRITE(iounit) bknew(:)
863
864            dim1dids = (/numclass_dimid/)
865
866            ncret = nf90_def_var(ncid, 'z0', NF90_FLOAT, &
867&                                       dim1dids, ncvarid)
868            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
869&                                        shuffle=0,     &
870&                                        deflate=1,     &
871&                                        deflate_level=DEF_LEVEL)
872            ncret = nf90_put_var(ncid, ncvarid, &
873&                                z0(1:numclass))
874
875
876            dim1dids = (/nwzmax_dimid/)
877
878            ncret = nf90_def_var(ncid, 'akm', NF90_FLOAT, &
879&                                       dim1dids, ncvarid)
880            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
881&                                        shuffle=0,     &
882&                                        deflate=1,     &
883&                                        deflate_level=DEF_LEVEL)
884            ncret = nf90_put_var(ncid, ncvarid, &
885&                                akm(1:nwzmax))
886
887            ncret = nf90_def_var(ncid, 'bkm', NF90_FLOAT, &
888&                                       dim1dids, ncvarid)
889            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
890&                                        shuffle=0,     &
891&                                        deflate=1,     &
892&                                        deflate_level=DEF_LEVEL)
893            ncret = nf90_put_var(ncid, ncvarid, &
894&                                bkm(1:nwzmax))
895
896
897            dim1dids = (/nuvzmax_dimid/)
898
899            ncret = nf90_def_var(ncid, 'akz', NF90_FLOAT, &
900&                                       dim1dids, ncvarid)
901            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
902&                                        shuffle=0,     &
903&                                        deflate=1,     &
904&                                        deflate_level=DEF_LEVEL)
905            ncret = nf90_put_var(ncid, ncvarid, &
906&                                akz(1:nuvzmax))
907
908            ncret = nf90_def_var(ncid, 'bkz', NF90_FLOAT, &
909&                                       dim1dids, ncvarid)
910            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
911&                                        shuffle=0,     &
912&                                        deflate=1,     &
913&                                        deflate_level=DEF_LEVEL)
914            ncret = nf90_put_var(ncid, ncvarid, &
915&                                bkz(1:nuvzmax))
916
917
918            dim1dids = (/nzmax_dimid/)
919
920            ncret = nf90_def_var(ncid, 'aknew', NF90_FLOAT, &
921&                                       dim1dids, ncvarid)
922            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
923&                                        shuffle=0,     &
924&                                        deflate=1,     &
925&                                        deflate_level=DEF_LEVEL)
926            ncret = nf90_put_var(ncid, ncvarid, &
927&                                aknew(1:nzmax))
928
929            ncret = nf90_def_var(ncid, 'bknew', NF90_FLOAT, &
930&                                       dim1dids, ncvarid)
931            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
932&                                        shuffle=0,     &
933&                                        deflate=1,     &
934&                                        deflate_level=DEF_LEVEL)
935            ncret = nf90_put_var(ncid, ncvarid, &
936&                                bknew(1:nzmax))
937
938
939            PRINT *, 'SUM(bknew(1:nzmax)): ', &
940&                                        SUM(bknew(1:nzmax))
941
942
943
944            ! Getting ready to add in nested code
945
946            ! These are compiled max dimensions from par_mod - these are
947            ! not meant to be reassigned during a LOAD, but used as "header"
948            ! information to provide the structure of arrays
949            ncret = nf90_def_dim(ncid, 'maxnests', maxnests, maxnests_dimid)
950            ncret = nf90_def_dim(ncid, 'nxmaxn', nxmaxn, nxmaxn_dimid)
951            ncret = nf90_def_dim(ncid, 'nymaxn', nymaxn, nymaxn_dimid)
952
953            WRITE(iounit) nxn(:)
954            WRITE(iounit) nyn(:)
955            WRITE(iounit) dxn(:)
956            WRITE(iounit) dyn(:)
957            WRITE(iounit) xlon0n(:)
958            WRITE(iounit) ylat0n(:)
959
960            ! Nested, scalar values (for each nest)
961
962            dim1dids = (/maxnests_dimid/)
963
964            ncret = nf90_def_var(ncid, 'nxn', NF90_INT, &
965&                                       dim1dids, ncvarid)
966            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
967&                                        shuffle=0,     &
968&                                        deflate=1,     &
969&                                        deflate_level=DEF_LEVEL)
970            ncret = nf90_put_var(ncid, ncvarid, &
971&                                nxn(1:maxnests))
972
973            ncret = nf90_def_var(ncid, 'nyn', NF90_INT, &
974&                                       dim1dids, ncvarid)
975            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
976&                                        shuffle=0,     &
977&                                        deflate=1,     &
978&                                        deflate_level=DEF_LEVEL)
979            ncret = nf90_put_var(ncid, ncvarid, &
980&                                nyn(1:maxnests))
981
982            ncret = nf90_def_var(ncid, 'dxn', NF90_FLOAT, &
983&                                       dim1dids, ncvarid)
984            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
985&                                        shuffle=0,     &
986&                                        deflate=1,     &
987&                                        deflate_level=DEF_LEVEL)
988            ncret = nf90_put_var(ncid, ncvarid, &
989&                                dxn(1:maxnests))
990
991            ncret = nf90_def_var(ncid, 'dyn', NF90_FLOAT, &
992&                                       dim1dids, ncvarid)
993            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
994&                                        shuffle=0,     &
995&                                        deflate=1,     &
996&                                        deflate_level=DEF_LEVEL)
997            ncret = nf90_put_var(ncid, ncvarid, &
998&                                dyn(1:maxnests))
999
1000            ncret = nf90_def_var(ncid, 'xlon0n', NF90_FLOAT, &
1001&                                       dim1dids, ncvarid)
1002            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1003&                                        shuffle=0,     &
1004&                                        deflate=1,     &
1005&                                        deflate_level=DEF_LEVEL)
1006            ncret = nf90_put_var(ncid, ncvarid, &
1007&                                xlon0n(1:maxnests))
1008
1009            ncret = nf90_def_var(ncid, 'ylat0n', NF90_FLOAT, &
1010&                                       dim1dids, ncvarid)
1011            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1012&                                        shuffle=0,     &
1013&                                        deflate=1,     &
1014&                                        deflate_level=DEF_LEVEL)
1015            ncret = nf90_put_var(ncid, ncvarid, &
1016&                                ylat0n(1:maxnests))
1017
1018
1019
1020
1021            ! Nested fields, static over time
1022            WRITE(iounit) oron, excessoron, lsmn, xlandusen
1023
1024            dim3dids = (/nxmaxn_dimid, nymaxn_dimid, maxnests_dimid/)
1025
1026            ncret = nf90_def_var(ncid, 'oron', NF90_FLOAT, &
1027&                                       dim3dids, ncvarid)
1028            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1029&                                        shuffle=0,     &
1030&                                        deflate=1,     &
1031&                                        deflate_level=DEF_LEVEL)
1032            ncret = nf90_put_var(ncid, ncvarid, &
1033&                                oron(0:nxmaxn-1, 0:nymaxn-1, 1:maxnests))
1034
1035            ncret = nf90_def_var(ncid, 'excessoron', NF90_FLOAT, &
1036&                                       dim3dids, ncvarid)
1037            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1038&                                        shuffle=0,     &
1039&                                        deflate=1,     &
1040&                                        deflate_level=DEF_LEVEL)
1041            ncret = nf90_put_var(ncid, ncvarid, &
1042&                                excessoron(0:nxmaxn-1, 0:nymaxn-1, 1:maxnests))
1043
1044            ncret = nf90_def_var(ncid, 'lsmn', NF90_FLOAT, &
1045&                                       dim3dids, ncvarid)
1046            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1047&                                        shuffle=0,     &
1048&                                        deflate=1,     &
1049&                                        deflate_level=DEF_LEVEL)
1050            ncret = nf90_put_var(ncid, ncvarid, &
1051&                                lsmn(0:nxmaxn-1, 0:nymaxn-1, 1:maxnests))
1052
1053            dim4dids = (/nxmaxn_dimid, nymaxn_dimid, numclass_dimid, maxnests_dimid/)
1054
1055            ncret = nf90_def_var(ncid, 'xlandusen', NF90_FLOAT, &
1056&                                       dim4dids, ncvarid)
1057            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1058&                                        shuffle=0,     &
1059&                                        deflate=1,     &
1060&                                        deflate_level=DEF_LEVEL)
1061            ncret = nf90_put_var(ncid, ncvarid, &
1062&                                xlandusen(0:nxmaxn-1, 0:nymaxn-1, 1:numclass, 1:maxnests))
1063
1064            PRINT *, 'SUM(xlandusen): ', &
1065&                                        SUM(xlandusen)
1066
1067
1068
1069            ! 3d nested fields
1070            WRITE(iounit) uun(:,:,:,cm_index,:)
1071            WRITE(iounit) vvn(:,:,:,cm_index,:)
1072            WRITE(iounit) wwn(:,:,:,cm_index,:)
1073            WRITE(iounit) ttn(:,:,:,cm_index,:)
1074            WRITE(iounit) qvn(:,:,:,cm_index,:)
1075            WRITE(iounit) pvn(:,:,:,cm_index,:)
1076            WRITE(iounit) cloudsn(:,:,:,cm_index,:)
1077            WRITE(iounit) cloudsnh(:,:,cm_index,:)
1078            WRITE(iounit) rhon(:,:,:,cm_index,:)
1079            WRITE(iounit) drhodzn(:,:,:,cm_index,:)
1080            WRITE(iounit) tthn(:,:,:,cm_index,:)
1081            WRITE(iounit) qvhn(:,:,:,cm_index,:)
1082
1083            ! 2d nested fields
1084            WRITE(iounit) psn(:,:,:,cm_index,:)
1085            WRITE(iounit) sdn(:,:,:,cm_index,:)
1086            WRITE(iounit) msln(:,:,:,cm_index,:)
1087            WRITE(iounit) tccn(:,:,:,cm_index,:)
1088            WRITE(iounit) u10n(:,:,:,cm_index,:)
1089            WRITE(iounit) v10n(:,:,:,cm_index,:)
1090            WRITE(iounit) tt2n(:,:,:,cm_index,:)
1091            WRITE(iounit) td2n(:,:,:,cm_index,:)
1092            WRITE(iounit) lsprecn(:,:,:,cm_index,:)
1093            WRITE(iounit) convprecn(:,:,:,cm_index,:)
1094            WRITE(iounit) sshfn(:,:,:,cm_index,:)
1095            WRITE(iounit) ssrn(:,:,:,cm_index,:)
1096            WRITE(iounit) surfstrn(:,:,:,cm_index,:)
1097            WRITE(iounit) ustarn(:,:,:,cm_index,:)
1098            WRITE(iounit) wstarn(:,:,:,cm_index,:)
1099            WRITE(iounit) hmixn(:,:,:,cm_index,:)
1100            WRITE(iounit) tropopausen(:,:,:,cm_index,:)
1101            WRITE(iounit) olin(:,:,:,cm_index,:)
1102            WRITE(iounit) diffkn(:,:,:,cm_index,:)
1103            WRITE(iounit) vdepn(:,:,:,cm_index,:)
1104
1105            ! Auxiliary variables for nests
1106            WRITE(iounit) xresoln(:)
1107            WRITE(iounit) yresoln(:)
1108            WRITE(iounit) xln(:)
1109            WRITE(iounit) yln(:)
1110            WRITE(iounit) xrn(:)
1111            WRITE(iounit) yrn(:)
1112
1113            ! Variables for polar stereographic projection
1114            WRITE(iounit) xglobal, sglobal, nglobal
1115            WRITE(iounit) switchnorthg, switchsouthg
1116            WRITE(iounit) southpolemap(:)
1117            WRITE(iounit) northpolemap(:)
1118
1119            ! Variables declared in conv_mod (convection)
1120            WRITE(iounit) pconv(:)
1121            WRITE(iounit) phconv(:)
1122            WRITE(iounit) dpr(:)
1123            WRITE(iounit) pconv_hpa(:)
1124            WRITE(iounit) phconv_hpa(:)
1125            WRITE(iounit) ft(:)
1126            WRITE(iounit) fq(:)
1127            WRITE(iounit) fmass(:,:)
1128            WRITE(iounit) sub(:)
1129            WRITE(iounit) fmassfrac(:,:)
1130            WRITE(iounit) cbaseflux(:,:)
1131            WRITE(iounit) cbasefluxn(:,:,:)
1132            WRITE(iounit) tconv(:)
1133            WRITE(iounit) qconv(:)
1134            WRITE(iounit) qsconv(:)
1135            WRITE(iounit) psconv, tt2conv, td2conv
1136            WRITE(iounit) nconvlev, nconvtop
1137
1138        ELSE IF (op == 'LOAD') THEN
1139
1140            ! Read the preprocessed format version string and insure it
1141            ! matches this version
1142            READ (iounit) temp_preproc_format_version_str
1143            PRINT *, 'Reading preprocessed file format version: ', &
1144&                    temp_preproc_format_version_str
1145
1146            IF (TRIM(temp_preproc_format_version_str) == &
1147&                                        TRIM(PREPROC_FORMAT_VERSION_STR)) THEN
1148                CONTINUE
1149            ELSE
1150                ! PRINT *, ''  GK: causes relocation truncated to fit: R_X86_64_32
1151                PRINT *, 'Inconsistent preprocessing format version'
1152                PRINT *, 'Expected Version: ', PREPROC_FORMAT_VERSION_STR
1153                PRINT *, 'Detected Version: ', temp_preproc_format_version_str
1154                ! PRINT *, ''
1155                STOP
1156            END IF
1157
1158            ! Read the compiled max dimensions that were dumped from par_mod
1159            ! when creating the fp file, so that we can compare against
1160            ! current FLEXPART dimensions - they need to be the same, or else
1161            ! we abort.
1162            READ (iounit) temp_nxmax, temp_nymax, temp_nzmax, &
1163&                         temp_nuvzmax, temp_nwzmax
1164
1165            ! Get dimensions
1166            ncret = nf90_inq_dimid(ncid, 'nxmax', nxmax_dimid)
1167            ncret = nf90_inquire_dimension(ncid, nxmax_dimid, nxmax_dimname, &
1168&                                                temp_nxmax)
1169            PRINT *, 'temp_nxmax: ', temp_nxmax
1170
1171            ncret = nf90_inq_dimid(ncid, 'nymax', nymax_dimid)
1172            ncret = nf90_inquire_dimension(ncid, nymax_dimid, nymax_dimname, &
1173&                                                temp_nymax)
1174            PRINT *, 'temp_nymax: ', temp_nymax
1175
1176            ncret = nf90_inq_dimid(ncid, 'nzmax', nzmax_dimid)
1177            ncret = nf90_inquire_dimension(ncid, nzmax_dimid, nzmax_dimname, &
1178&                                                temp_nzmax)
1179            PRINT *, 'temp_nzmax: ', temp_nzmax
1180
1181            ncret = nf90_inq_dimid(ncid, 'nuvzmax', nuvzmax_dimid)
1182            ncret = nf90_inquire_dimension(ncid, nuvzmax_dimid, nuvzmax_dimname, &
1183&                                                temp_nuvzmax)
1184            PRINT *, 'temp_nuvzmax: ', temp_nuvzmax
1185
1186            ncret = nf90_inq_dimid(ncid, 'nwzmax', nwzmax_dimid)
1187            ncret = nf90_inquire_dimension(ncid, nwzmax_dimid, nwzmax_dimname, &
1188&                                                temp_nwzmax)
1189            PRINT *, 'temp_nwzmax: ', temp_nwzmax
1190
1191            ncret = nf90_inq_dimid(ncid, 'numclass', numclass_dimid)
1192            ncret = nf90_inquire_dimension(ncid, numclass_dimid, numclass_dimname, &
1193&                                                temp_numclass)
1194            PRINT *, 'temp_numclass: ', temp_numclass
1195
1196            ncret = nf90_inq_dimid(ncid, 'maxspec', maxspec_dimid)
1197            ncret = nf90_inquire_dimension(ncid, maxspec_dimid, maxspec_dimname, &
1198&                                                temp_maxspec)
1199            PRINT *, 'temp_maxspec: ', temp_maxspec
1200
1201
1202
1203            IF ( (temp_nxmax == nxmax) .AND. (temp_nymax == nymax) .AND. &
1204&                   (temp_nzmax == nzmax) .AND. &
1205&                   (temp_nuvzmax == nuvzmax) .AND. &
1206&                   (temp_nwzmax == nwzmax) .AND. &
1207&                   (temp_numclass == numclass) .AND. &
1208&                   (temp_maxspec == maxspec) ) THEN
1209                CONTINUE
1210            ELSE
1211                PRINT *, 'Incompatible dimensions between fp file and current FLEXPART!'
1212                ! PRINT *, ''
1213                PRINT *, '                  FP file     Compiled FP'
1214                PRINT *, 'nxmax:     ', temp_nxmax, '    ', nxmax
1215                PRINT *, 'nymax:     ', temp_nymax, '    ', nymax
1216                PRINT *, 'nzmax:     ', temp_nzmax, '    ', nzmax
1217                PRINT *, 'nuvzmax:     ', temp_nuvzmax, '    ', nuvzmax
1218                PRINT *, 'nwzmax:     ', temp_nwzmax, '    ', nwzmax
1219                PRINT *, 'numclass:     ', temp_numclass, '    ', numclass
1220                PRINT *, 'maxspec:     ', temp_maxspec, '    ', maxspec
1221                ! PRINT *, ''
1222                STOP
1223            END IF
1224
1225
1226
1227
1228            ! Scalar values
1229            READ(iounit) nx, ny, nxmin1, nymin1, nxfield
1230            READ(iounit) nuvz, nwz, nz, nmixz, nlev_ec
1231            READ(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst
1232
1233
1234
1235            ! Get the varid , then read into scalar variable
1236            ncret = nf90_inq_varid(ncid, 'nx', ncvarid)
1237            ncret = nf90_get_var(ncid, ncvarid, nx)
1238
1239            ncret = nf90_inq_varid(ncid, 'ny', ncvarid)
1240            ncret = nf90_get_var(ncid, ncvarid, ny)
1241
1242            ncret = nf90_inq_varid(ncid, 'nxmin1', ncvarid)
1243            ncret = nf90_get_var(ncid, ncvarid, nxmin1)
1244
1245            ncret = nf90_inq_varid(ncid, 'nymin1', ncvarid)
1246            ncret = nf90_get_var(ncid, ncvarid, nymin1)
1247
1248            ncret = nf90_inq_varid(ncid, 'nxfield', ncvarid)
1249            ncret = nf90_get_var(ncid, ncvarid, nxfield)
1250
1251            ncret = nf90_inq_varid(ncid, 'nuvz', ncvarid)
1252            ncret = nf90_get_var(ncid, ncvarid, nuvz)
1253
1254            ncret = nf90_inq_varid(ncid, 'nwz', ncvarid)
1255            ncret = nf90_get_var(ncid, ncvarid, nwz)
1256
1257            ncret = nf90_inq_varid(ncid, 'nz', ncvarid)
1258            ncret = nf90_get_var(ncid, ncvarid, nz)
1259
1260            ncret = nf90_inq_varid(ncid, 'nmixz', ncvarid)
1261            ncret = nf90_get_var(ncid, ncvarid, nmixz)
1262
1263            ncret = nf90_inq_varid(ncid, 'nlev_ec', ncvarid)
1264            ncret = nf90_get_var(ncid, ncvarid, nlev_ec)
1265
1266            ncret = nf90_inq_varid(ncid, 'dx', ncvarid)
1267            ncret = nf90_get_var(ncid, ncvarid, dx)
1268
1269            ncret = nf90_inq_varid(ncid, 'dy', ncvarid)
1270            ncret = nf90_get_var(ncid, ncvarid, dy)
1271
1272            ncret = nf90_inq_varid(ncid, 'xlon0', ncvarid)
1273            ncret = nf90_get_var(ncid, ncvarid, xlon0)
1274
1275            ncret = nf90_inq_varid(ncid, 'ylat0', ncvarid)
1276            ncret = nf90_get_var(ncid, ncvarid, ylat0)
1277
1278            ncret = nf90_inq_varid(ncid, 'dxconst', ncvarid)
1279            ncret = nf90_get_var(ncid, ncvarid, dxconst)
1280
1281            ncret = nf90_inq_varid(ncid, 'dyconst', ncvarid)
1282            ncret = nf90_get_var(ncid, ncvarid, dyconst)
1283
1284
1285
1286
1287
1288
1289            ! Fixed fields, static in time
1290            READ(iounit) oro, excessoro, lsm, xlanduse, height
1291
1292            ncret = nf90_inq_varid(ncid, 'oro', ncvarid)
1293            ncret = nf90_get_var(ncid, ncvarid, oro(0:nxmax-1,0:nymax-1))
1294
1295            ncret = nf90_inq_varid(ncid, 'excessoro', ncvarid)
1296            ncret = nf90_get_var(ncid, ncvarid, excessoro(0:nxmax-1,0:nymax-1))
1297
1298            ncret = nf90_inq_varid(ncid, 'lsm', ncvarid)
1299            ncret = nf90_get_var(ncid, ncvarid, lsm(0:nxmax-1,0:nymax-1))
1300
1301            ncret = nf90_inq_varid(ncid, 'xlanduse', ncvarid)
1302            ncret = nf90_get_var(ncid, ncvarid, xlanduse(0:nxmax-1,0:nymax-1, 1:numclass))
1303
1304            ncret = nf90_inq_varid(ncid, 'height', ncvarid)
1305            ncret = nf90_get_var(ncid, ncvarid, height(1:nzmax))
1306
1307
1308
1309
1310            ! 3d fields
1311            READ(iounit) uu(:,:,:,cm_index)
1312            READ(iounit) vv(:,:,:,cm_index)
1313            READ(iounit) uupol(:,:,:,cm_index)
1314            READ(iounit) vvpol(:,:,:,cm_index)
1315            READ(iounit) ww(:,:,:,cm_index)
1316            READ(iounit) tt(:,:,:,cm_index)
1317            READ(iounit) qv(:,:,:,cm_index)
1318            READ(iounit) pv(:,:,:,cm_index)
1319            READ(iounit) rho(:,:,:,cm_index)
1320            READ(iounit) drhodz(:,:,:,cm_index)
1321            READ(iounit) tth(:,:,:,cm_index)
1322            READ(iounit) qvh(:,:,:,cm_index)
1323            READ(iounit) pplev(:,:,:,cm_index)
1324            READ(iounit) clouds(:,:,:,cm_index)
1325            READ(iounit) cloudsh(:,:,cm_index)
1326
1327
1328
1329
1330            ! Get the varid and read the variable into the array
1331            ncret = nf90_inq_varid(ncid, 'uu', ncvarid)
1332            ncret = nf90_get_var(ncid, ncvarid, uu(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1333
1334            ncret = nf90_inq_varid(ncid, 'vv', ncvarid)
1335            ncret = nf90_get_var(ncid, ncvarid, vv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1336
1337            ncret = nf90_inq_varid(ncid, 'uupol', ncvarid)
1338            ncret = nf90_get_var(ncid, ncvarid, uupol(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1339
1340            ncret = nf90_inq_varid(ncid, 'vvpol', ncvarid)
1341            ncret = nf90_get_var(ncid, ncvarid, vvpol(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1342
1343            ncret = nf90_inq_varid(ncid, 'ww', ncvarid)
1344            ncret = nf90_get_var(ncid, ncvarid, ww(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1345
1346            ncret = nf90_inq_varid(ncid, 'tt', ncvarid)
1347            ncret = nf90_get_var(ncid, ncvarid, tt(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1348
1349            ncret = nf90_inq_varid(ncid, 'qv', ncvarid)
1350            ncret = nf90_get_var(ncid, ncvarid, qv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1351
1352            ncret = nf90_inq_varid(ncid, 'pv', ncvarid)
1353            ncret = nf90_get_var(ncid, ncvarid, pv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1354
1355            ncret = nf90_inq_varid(ncid, 'rho', ncvarid)
1356            ncret = nf90_get_var(ncid, ncvarid, rho(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1357
1358            ncret = nf90_inq_varid(ncid, 'drhodz', ncvarid)
1359            ncret = nf90_get_var(ncid, ncvarid, drhodz(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1360
1361            ncret = nf90_inq_varid(ncid, 'clouds', ncvarid)
1362            ncret = nf90_get_var(ncid, ncvarid, clouds(0:nxmax-1,0:nzmax-1,1:nzmax,cm_index))
1363
1364            ncret = nf90_inq_varid(ncid, 'tth', ncvarid)
1365            ncret = nf90_get_var(ncid, ncvarid, tth(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
1366
1367            ncret = nf90_inq_varid(ncid, 'qvh', ncvarid)
1368            ncret = nf90_get_var(ncid, ncvarid, qvh(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
1369
1370            ncret = nf90_inq_varid(ncid, 'pplev', ncvarid)
1371            ncret = nf90_get_var(ncid, ncvarid, pplev(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
1372
1373            ncret = nf90_inq_varid(ncid, 'cloudsh', ncvarid)
1374            ncret = nf90_get_var(ncid, ncvarid, cloudsh(0:nxmax-1,0:nymax-1,cm_index))
1375
1376
1377
1378
1379            PRINT *, 'SUM(tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
1380&                                        SUM(tt(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
1381
1382
1383            PRINT *, 'SUM(clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
1384&                                        SUM(clouds(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
1385
1386
1387
1388            ! 2d fields
1389            READ(iounit) ps(:,:,:,cm_index)
1390            READ(iounit) sd(:,:,:,cm_index)
1391            READ(iounit) msl(:,:,:,cm_index)
1392            READ(iounit) tcc(:,:,:,cm_index)
1393            READ(iounit) u10(:,:,:,cm_index)
1394            READ(iounit) v10(:,:,:,cm_index)
1395            READ(iounit) tt2(:,:,:,cm_index)
1396            READ(iounit) td2(:,:,:,cm_index)
1397            READ(iounit) lsprec(:,:,:,cm_index)
1398            READ(iounit) convprec(:,:,:,cm_index)
1399            READ(iounit) sshf(:,:,:,cm_index)
1400            READ(iounit) ssr(:,:,:,cm_index)
1401            READ(iounit) surfstr(:,:,:,cm_index)
1402            READ(iounit) ustar(:,:,:,cm_index)
1403            READ(iounit) wstar(:,:,:,cm_index)
1404            READ(iounit) hmix(:,:,:,cm_index)
1405            READ(iounit) tropopause(:,:,:,cm_index)
1406            READ(iounit) oli(:,:,:,cm_index)
1407            READ(iounit) diffk(:,:,:,cm_index)
1408            READ(iounit) vdep(:,:,:,cm_index)
1409
1410            ! Get the varid and read the variable into the array
1411            ncret = nf90_inq_varid(ncid, 'ps', ncvarid)
1412            ncret = nf90_get_var(ncid, ncvarid, ps(0:nxmax-1,0:nymax-1,1,cm_index))
1413
1414            ncret = nf90_inq_varid(ncid, 'sd', ncvarid)
1415            ncret = nf90_get_var(ncid, ncvarid, sd(0:nxmax-1,0:nymax-1,1,cm_index))
1416
1417            ncret = nf90_inq_varid(ncid, 'msl', ncvarid)
1418            ncret = nf90_get_var(ncid, ncvarid, msl(0:nxmax-1,0:nymax-1,1,cm_index))
1419
1420            ncret = nf90_inq_varid(ncid, 'tcc', ncvarid)
1421            ncret = nf90_get_var(ncid, ncvarid, tcc(0:nxmax-1,0:nymax-1,1,cm_index))
1422
1423            ncret = nf90_inq_varid(ncid, 'u10', ncvarid)
1424            ncret = nf90_get_var(ncid, ncvarid, u10(0:nxmax-1,0:nymax-1,1,cm_index))
1425
1426            ncret = nf90_inq_varid(ncid, 'v10', ncvarid)
1427            ncret = nf90_get_var(ncid, ncvarid, v10(0:nxmax-1,0:nymax-1,1,cm_index))
1428
1429            ncret = nf90_inq_varid(ncid, 'tt2', ncvarid)
1430            ncret = nf90_get_var(ncid, ncvarid, tt2(0:nxmax-1,0:nymax-1,1,cm_index))
1431
1432            ncret = nf90_inq_varid(ncid, 'td2', ncvarid)
1433            ncret = nf90_get_var(ncid, ncvarid, td2(0:nxmax-1,0:nymax-1,1,cm_index))
1434
1435            ncret = nf90_inq_varid(ncid, 'lsprec', ncvarid)
1436            ncret = nf90_get_var(ncid, ncvarid, lsprec(0:nxmax-1,0:nymax-1,1,cm_index))
1437
1438            ncret = nf90_inq_varid(ncid, 'convprec', ncvarid)
1439            ncret = nf90_get_var(ncid, ncvarid, convprec(0:nxmax-1,0:nymax-1,1,cm_index))
1440
1441            ncret = nf90_inq_varid(ncid, 'sshf', ncvarid)
1442            ncret = nf90_get_var(ncid, ncvarid, sshf(0:nxmax-1,0:nymax-1,1,cm_index))
1443
1444            ncret = nf90_inq_varid(ncid, 'ssr', ncvarid)
1445            ncret = nf90_get_var(ncid, ncvarid, ssr(0:nxmax-1,0:nymax-1,1,cm_index))
1446
1447            ncret = nf90_inq_varid(ncid, 'surfstr', ncvarid)
1448            ncret = nf90_get_var(ncid, ncvarid, surfstr(0:nxmax-1,0:nymax-1,1,cm_index))
1449
1450            ncret = nf90_inq_varid(ncid, 'ustar', ncvarid)
1451            ncret = nf90_get_var(ncid, ncvarid, ustar(0:nxmax-1,0:nymax-1,1,cm_index))
1452
1453            ncret = nf90_inq_varid(ncid, 'wstar', ncvarid)
1454            ncret = nf90_get_var(ncid, ncvarid, wstar(0:nxmax-1,0:nymax-1,1,cm_index))
1455
1456            ncret = nf90_inq_varid(ncid, 'hmix', ncvarid)
1457            ncret = nf90_get_var(ncid, ncvarid, hmix(0:nxmax-1,0:nymax-1,1,cm_index))
1458
1459            ncret = nf90_inq_varid(ncid, 'tropopause', ncvarid)
1460            ncret = nf90_get_var(ncid, ncvarid, tropopause(0:nxmax-1,0:nymax-1,1,cm_index))
1461
1462            ncret = nf90_inq_varid(ncid, 'oli', ncvarid)
1463            ncret = nf90_get_var(ncid, ncvarid, oli(0:nxmax-1,0:nymax-1,1,cm_index))
1464
1465            ncret = nf90_inq_varid(ncid, 'diffk', ncvarid)
1466            ncret = nf90_get_var(ncid, ncvarid, diffk(0:nxmax-1,0:nymax-1,1,cm_index))
1467
1468
1469
1470
1471            PRINT *, 'SUM(ps(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
1472&                                        SUM(ps(0:nxmax-1,0:nymax-1,1, cm_index))
1473
1474            PRINT *, 'SUM(wstar(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
1475&                                        SUM(wstar(0:nxmax-1,0:nymax-1,1, cm_index))
1476
1477
1478            ncret = nf90_inq_varid(ncid, 'vdep', ncvarid)
1479            ncret = nf90_get_var(ncid, ncvarid, vdep(0:nxmax-1,0:nymax-1,1:maxspec, cm_index))
1480
1481
1482
1483
1484
1485            ! 1d fields
1486            READ(iounit) z0(:)
1487            READ(iounit) akm(:)
1488            READ(iounit) bkm(:)
1489            READ(iounit) akz(:)
1490            READ(iounit) bkz(:)
1491            READ(iounit) aknew(:)
1492            READ(iounit) bknew(:)
1493
1494
1495            ncret = nf90_inq_varid(ncid, 'z0', ncvarid)
1496            ncret = nf90_get_var(ncid, ncvarid, z0(1:numclass))
1497
1498            ncret = nf90_inq_varid(ncid, 'akm', ncvarid)
1499            ncret = nf90_get_var(ncid, ncvarid, akm(1:nwzmax))
1500
1501            ncret = nf90_inq_varid(ncid, 'bkm', ncvarid)
1502            ncret = nf90_get_var(ncid, ncvarid, bkm(1:nwzmax))
1503
1504            ncret = nf90_inq_varid(ncid, 'akz', ncvarid)
1505            ncret = nf90_get_var(ncid, ncvarid, akz(1:nuvzmax))
1506
1507            ncret = nf90_inq_varid(ncid, 'bkz', ncvarid)
1508            ncret = nf90_get_var(ncid, ncvarid, bkz(1:nuvzmax))
1509
1510            ncret = nf90_inq_varid(ncid, 'aknew', ncvarid)
1511            ncret = nf90_get_var(ncid, ncvarid, aknew(1:nzmax))
1512
1513            ncret = nf90_inq_varid(ncid, 'bknew', ncvarid)
1514            ncret = nf90_get_var(ncid, ncvarid, bknew(1:nzmax))
1515
1516
1517
1518            PRINT *, 'SUM(bknew(1:nzmax)): ', &
1519&                                        SUM(bknew(1:nzmax))
1520
1521
1522
1523
1524            ! Now the nested input grid variables
1525            ! Get the compiled values that were written into the FP file, and
1526            ! make sure they are equal to the current compiled values, to make
1527            ! sure we are working with consistent arrays
1528            ncret = nf90_inq_dimid(ncid, 'maxnests', maxnests_dimid)
1529            ncret = nf90_inquire_dimension(ncid, maxnests_dimid, maxnests_dimname, &
1530&                                                temp_maxnests)
1531            PRINT *, 'temp_maxnests: ', temp_maxnests
1532
1533            ncret = nf90_inq_dimid(ncid, 'nxmaxn', nxmaxn_dimid)
1534            ncret = nf90_inquire_dimension(ncid, nxmaxn_dimid, nxmaxn_dimname, &
1535&                                                temp_nxmaxn)
1536            PRINT *, 'temp_nxmaxn: ', temp_nxmaxn
1537
1538            ncret = nf90_inq_dimid(ncid, 'nymaxn', nymaxn_dimid)
1539            ncret = nf90_inquire_dimension(ncid, nymaxn_dimid, nymaxn_dimname, &
1540&                                                temp_nymaxn)
1541            PRINT *, 'temp_nymaxn: ', temp_nymaxn
1542
1543            ! Note that maxspec_dimid and numclass_dimid were checked above
1544
1545
1546
1547
1548            IF ( (temp_nxmaxn == nxmaxn) .AND. (temp_nymaxn == nymaxn) .AND. &
1549&                   (temp_maxnests == maxnests) ) THEN
1550                CONTINUE
1551            ELSE
1552                PRINT *, 'Incompatible dimensions between fp file and current FLEXPART!'
1553                ! PRINT *, ''
1554                PRINT *, '                  FP file     Compiled FP'
1555                PRINT *, 'nxmaxn:     ', temp_nxmaxn, '    ', nxmaxn
1556                PRINT *, 'nymaxn:     ', temp_nymaxn, '    ', nymaxn
1557                PRINT *, 'maxnests:     ', temp_maxnests, '    ', maxnests
1558                STOP
1559            END IF
1560
1561
1562
1563            ! Nested, scalar values (for each nest)
1564            READ(iounit) nxn(:)
1565            READ(iounit) nyn(:)
1566            READ(iounit) dxn(:)
1567            READ(iounit) dyn(:)
1568            READ(iounit) xlon0n(:)
1569            READ(iounit) ylat0n(:)
1570
1571
1572            ncret = nf90_inq_varid(ncid, 'nxn', ncvarid)
1573            ncret = nf90_get_var(ncid, ncvarid, nxn(1:maxnests))
1574
1575            ncret = nf90_inq_varid(ncid, 'nyn', ncvarid)
1576            ncret = nf90_get_var(ncid, ncvarid, nyn(1:maxnests))
1577
1578            ncret = nf90_inq_varid(ncid, 'dxn', ncvarid)
1579            ncret = nf90_get_var(ncid, ncvarid, dxn(1:maxnests))
1580
1581            ncret = nf90_inq_varid(ncid, 'dyn', ncvarid)
1582            ncret = nf90_get_var(ncid, ncvarid, dyn(1:maxnests))
1583
1584            ncret = nf90_inq_varid(ncid, 'xlon0n', ncvarid)
1585            ncret = nf90_get_var(ncid, ncvarid, xlon0n(1:maxnests))
1586
1587            ncret = nf90_inq_varid(ncid, 'ylat0n', ncvarid)
1588            ncret = nf90_get_var(ncid, ncvarid, ylat0n(1:maxnests))
1589
1590
1591
1592            ! Nested fields, static over time
1593            READ(iounit) oron, excessoron, lsmn, xlandusen
1594
1595            ncret = nf90_inq_varid(ncid, 'oron', ncvarid)
1596            ncret = nf90_get_var(ncid, ncvarid, oron(0:nxmaxn-1,0:nymaxn-1,1:maxnests))
1597
1598            ncret = nf90_inq_varid(ncid, 'excessoron', ncvarid)
1599            ncret = nf90_get_var(ncid, ncvarid, excessoron(0:nxmaxn-1,0:nymaxn-1,1:maxnests))
1600
1601            ncret = nf90_inq_varid(ncid, 'lsmn', ncvarid)
1602            ncret = nf90_get_var(ncid, ncvarid, lsmn(0:nxmaxn-1,0:nymaxn-1,1:maxnests))
1603
1604            ncret = nf90_inq_varid(ncid, 'xlandusen', ncvarid)
1605            ncret = nf90_get_var(ncid, ncvarid, xlandusen(0:nxmaxn-1,0:nymaxn-1,1:numclass,1:maxnests))
1606
1607
1608            PRINT *, 'SUM(xlandusen): ', &
1609&                                        SUM(xlandusen)
1610
1611
1612
1613
1614            ! 3d nested fields
1615            READ(iounit) uun(:,:,:,cm_index,:)
1616            READ(iounit) vvn(:,:,:,cm_index,:)
1617            READ(iounit) wwn(:,:,:,cm_index,:)
1618            READ(iounit) ttn(:,:,:,cm_index,:)
1619            READ(iounit) qvn(:,:,:,cm_index,:)
1620            READ(iounit) pvn(:,:,:,cm_index,:)
1621            READ(iounit) cloudsn(:,:,:,cm_index,:)
1622            READ(iounit) cloudsnh(:,:,cm_index,:)
1623            READ(iounit) rhon(:,:,:,cm_index,:)
1624            READ(iounit) drhodzn(:,:,:,cm_index,:)
1625            READ(iounit) tthn(:,:,:,cm_index,:)
1626            READ(iounit) qvhn(:,:,:,cm_index,:)
1627
1628            ! 2d nested fields
1629            READ(iounit) psn(:,:,:,cm_index,:)
1630            READ(iounit) sdn(:,:,:,cm_index,:)
1631            READ(iounit) msln(:,:,:,cm_index,:)
1632            READ(iounit) tccn(:,:,:,cm_index,:)
1633            READ(iounit) u10n(:,:,:,cm_index,:)
1634            READ(iounit) v10n(:,:,:,cm_index,:)
1635            READ(iounit) tt2n(:,:,:,cm_index,:)
1636            READ(iounit) td2n(:,:,:,cm_index,:)
1637            READ(iounit) lsprecn(:,:,:,cm_index,:)
1638            READ(iounit) convprecn(:,:,:,cm_index,:)
1639            READ(iounit) sshfn(:,:,:,cm_index,:)
1640            READ(iounit) ssrn(:,:,:,cm_index,:)
1641            READ(iounit) surfstrn(:,:,:,cm_index,:)
1642            READ(iounit) ustarn(:,:,:,cm_index,:)
1643            READ(iounit) wstarn(:,:,:,cm_index,:)
1644            READ(iounit) hmixn(:,:,:,cm_index,:)
1645            READ(iounit) tropopausen(:,:,:,cm_index,:)
1646            READ(iounit) olin(:,:,:,cm_index,:)
1647            READ(iounit) diffkn(:,:,:,cm_index,:)
1648            READ(iounit) vdepn(:,:,:,cm_index,:)
1649
1650            ! Auxiliary variables for nests
1651            READ(iounit) xresoln(:)
1652            READ(iounit) yresoln(:)
1653            READ(iounit) xln(:)
1654            READ(iounit) yln(:)
1655            READ(iounit) xrn(:)
1656            READ(iounit) yrn(:)
1657
1658            ! Variables for polar stereographic projection
1659            READ(iounit) xglobal, sglobal, nglobal
1660            READ(iounit) switchnorthg, switchsouthg
1661            READ(iounit) southpolemap(:)
1662            READ(iounit) northpolemap(:)
1663
1664            ! Variables declared in conv_mod (convection)
1665            READ(iounit) pconv(:)
1666            READ(iounit) phconv(:)
1667            READ(iounit) dpr(:)
1668            READ(iounit) pconv_hpa(:)
1669            READ(iounit) phconv_hpa(:)
1670            READ(iounit) ft(:)
1671            READ(iounit) fq(:)
1672            READ(iounit) fmass(:,:)
1673            READ(iounit) sub(:)
1674            READ(iounit) fmassfrac(:,:)
1675            READ(iounit) cbaseflux(:,:)
1676            READ(iounit) cbasefluxn(:,:,:)
1677            READ(iounit) tconv(:)
1678            READ(iounit) qconv(:)
1679            READ(iounit) qsconv(:)
1680            READ(iounit) psconv, tt2conv, td2conv
1681            READ(iounit) nconvlev, nconvtop
1682
1683        ELSE
1684            STOP 'fpio(): Illegal operation'
1685           
1686        ENDIF
1687    END SUBROUTINE fpio
1688
1689    SUBROUTINE fpmetbinary_filetext(filename, cm_index)
1690
1691        ! This is a utility subroutine meant to be used for testing purposes.
1692        ! It facilitates the text output of variables read in from the
1693        ! specified .fp file.  This routine will easily cause the program
1694        ! to crash due memory allocation issues, particularly when you are
1695        ! trying to text print 3d arrays in a single formatted statetment.
1696       
1697        CHARACTER(LEN=*), INTENT(IN) :: filename
1698        INTEGER, INTENT(IN) :: cm_index           ! Index of last dimension in
1699                                                  ! most com_mod variables.
1700                                                  ! Should be 1 or 2
1701
1702        !OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', status='REPLACE', &
1703        !    form="formatted", access="stream")
1704        OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', &
1705             form="formatted", access="APPEND")
1706
1707        WRITE(IOUNIT_TEXTOUT, *) 'oro: ', oro
1708        WRITE(IOUNIT_TEXTOUT, *) 'excessoro: ', excessoro
1709        WRITE(IOUNIT_TEXTOUT, *) 'lsm: ', lsm
1710        WRITE(IOUNIT_TEXTOUT, *) 'xlanduse: ', xlanduse
1711        WRITE(IOUNIT_TEXTOUT, *) 'height: ', height
1712
1713        WRITE(IOUNIT_TEXTOUT, *) 'uu: ', uu(:,:,:,cm_index)
1714        WRITE(IOUNIT_TEXTOUT, *) 'vv: ', vv(:,:,:,cm_index)
1715        WRITE(IOUNIT_TEXTOUT, *) 'uupol: ', uupol(:,:,:,cm_index)
1716        WRITE(IOUNIT_TEXTOUT, *) 'vvpol: ', vvpol(:,:,:,cm_index)
1717        WRITE(IOUNIT_TEXTOUT, *) 'ww: ', ww(:,:,:,cm_index)
1718        WRITE(IOUNIT_TEXTOUT, *) 'tt: ', tt(:,:,:,cm_index)
1719        WRITE(IOUNIT_TEXTOUT, *) 'qv: ', qv(:,:,:,cm_index)
1720        WRITE(IOUNIT_TEXTOUT, *) 'pv: ', pv(:,:,:,cm_index)
1721        WRITE(IOUNIT_TEXTOUT, *) 'rho: ', rho(:,:,:,cm_index)
1722        WRITE(IOUNIT_TEXTOUT, *) 'drhodz: ', drhodz(:,:,:,cm_index)
1723        WRITE(IOUNIT_TEXTOUT, *) 'tth: ', tth(:,:,:,cm_index)
1724        WRITE(IOUNIT_TEXTOUT, *) 'qvh: ', qvh(:,:,:,cm_index)
1725        WRITE(IOUNIT_TEXTOUT, *) 'pplev: ', pplev(:,:,:,cm_index)
1726        WRITE(IOUNIT_TEXTOUT, *) 'clouds: ', clouds(:,:,:,cm_index)
1727        WRITE(IOUNIT_TEXTOUT, *) 'cloudsh: ', cloudsh(:,:,cm_index)
1728
1729
1730
1731
1732        CLOSE(IOUNIT_TEXTOUT)
1733    END SUBROUTINE fpmetbinary_filetext
1734
1735
1736END MODULE fpmetbinary_mod
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG