source: flexpart.git/flexpart_code/fpmetbinary_mod.F90 @ 8413bc9

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

Continuing work on NC4 implementation of intermediate format

  • Property mode set to 100644
File size: 103.1 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&                  zero_to_nzmax_dimid
312
313
314        INTEGER, DIMENSION(1) :: dim1dids    ! Dimension IDs for 1D arrays
315        INTEGER, DIMENSION(2) :: dim2dids    ! Dimension IDs for 2D arrays
316        INTEGER, DIMENSION(3) :: dim3dids    ! Dimension IDs for 3D arrays
317        INTEGER, DIMENSION(4) :: dim4dids    ! Dimension IDs for 4D arrays
318        INTEGER, DIMENSION(5) :: dim5dids    ! Dimension IDs for 5D arrays
319
320
321
322
323        ! These are used when loading in dimensions from NC file
324        CHARACTER(LEN=NF90_MAX_NAME) :: nxmax_dimname, nymax_dimname, nzmax_dimname, &
325&                                       nuvzmax_dimname, nwzmax_dimname,&
326&                                       maxspec_dimname, numclass_dimname,&
327&                                       maxnests_dimname, nxmaxn_dimname, nymaxn_dimname, &
328&                                       zero_to_nzmax_dimname
329
330        ! These are temporary variables, used in the LOAD option, for
331        ! comparing against the current values in FLEXPART of nxmax, nymax, ...
332        INTEGER :: temp_nxmax, temp_nymax, temp_nzmax, &
333&                  temp_nuvzmax, temp_nwzmax, &
334&                  temp_maxspec, temp_numclass,&
335&                  temp_maxnests, temp_nxmaxn, temp_nymaxn
336
337        CHARACTER(LEN=12) :: temp_preproc_format_version_str
338
339        CHARACTER(LEN=128) :: errmesg
340
341        INTEGER, PARAMETER :: DEF_LEVEL = 3
342
343        if (op == 'DUMP') THEN
344
345
346            ! Write the preprocessing format version string
347            WRITE (iounit) PREPROC_FORMAT_VERSION_STR
348
349            ! Write the compiled max dimensions from par_mod - these are
350            ! not meant to be reassigned during a LOAD, but used as "header"
351            ! information to provide the structure of arrays
352            WRITE (iounit) nxmax, nymax, nzmax, nuvzmax, nwzmax
353
354            ncret = nf90_def_dim(ncid, 'nxmax', nxmax, nxmax_dimid)
355            call handle_nf90_err(ncret)
356            ncret = nf90_def_dim(ncid, 'nymax', nymax, nymax_dimid)
357            call handle_nf90_err(ncret)
358            ncret = nf90_def_dim(ncid, 'nzmax', nzmax, nzmax_dimid)
359            call handle_nf90_err(ncret)
360            ncret = nf90_def_dim(ncid, 'nuvzmax', nuvzmax, nuvzmax_dimid)
361            call handle_nf90_err(ncret)
362            ncret = nf90_def_dim(ncid, 'nwzmax', nwzmax, nwzmax_dimid)
363            call handle_nf90_err(ncret)
364            ncret = nf90_def_dim(ncid, 'maxspec', maxspec, maxspec_dimid)
365            call handle_nf90_err(ncret)
366            ncret = nf90_def_dim(ncid, 'numclass', numclass, numclass_dimid)
367            call handle_nf90_err(ncret)
368            ncret = nf90_def_dim(ncid, 'zero_to_nzmax', nzmax+1, zero_to_nzmax_dimid)
369            call handle_nf90_err(ncret)
370
371
372            ! Scalar values
373            WRITE(iounit) nx, ny, nxmin1, nymin1, nxfield
374            WRITE(iounit) nuvz, nwz, nz, nmixz, nlev_ec
375            WRITE(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst
376
377            ncret = nf90_def_var(ncid, 'nx', NF90_INT, ncvarid)
378            call handle_nf90_err(ncret)
379            ncret = nf90_put_var(ncid, ncvarid, nx)
380            call handle_nf90_err(ncret)
381
382            ncret = nf90_def_var(ncid, 'ny', NF90_INT, ncvarid)
383            call handle_nf90_err(ncret)
384            ncret = nf90_put_var(ncid, ncvarid, ny)
385            call handle_nf90_err(ncret)
386
387            ncret = nf90_def_var(ncid, 'nxmin1', NF90_INT, ncvarid)
388            call handle_nf90_err(ncret)
389            ncret = nf90_put_var(ncid, ncvarid, nxmin1)
390            call handle_nf90_err(ncret)
391
392            ncret = nf90_def_var(ncid, 'nymin1', NF90_INT, ncvarid)
393            call handle_nf90_err(ncret)
394            ncret = nf90_put_var(ncid, ncvarid, nymin1)
395            call handle_nf90_err(ncret)
396
397            ncret = nf90_def_var(ncid, 'nxfield', NF90_INT, ncvarid)
398            call handle_nf90_err(ncret)
399            ncret = nf90_put_var(ncid, ncvarid, nxfield)
400            call handle_nf90_err(ncret)
401
402            ncret = nf90_def_var(ncid, 'nuvz', NF90_INT, ncvarid)
403            call handle_nf90_err(ncret)
404            ncret = nf90_put_var(ncid, ncvarid, nuvz)
405            call handle_nf90_err(ncret)
406
407            ncret = nf90_def_var(ncid, 'nwz', NF90_INT, ncvarid)
408            call handle_nf90_err(ncret)
409            ncret = nf90_put_var(ncid, ncvarid, nwz)
410            call handle_nf90_err(ncret)
411
412            ncret = nf90_def_var(ncid, 'nz', NF90_INT, ncvarid)
413            call handle_nf90_err(ncret)
414            ncret = nf90_put_var(ncid, ncvarid, nz)
415            call handle_nf90_err(ncret)
416
417            ncret = nf90_def_var(ncid, 'nmixz', NF90_INT, ncvarid)
418            call handle_nf90_err(ncret)
419            ncret = nf90_put_var(ncid, ncvarid, nmixz)
420            call handle_nf90_err(ncret)
421
422            ncret = nf90_def_var(ncid, 'nlev_ec', NF90_INT, ncvarid)
423            call handle_nf90_err(ncret)
424            ncret = nf90_put_var(ncid, ncvarid, nlev_ec)
425            call handle_nf90_err(ncret)
426
427            ncret = nf90_def_var(ncid, 'dx', NF90_FLOAT, ncvarid)
428            call handle_nf90_err(ncret)
429            ncret = nf90_put_var(ncid, ncvarid, dx)
430            call handle_nf90_err(ncret)
431
432            ncret = nf90_def_var(ncid, 'dy', NF90_FLOAT, ncvarid)
433            call handle_nf90_err(ncret)
434            ncret = nf90_put_var(ncid, ncvarid, dy)
435            call handle_nf90_err(ncret)
436
437            ncret = nf90_def_var(ncid, 'xlon0', NF90_FLOAT, ncvarid)
438            call handle_nf90_err(ncret)
439            ncret = nf90_put_var(ncid, ncvarid, xlon0)
440            call handle_nf90_err(ncret)
441
442            ncret = nf90_def_var(ncid, 'ylat0', NF90_FLOAT, ncvarid)
443            call handle_nf90_err(ncret)
444            ncret = nf90_put_var(ncid, ncvarid, ylat0)
445            call handle_nf90_err(ncret)
446
447            ncret = nf90_def_var(ncid, 'dxconst', NF90_FLOAT, ncvarid)
448            call handle_nf90_err(ncret)
449            ncret = nf90_put_var(ncid, ncvarid, dxconst)
450            call handle_nf90_err(ncret)
451
452            ncret = nf90_def_var(ncid, 'dyconst', NF90_FLOAT, ncvarid)
453            call handle_nf90_err(ncret)
454            ncret = nf90_put_var(ncid, ncvarid, dyconst)
455            call handle_nf90_err(ncret)
456
457
458
459            ! Fixed fields, static in time
460            WRITE(iounit) oro, excessoro, lsm, xlanduse, height
461
462            dim2dids = (/nxmax_dimid, nymax_dimid/)
463
464            ncret = nf90_def_var(ncid, 'oro', NF90_FLOAT, &
465&                                       dim2dids, ncvarid)
466            call handle_nf90_err(ncret)
467            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
468&                                        shuffle=0,     &
469&                                        deflate=1,     &
470&                                        deflate_level=DEF_LEVEL)
471            call handle_nf90_err(ncret)
472            ncret = nf90_put_var(ncid, ncvarid, &
473&                                oro(0:nxmax-1, 0:nymax-1))
474            call handle_nf90_err(ncret)
475
476            ncret = nf90_def_var(ncid, 'excessoro', NF90_FLOAT, &
477&                                       dim2dids, ncvarid)
478            call handle_nf90_err(ncret)
479            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
480&                                        shuffle=0,     &
481&                                        deflate=1,     &
482&                                        deflate_level=DEF_LEVEL)
483            call handle_nf90_err(ncret)
484            ncret = nf90_put_var(ncid, ncvarid, &
485&                                excessoro(0:nxmax-1, 0:nymax-1))
486            call handle_nf90_err(ncret)
487
488            ncret = nf90_def_var(ncid, 'lsm', NF90_FLOAT, &
489&                                       dim2dids, ncvarid)
490            call handle_nf90_err(ncret)
491            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
492&                                        shuffle=0,     &
493&                                        deflate=1,     &
494&                                        deflate_level=DEF_LEVEL)
495            call handle_nf90_err(ncret)
496            ncret = nf90_put_var(ncid, ncvarid, &
497&                                lsm(0:nxmax-1, 0:nymax-1))
498            call handle_nf90_err(ncret)
499
500            dim3dids = (/nxmax_dimid, nymax_dimid, numclass_dimid/)
501            ! numclass comes from par_mod - number of land use classes
502            ncret = nf90_def_var(ncid, 'xlanduse', NF90_FLOAT, &
503&                                       dim3dids, ncvarid)
504            call handle_nf90_err(ncret)
505            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
506&                                        shuffle=0,     &
507&                                        deflate=1,     &
508&                                        deflate_level=DEF_LEVEL)
509            call handle_nf90_err(ncret)
510            ncret = nf90_put_var(ncid, ncvarid, &
511&                                xlanduse(0:nxmax-1, 0:nymax-1, 1:numclass))
512            call handle_nf90_err(ncret)
513
514            dim1dids = (/nzmax_dimid/)
515            ncret = nf90_def_var(ncid, 'height', NF90_FLOAT, &
516&                                       dim1dids, ncvarid)
517            call handle_nf90_err(ncret)
518            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
519&                                        shuffle=0,     &
520&                                        deflate=1,     &
521&                                        deflate_level=DEF_LEVEL)
522            call handle_nf90_err(ncret)
523            ncret = nf90_put_var(ncid, ncvarid, &
524&                                height(1:nzmax))
525            call handle_nf90_err(ncret)
526
527
528
529
530            ! 3d fields
531            WRITE(iounit) uu(:,:,:,cm_index)
532            WRITE(iounit) vv(:,:,:,cm_index)
533            WRITE(iounit) uupol(:,:,:,cm_index)
534            WRITE(iounit) vvpol(:,:,:,cm_index)
535            WRITE(iounit) ww(:,:,:,cm_index)
536            WRITE(iounit) tt(:,:,:,cm_index)
537            WRITE(iounit) qv(:,:,:,cm_index)
538            WRITE(iounit) pv(:,:,:,cm_index)
539            WRITE(iounit) rho(:,:,:,cm_index)
540            WRITE(iounit) drhodz(:,:,:,cm_index)
541            WRITE(iounit) tth(:,:,:,cm_index)
542            WRITE(iounit) qvh(:,:,:,cm_index)
543            WRITE(iounit) pplev(:,:,:,cm_index)
544            WRITE(iounit) clouds(:,:,:,cm_index)
545            WRITE(iounit) cloudsh(:,:,cm_index)
546
547            dim3dids = (/nxmax_dimid, nymax_dimid, nzmax_dimid/)
548
549            ncret = nf90_def_var(ncid, 'uu', NF90_FLOAT, &
550&                                       dim3dids, ncvarid)
551            call handle_nf90_err(ncret)
552            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
553&                                        shuffle=0,     &
554&                                        deflate=1,     &
555&                                        deflate_level=DEF_LEVEL)
556            call handle_nf90_err(ncret)
557            ncret = nf90_put_var(ncid, ncvarid, &
558&                                uu(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
559            call handle_nf90_err(ncret)
560
561            ncret = nf90_def_var(ncid, 'vv', NF90_FLOAT, &
562&                                       dim3dids, ncvarid)
563            call handle_nf90_err(ncret)
564            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
565&                                        shuffle=0,     &
566&                                        deflate=1,     &
567&                                        deflate_level=DEF_LEVEL)
568            call handle_nf90_err(ncret)
569            ncret = nf90_put_var(ncid, ncvarid, &
570&                                vv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
571            call handle_nf90_err(ncret)
572
573            ncret = nf90_def_var(ncid, 'uupol', NF90_FLOAT, &
574&                                       dim3dids, ncvarid)
575            call handle_nf90_err(ncret)
576            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
577&                                        shuffle=0,     &
578&                                        deflate=1,     &
579&                                        deflate_level=DEF_LEVEL)
580            call handle_nf90_err(ncret)
581            ncret = nf90_put_var(ncid, ncvarid, &
582&                                uupol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
583            call handle_nf90_err(ncret)
584
585            ncret = nf90_def_var(ncid, 'vvpol', NF90_FLOAT, &
586&                                       dim3dids, ncvarid)
587            call handle_nf90_err(ncret)
588            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
589&                                        shuffle=0,     &
590&                                        deflate=1,     &
591&                                        deflate_level=DEF_LEVEL)
592            call handle_nf90_err(ncret)
593            ncret = nf90_put_var(ncid, ncvarid, &
594&                                vvpol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
595            call handle_nf90_err(ncret)
596
597            ncret = nf90_def_var(ncid, 'ww', NF90_FLOAT, &
598&                                       dim3dids, ncvarid)
599            call handle_nf90_err(ncret)
600            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
601&                                        shuffle=0,     &
602&                                        deflate=1,     &
603&                                        deflate_level=DEF_LEVEL)
604            call handle_nf90_err(ncret)
605            ncret = nf90_put_var(ncid, ncvarid, &
606&                                ww(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
607            call handle_nf90_err(ncret)
608
609            ncret = nf90_def_var(ncid, 'tt', NF90_FLOAT, &
610&                                       dim3dids, ncvarid)
611            call handle_nf90_err(ncret)
612            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
613&                                        shuffle=0,     &
614&                                        deflate=1,     &
615&                                        deflate_level=DEF_LEVEL)
616            call handle_nf90_err(ncret)
617            ncret = nf90_put_var(ncid, ncvarid, &
618&                                tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
619            call handle_nf90_err(ncret)
620
621            ncret = nf90_def_var(ncid, 'qv', NF90_FLOAT, &
622&                                       dim3dids, ncvarid)
623            call handle_nf90_err(ncret)
624            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
625&                                        shuffle=0,     &
626&                                        deflate=1,     &
627&                                        deflate_level=DEF_LEVEL)
628            call handle_nf90_err(ncret)
629            ncret = nf90_put_var(ncid, ncvarid, &
630&                                qv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
631            call handle_nf90_err(ncret)
632
633            ncret = nf90_def_var(ncid, 'pv', NF90_FLOAT, &
634&                                       dim3dids, ncvarid)
635            call handle_nf90_err(ncret)
636            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
637&                                        shuffle=0,     &
638&                                        deflate=1,     &
639&                                        deflate_level=DEF_LEVEL)
640            call handle_nf90_err(ncret)
641            ncret = nf90_put_var(ncid, ncvarid, &
642&                                pv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
643            call handle_nf90_err(ncret)
644
645            ncret = nf90_def_var(ncid, 'rho', NF90_FLOAT, &
646&                                       dim3dids, ncvarid)
647            call handle_nf90_err(ncret)
648            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
649&                                        shuffle=0,     &
650&                                        deflate=1,     &
651&                                        deflate_level=DEF_LEVEL)
652            call handle_nf90_err(ncret)
653            ncret = nf90_put_var(ncid, ncvarid, &
654&                                rho(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
655            call handle_nf90_err(ncret)
656
657            ncret = nf90_def_var(ncid, 'drhodz', NF90_FLOAT, &
658&                                       dim3dids, ncvarid)
659            call handle_nf90_err(ncret)
660            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
661&                                        shuffle=0,     &
662&                                        deflate=1,     &
663&                                        deflate_level=DEF_LEVEL)
664            call handle_nf90_err(ncret)
665            ncret = nf90_put_var(ncid, ncvarid, &
666&                                drhodz(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
667            call handle_nf90_err(ncret)
668
669            ncret = nf90_def_var(ncid, 'clouds', NF90_BYTE, &
670&                                       dim3dids, ncvarid)
671            call handle_nf90_err(ncret)
672            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
673&                                        shuffle=0,     &
674&                                        deflate=1,     &
675&                                        deflate_level=DEF_LEVEL)
676            call handle_nf90_err(ncret)
677            ncret = nf90_put_var(ncid, ncvarid, &
678&                                clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
679            call handle_nf90_err(ncret)
680
681
682
683            ! Note the change in z dimension for the following
684            dim3dids = (/nxmax_dimid, nymax_dimid, nuvzmax_dimid/)
685
686            ncret = nf90_def_var(ncid, 'tth', NF90_FLOAT, &
687&                                       dim3dids, ncvarid)
688            call handle_nf90_err(ncret)
689            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
690&                                        shuffle=0,     &
691&                                        deflate=1,     &
692&                                        deflate_level=DEF_LEVEL)
693            call handle_nf90_err(ncret)
694            ncret = nf90_put_var(ncid, ncvarid, &
695&                                tth(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
696            call handle_nf90_err(ncret)
697
698            ncret = nf90_def_var(ncid, 'qvh', NF90_FLOAT, &
699&                                       dim3dids, ncvarid)
700            call handle_nf90_err(ncret)
701            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
702&                                        shuffle=0,     &
703&                                        deflate=1,     &
704&                                        deflate_level=DEF_LEVEL)
705            call handle_nf90_err(ncret)
706            ncret = nf90_put_var(ncid, ncvarid, &
707&                                qvh(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
708            call handle_nf90_err(ncret)
709
710            ncret = nf90_def_var(ncid, 'pplev', NF90_FLOAT, &
711&                                       dim3dids, ncvarid)
712            call handle_nf90_err(ncret)
713            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
714&                                        shuffle=0,     &
715&                                        deflate=1,     &
716&                                        deflate_level=DEF_LEVEL)
717            call handle_nf90_err(ncret)
718            ncret = nf90_put_var(ncid, ncvarid, &
719&                                pplev(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
720            call handle_nf90_err(ncret)
721
722
723            dim2dids = (/nxmax_dimid, nymax_dimid/)
724            ncret = nf90_def_var(ncid, 'cloudsh', NF90_INT, &
725&                                       dim2dids, ncvarid)
726            call handle_nf90_err(ncret)
727            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
728&                                        shuffle=0,     &
729&                                        deflate=1,     &
730&                                        deflate_level=DEF_LEVEL)
731            call handle_nf90_err(ncret)
732            ncret = nf90_put_var(ncid, ncvarid, &
733&                                cloudsh(0:nxmax-1, 0:nymax-1, cm_index))
734            call handle_nf90_err(ncret)
735
736
737
738            PRINT *, 'SUM(tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
739&                                        SUM(tt(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
740
741            PRINT *, 'SUM(clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
742&                                        SUM(clouds(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
743
744
745
746            ! 2d fields
747            WRITE(iounit) ps(:,:,:,cm_index)
748            WRITE(iounit) sd(:,:,:,cm_index)
749            WRITE(iounit) msl(:,:,:,cm_index)
750            WRITE(iounit) tcc(:,:,:,cm_index)
751            WRITE(iounit) u10(:,:,:,cm_index)
752            WRITE(iounit) v10(:,:,:,cm_index)
753            WRITE(iounit) tt2(:,:,:,cm_index)
754            WRITE(iounit) td2(:,:,:,cm_index)
755            WRITE(iounit) lsprec(:,:,:,cm_index)
756            WRITE(iounit) convprec(:,:,:,cm_index)
757            WRITE(iounit) sshf(:,:,:,cm_index)
758            WRITE(iounit) ssr(:,:,:,cm_index)
759            WRITE(iounit) surfstr(:,:,:,cm_index)
760            WRITE(iounit) ustar(:,:,:,cm_index)
761            WRITE(iounit) wstar(:,:,:,cm_index)
762            WRITE(iounit) hmix(:,:,:,cm_index)
763            WRITE(iounit) tropopause(:,:,:,cm_index)
764            WRITE(iounit) oli(:,:,:,cm_index)
765            WRITE(iounit) diffk(:,:,:,cm_index)
766            WRITE(iounit) vdep(:,:,:,cm_index)
767
768            dim2dids = (/nxmax_dimid, nymax_dimid/)
769
770            ncret = nf90_def_var(ncid, 'ps', NF90_FLOAT, &
771&                                       dim2dids, ncvarid)
772            call handle_nf90_err(ncret)
773            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
774&                                        shuffle=0,     &
775&                                        deflate=1,     &
776&                                        deflate_level=DEF_LEVEL)
777            call handle_nf90_err(ncret)
778            ncret = nf90_put_var(ncid, ncvarid, &
779&                                ps(0:nxmax-1, 0:nymax-1, 1, cm_index))
780            call handle_nf90_err(ncret)
781
782            ncret = nf90_def_var(ncid, 'sd', NF90_FLOAT, &
783&                                       dim2dids, ncvarid)
784            call handle_nf90_err(ncret)
785            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
786&                                        shuffle=0,     &
787&                                        deflate=1,     &
788&                                        deflate_level=DEF_LEVEL)
789            call handle_nf90_err(ncret)
790            ncret = nf90_put_var(ncid, ncvarid, &
791&                                sd(0:nxmax-1, 0:nymax-1, 1, cm_index))
792            call handle_nf90_err(ncret)
793
794            ncret = nf90_def_var(ncid, 'msl', NF90_FLOAT, &
795&                                       dim2dids, ncvarid)
796            call handle_nf90_err(ncret)
797            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
798&                                        shuffle=0,     &
799&                                        deflate=1,     &
800&                                        deflate_level=DEF_LEVEL)
801            call handle_nf90_err(ncret)
802            ncret = nf90_put_var(ncid, ncvarid, &
803&                                msl(0:nxmax-1, 0:nymax-1, 1, cm_index))
804            call handle_nf90_err(ncret)
805
806            ncret = nf90_def_var(ncid, 'tcc', NF90_FLOAT, &
807&                                       dim2dids, ncvarid)
808            call handle_nf90_err(ncret)
809            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
810&                                        shuffle=0,     &
811&                                        deflate=1,     &
812&                                        deflate_level=DEF_LEVEL)
813            call handle_nf90_err(ncret)
814            ncret = nf90_put_var(ncid, ncvarid, &
815&                                tcc(0:nxmax-1, 0:nymax-1, 1, cm_index))
816            call handle_nf90_err(ncret)
817
818            ncret = nf90_def_var(ncid, 'u10', NF90_FLOAT, &
819&                                       dim2dids, ncvarid)
820            call handle_nf90_err(ncret)
821            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
822&                                        shuffle=0,     &
823&                                        deflate=1,     &
824&                                        deflate_level=DEF_LEVEL)
825            call handle_nf90_err(ncret)
826            ncret = nf90_put_var(ncid, ncvarid, &
827&                                u10(0:nxmax-1, 0:nymax-1, 1, cm_index))
828            call handle_nf90_err(ncret)
829
830            ncret = nf90_def_var(ncid, 'v10', NF90_FLOAT, &
831&                                       dim2dids, ncvarid)
832            call handle_nf90_err(ncret)
833            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
834&                                        shuffle=0,     &
835&                                        deflate=1,     &
836&                                        deflate_level=DEF_LEVEL)
837            call handle_nf90_err(ncret)
838            ncret = nf90_put_var(ncid, ncvarid, &
839&                                v10(0:nxmax-1, 0:nymax-1, 1, cm_index))
840            call handle_nf90_err(ncret)
841
842            ncret = nf90_def_var(ncid, 'tt2', NF90_FLOAT, &
843&                                       dim2dids, ncvarid)
844            call handle_nf90_err(ncret)
845            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
846&                                        shuffle=0,     &
847&                                        deflate=1,     &
848&                                        deflate_level=DEF_LEVEL)
849            call handle_nf90_err(ncret)
850            ncret = nf90_put_var(ncid, ncvarid, &
851&                                tt2(0:nxmax-1, 0:nymax-1, 1, cm_index))
852            call handle_nf90_err(ncret)
853
854            ncret = nf90_def_var(ncid, 'td2', NF90_FLOAT, &
855&                                       dim2dids, ncvarid)
856            call handle_nf90_err(ncret)
857            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
858&                                        shuffle=0,     &
859&                                        deflate=1,     &
860&                                        deflate_level=DEF_LEVEL)
861            call handle_nf90_err(ncret)
862            ncret = nf90_put_var(ncid, ncvarid, &
863&                                td2(0:nxmax-1, 0:nymax-1, 1, cm_index))
864            call handle_nf90_err(ncret)
865
866            ncret = nf90_def_var(ncid, 'lsprec', NF90_FLOAT, &
867&                                       dim2dids, ncvarid)
868            call handle_nf90_err(ncret)
869            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
870&                                        shuffle=0,     &
871&                                        deflate=1,     &
872&                                        deflate_level=DEF_LEVEL)
873            call handle_nf90_err(ncret)
874            ncret = nf90_put_var(ncid, ncvarid, &
875&                                lsprec(0:nxmax-1, 0:nymax-1, 1, cm_index))
876            call handle_nf90_err(ncret)
877
878            ncret = nf90_def_var(ncid, 'convprec', NF90_FLOAT, &
879&                                       dim2dids, ncvarid)
880            call handle_nf90_err(ncret)
881            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
882&                                        shuffle=0,     &
883&                                        deflate=1,     &
884&                                        deflate_level=DEF_LEVEL)
885            call handle_nf90_err(ncret)
886            ncret = nf90_put_var(ncid, ncvarid, &
887&                                convprec(0:nxmax-1, 0:nymax-1, 1, cm_index))
888            call handle_nf90_err(ncret)
889
890            ncret = nf90_def_var(ncid, 'sshf', NF90_FLOAT, &
891&                                       dim2dids, ncvarid)
892            call handle_nf90_err(ncret)
893            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
894&                                        shuffle=0,     &
895&                                        deflate=1,     &
896&                                        deflate_level=DEF_LEVEL)
897            call handle_nf90_err(ncret)
898            ncret = nf90_put_var(ncid, ncvarid, &
899&                                sshf(0:nxmax-1, 0:nymax-1, 1, cm_index))
900            call handle_nf90_err(ncret)
901
902            ncret = nf90_def_var(ncid, 'ssr', NF90_FLOAT, &
903&                                       dim2dids, ncvarid)
904            call handle_nf90_err(ncret)
905            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
906&                                        shuffle=0,     &
907&                                        deflate=1,     &
908&                                        deflate_level=DEF_LEVEL)
909            call handle_nf90_err(ncret)
910            ncret = nf90_put_var(ncid, ncvarid, &
911&                                ssr(0:nxmax-1, 0:nymax-1, 1, cm_index))
912            call handle_nf90_err(ncret)
913
914            ncret = nf90_def_var(ncid, 'surfstr', NF90_FLOAT, &
915&                                       dim2dids, ncvarid)
916            call handle_nf90_err(ncret)
917            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
918&                                        shuffle=0,     &
919&                                        deflate=1,     &
920&                                        deflate_level=DEF_LEVEL)
921            call handle_nf90_err(ncret)
922            ncret = nf90_put_var(ncid, ncvarid, &
923&                                surfstr(0:nxmax-1, 0:nymax-1, 1, cm_index))
924            call handle_nf90_err(ncret)
925
926            ncret = nf90_def_var(ncid, 'ustar', NF90_FLOAT, &
927&                                       dim2dids, ncvarid)
928            call handle_nf90_err(ncret)
929            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
930&                                        shuffle=0,     &
931&                                        deflate=1,     &
932&                                        deflate_level=DEF_LEVEL)
933            call handle_nf90_err(ncret)
934            ncret = nf90_put_var(ncid, ncvarid, &
935&                                ustar(0:nxmax-1, 0:nymax-1, 1, cm_index))
936            call handle_nf90_err(ncret)
937
938            ncret = nf90_def_var(ncid, 'wstar', NF90_FLOAT, &
939&                                       dim2dids, ncvarid)
940            call handle_nf90_err(ncret)
941            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
942&                                        shuffle=0,     &
943&                                        deflate=1,     &
944&                                        deflate_level=DEF_LEVEL)
945            call handle_nf90_err(ncret)
946            ncret = nf90_put_var(ncid, ncvarid, &
947&                                wstar(0:nxmax-1, 0:nymax-1, 1, cm_index))
948            call handle_nf90_err(ncret)
949
950            ncret = nf90_def_var(ncid, 'hmix', NF90_FLOAT, &
951&                                       dim2dids, ncvarid)
952            call handle_nf90_err(ncret)
953            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
954&                                        shuffle=0,     &
955&                                        deflate=1,     &
956&                                        deflate_level=DEF_LEVEL)
957            call handle_nf90_err(ncret)
958            ncret = nf90_put_var(ncid, ncvarid, &
959&                                hmix(0:nxmax-1, 0:nymax-1, 1, cm_index))
960            call handle_nf90_err(ncret)
961
962            ncret = nf90_def_var(ncid, 'tropopause', NF90_FLOAT, &
963&                                       dim2dids, ncvarid)
964            call handle_nf90_err(ncret)
965            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
966&                                        shuffle=0,     &
967&                                        deflate=1,     &
968&                                        deflate_level=DEF_LEVEL)
969            call handle_nf90_err(ncret)
970            ncret = nf90_put_var(ncid, ncvarid, &
971&                                tropopause(0:nxmax-1, 0:nymax-1, 1, cm_index))
972            call handle_nf90_err(ncret)
973
974            ncret = nf90_def_var(ncid, 'oli', NF90_FLOAT, &
975&                                       dim2dids, ncvarid)
976            call handle_nf90_err(ncret)
977            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
978&                                        shuffle=0,     &
979&                                        deflate=1,     &
980&                                        deflate_level=DEF_LEVEL)
981            call handle_nf90_err(ncret)
982            ncret = nf90_put_var(ncid, ncvarid, &
983&                                oli(0:nxmax-1, 0:nymax-1, 1, cm_index))
984            call handle_nf90_err(ncret)
985
986            ncret = nf90_def_var(ncid, 'diffk', NF90_FLOAT, &
987&                                       dim2dids, ncvarid)
988            call handle_nf90_err(ncret)
989            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
990&                                        shuffle=0,     &
991&                                        deflate=1,     &
992&                                        deflate_level=DEF_LEVEL)
993            call handle_nf90_err(ncret)
994            ncret = nf90_put_var(ncid, ncvarid, &
995&                                diffk(0:nxmax-1, 0:nymax-1, 1, cm_index))
996            call handle_nf90_err(ncret)
997
998
999
1000            PRINT *, 'SUM(ps(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
1001&                                        SUM(ps(0:nxmax-1,0:nymax-1,1, cm_index))
1002
1003            PRINT *, 'SUM(wstar(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
1004&                                        SUM(wstar(0:nxmax-1,0:nymax-1,1, cm_index))
1005
1006
1007            dim3dids = (/nxmax_dimid, nymax_dimid, maxspec_dimid/)
1008
1009            ncret = nf90_def_var(ncid, 'vdep', NF90_FLOAT, &
1010&                                       dim3dids, ncvarid)
1011            call handle_nf90_err(ncret)
1012            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1013&                                        shuffle=0,     &
1014&                                        deflate=1,     &
1015&                                        deflate_level=DEF_LEVEL)
1016            call handle_nf90_err(ncret)
1017            ncret = nf90_put_var(ncid, ncvarid, &
1018&                                vdep(0:nxmax-1, 0:nymax-1, 1:maxspec, cm_index))
1019            call handle_nf90_err(ncret)
1020
1021
1022
1023            ! 1d fields
1024            WRITE(iounit) z0(:)
1025            WRITE(iounit) akm(:)
1026            WRITE(iounit) bkm(:)
1027            WRITE(iounit) akz(:)
1028            WRITE(iounit) bkz(:)
1029            WRITE(iounit) aknew(:)
1030            WRITE(iounit) bknew(:)
1031
1032            dim1dids = (/numclass_dimid/)
1033
1034            ncret = nf90_def_var(ncid, 'z0', NF90_FLOAT, &
1035&                                       dim1dids, ncvarid)
1036            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1037&                                        shuffle=0,     &
1038&                                        deflate=1,     &
1039&                                        deflate_level=DEF_LEVEL)
1040            ncret = nf90_put_var(ncid, ncvarid, &
1041&                                z0(1:numclass))
1042
1043
1044            dim1dids = (/nwzmax_dimid/)
1045
1046            ncret = nf90_def_var(ncid, 'akm', NF90_FLOAT, &
1047&                                       dim1dids, ncvarid)
1048            call handle_nf90_err(ncret)
1049            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1050&                                        shuffle=0,     &
1051&                                        deflate=1,     &
1052&                                        deflate_level=DEF_LEVEL)
1053            call handle_nf90_err(ncret)
1054            ncret = nf90_put_var(ncid, ncvarid, &
1055&                                akm(1:nwzmax))
1056            call handle_nf90_err(ncret)
1057
1058            ncret = nf90_def_var(ncid, 'bkm', NF90_FLOAT, &
1059&                                       dim1dids, ncvarid)
1060            call handle_nf90_err(ncret)
1061            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1062&                                        shuffle=0,     &
1063&                                        deflate=1,     &
1064&                                        deflate_level=DEF_LEVEL)
1065            call handle_nf90_err(ncret)
1066            ncret = nf90_put_var(ncid, ncvarid, &
1067&                                bkm(1:nwzmax))
1068            call handle_nf90_err(ncret)
1069
1070
1071            dim1dids = (/nuvzmax_dimid/)
1072
1073            ncret = nf90_def_var(ncid, 'akz', NF90_FLOAT, &
1074&                                       dim1dids, ncvarid)
1075            call handle_nf90_err(ncret)
1076            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1077&                                        shuffle=0,     &
1078&                                        deflate=1,     &
1079&                                        deflate_level=DEF_LEVEL)
1080            call handle_nf90_err(ncret)
1081            ncret = nf90_put_var(ncid, ncvarid, &
1082&                                akz(1:nuvzmax))
1083            call handle_nf90_err(ncret)
1084
1085            ncret = nf90_def_var(ncid, 'bkz', NF90_FLOAT, &
1086&                                       dim1dids, ncvarid)
1087            call handle_nf90_err(ncret)
1088            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1089&                                        shuffle=0,     &
1090&                                        deflate=1,     &
1091&                                        deflate_level=DEF_LEVEL)
1092            call handle_nf90_err(ncret)
1093            ncret = nf90_put_var(ncid, ncvarid, &
1094&                                bkz(1:nuvzmax))
1095            call handle_nf90_err(ncret)
1096
1097
1098            dim1dids = (/nzmax_dimid/)
1099
1100            ncret = nf90_def_var(ncid, 'aknew', NF90_FLOAT, &
1101&                                       dim1dids, ncvarid)
1102            call handle_nf90_err(ncret)
1103            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1104&                                        shuffle=0,     &
1105&                                        deflate=1,     &
1106&                                        deflate_level=DEF_LEVEL)
1107            call handle_nf90_err(ncret)
1108            ncret = nf90_put_var(ncid, ncvarid, &
1109&                                aknew(1:nzmax))
1110            call handle_nf90_err(ncret)
1111
1112            ncret = nf90_def_var(ncid, 'bknew', NF90_FLOAT, &
1113&                                       dim1dids, ncvarid)
1114            call handle_nf90_err(ncret)
1115            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1116&                                        shuffle=0,     &
1117&                                        deflate=1,     &
1118&                                        deflate_level=DEF_LEVEL)
1119            call handle_nf90_err(ncret)
1120            ncret = nf90_put_var(ncid, ncvarid, &
1121&                                bknew(1:nzmax))
1122            call handle_nf90_err(ncret)
1123
1124
1125            PRINT *, 'SUM(bknew(1:nzmax)): ', &
1126&                                        SUM(bknew(1:nzmax))
1127
1128
1129
1130            ! Getting ready to add in nested code
1131
1132            ! These are compiled max dimensions from par_mod - these are
1133            ! not meant to be reassigned during a LOAD, but used as "header"
1134            ! information to provide the structure of arrays
1135            ncret = nf90_def_dim(ncid, 'maxnests', maxnests, maxnests_dimid)
1136            call handle_nf90_err(ncret)
1137            ncret = nf90_def_dim(ncid, 'nxmaxn', nxmaxn, nxmaxn_dimid)
1138            call handle_nf90_err(ncret)
1139            ncret = nf90_def_dim(ncid, 'nymaxn', nymaxn, nymaxn_dimid)
1140            call handle_nf90_err(ncret)
1141
1142            WRITE(iounit) nxn(:)
1143            WRITE(iounit) nyn(:)
1144            WRITE(iounit) dxn(:)
1145            WRITE(iounit) dyn(:)
1146            WRITE(iounit) xlon0n(:)
1147            WRITE(iounit) ylat0n(:)
1148
1149            ! Nested, scalar values (for each nest)
1150
1151            dim1dids = (/maxnests_dimid/)
1152
1153            ncret = nf90_def_var(ncid, 'nxn', NF90_INT, &
1154&                                       dim1dids, ncvarid)
1155            call handle_nf90_err(ncret)
1156            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1157&                                        shuffle=0,     &
1158&                                        deflate=1,     &
1159&                                        deflate_level=DEF_LEVEL)
1160            call handle_nf90_err(ncret)
1161            ncret = nf90_put_var(ncid, ncvarid, &
1162&                                nxn(1:maxnests))
1163            call handle_nf90_err(ncret)
1164
1165            ncret = nf90_def_var(ncid, 'nyn', NF90_INT, &
1166&                                       dim1dids, ncvarid)
1167            call handle_nf90_err(ncret)
1168            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1169&                                        shuffle=0,     &
1170&                                        deflate=1,     &
1171&                                        deflate_level=DEF_LEVEL)
1172            call handle_nf90_err(ncret)
1173            ncret = nf90_put_var(ncid, ncvarid, &
1174&                                nyn(1:maxnests))
1175            call handle_nf90_err(ncret)
1176
1177            ncret = nf90_def_var(ncid, 'dxn', NF90_FLOAT, &
1178&                                       dim1dids, ncvarid)
1179            call handle_nf90_err(ncret)
1180            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1181&                                        shuffle=0,     &
1182&                                        deflate=1,     &
1183&                                        deflate_level=DEF_LEVEL)
1184            call handle_nf90_err(ncret)
1185            ncret = nf90_put_var(ncid, ncvarid, &
1186&                                dxn(1:maxnests))
1187            call handle_nf90_err(ncret)
1188
1189            ncret = nf90_def_var(ncid, 'dyn', NF90_FLOAT, &
1190&                                       dim1dids, ncvarid)
1191            call handle_nf90_err(ncret)
1192            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1193&                                        shuffle=0,     &
1194&                                        deflate=1,     &
1195&                                        deflate_level=DEF_LEVEL)
1196            call handle_nf90_err(ncret)
1197            ncret = nf90_put_var(ncid, ncvarid, &
1198&                                dyn(1:maxnests))
1199            call handle_nf90_err(ncret)
1200
1201            ncret = nf90_def_var(ncid, 'xlon0n', NF90_FLOAT, &
1202&                                       dim1dids, ncvarid)
1203            call handle_nf90_err(ncret)
1204            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1205&                                        shuffle=0,     &
1206&                                        deflate=1,     &
1207&                                        deflate_level=DEF_LEVEL)
1208            call handle_nf90_err(ncret)
1209            ncret = nf90_put_var(ncid, ncvarid, &
1210&                                xlon0n(1:maxnests))
1211            call handle_nf90_err(ncret)
1212
1213            ncret = nf90_def_var(ncid, 'ylat0n', NF90_FLOAT, &
1214&                                       dim1dids, ncvarid)
1215            call handle_nf90_err(ncret)
1216            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1217&                                        shuffle=0,     &
1218&                                        deflate=1,     &
1219&                                        deflate_level=DEF_LEVEL)
1220            call handle_nf90_err(ncret)
1221            ncret = nf90_put_var(ncid, ncvarid, &
1222&                                ylat0n(1:maxnests))
1223            call handle_nf90_err(ncret)
1224
1225
1226
1227
1228            ! Nested fields, static over time
1229            WRITE(iounit) oron, excessoron, lsmn, xlandusen
1230
1231            dim3dids = (/nxmaxn_dimid, nymaxn_dimid, maxnests_dimid/)
1232
1233            ncret = nf90_def_var(ncid, 'oron', NF90_FLOAT, &
1234&                                       dim3dids, ncvarid)
1235            call handle_nf90_err(ncret)
1236            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1237&                                        shuffle=0,     &
1238&                                        deflate=1,     &
1239&                                        deflate_level=DEF_LEVEL)
1240            call handle_nf90_err(ncret)
1241            ncret = nf90_put_var(ncid, ncvarid, &
1242&                                oron(0:nxmaxn-1, 0:nymaxn-1, 1:maxnests))
1243            call handle_nf90_err(ncret)
1244
1245            ncret = nf90_def_var(ncid, 'excessoron', NF90_FLOAT, &
1246&                                       dim3dids, ncvarid)
1247            call handle_nf90_err(ncret)
1248            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1249&                                        shuffle=0,     &
1250&                                        deflate=1,     &
1251&                                        deflate_level=DEF_LEVEL)
1252            call handle_nf90_err(ncret)
1253            ncret = nf90_put_var(ncid, ncvarid, &
1254&                                excessoron(0:nxmaxn-1, 0:nymaxn-1, 1:maxnests))
1255            call handle_nf90_err(ncret)
1256
1257            ncret = nf90_def_var(ncid, 'lsmn', NF90_FLOAT, &
1258&                                       dim3dids, ncvarid)
1259            call handle_nf90_err(ncret)
1260            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1261&                                        shuffle=0,     &
1262&                                        deflate=1,     &
1263&                                        deflate_level=DEF_LEVEL)
1264            call handle_nf90_err(ncret)
1265            ncret = nf90_put_var(ncid, ncvarid, &
1266&                                lsmn(0:nxmaxn-1, 0:nymaxn-1, 1:maxnests))
1267            call handle_nf90_err(ncret)
1268
1269            dim4dids = (/nxmaxn_dimid, nymaxn_dimid, numclass_dimid, maxnests_dimid/)
1270
1271            ncret = nf90_def_var(ncid, 'xlandusen', NF90_FLOAT, &
1272&                                       dim4dids, ncvarid)
1273            call handle_nf90_err(ncret)
1274            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1275&                                        shuffle=0,     &
1276&                                        deflate=1,     &
1277&                                        deflate_level=DEF_LEVEL)
1278            call handle_nf90_err(ncret)
1279            ncret = nf90_put_var(ncid, ncvarid, &
1280&                                xlandusen(0:nxmaxn-1, 0:nymaxn-1, 1:numclass, 1:maxnests))
1281            call handle_nf90_err(ncret)
1282
1283            PRINT *, 'SUM(oron): ', SUM(oron)
1284
1285
1286
1287            ! 3d nested fields
1288            WRITE(iounit) uun(:,:,:,cm_index,:)
1289            WRITE(iounit) vvn(:,:,:,cm_index,:)
1290            WRITE(iounit) wwn(:,:,:,cm_index,:)
1291            WRITE(iounit) ttn(:,:,:,cm_index,:)
1292            WRITE(iounit) qvn(:,:,:,cm_index,:)
1293            WRITE(iounit) pvn(:,:,:,cm_index,:)
1294            WRITE(iounit) cloudsn(:,:,:,cm_index,:)
1295            WRITE(iounit) cloudsnh(:,:,cm_index,:)
1296            WRITE(iounit) rhon(:,:,:,cm_index,:)
1297            WRITE(iounit) drhodzn(:,:,:,cm_index,:)
1298            WRITE(iounit) tthn(:,:,:,cm_index,:)
1299            WRITE(iounit) qvhn(:,:,:,cm_index,:)
1300
1301
1302            dim4dids = (/nxmaxn_dimid, nymaxn_dimid, nzmax_dimid, maxnests_dimid/)
1303
1304            ncret = nf90_def_var(ncid, 'uun', NF90_FLOAT, &
1305&                                       dim4dids, ncvarid)
1306            call handle_nf90_err(ncret)
1307            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1308&                                        shuffle=0,     &
1309&                                        deflate=1,     &
1310&                                        deflate_level=DEF_LEVEL)
1311            call handle_nf90_err(ncret)
1312            ncret = nf90_put_var(ncid, ncvarid, &
1313&                                uun(0:nxmaxn-1, 0:nymaxn-1, 1:nzmax, cm_index, 1:maxnests))
1314            call handle_nf90_err(ncret)
1315
1316            ncret = nf90_def_var(ncid, 'vvn', NF90_FLOAT, &
1317&                                       dim4dids, ncvarid)
1318            call handle_nf90_err(ncret)
1319            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1320&                                        shuffle=0,     &
1321&                                        deflate=1,     &
1322&                                        deflate_level=DEF_LEVEL)
1323            call handle_nf90_err(ncret)
1324            ncret = nf90_put_var(ncid, ncvarid, &
1325&                                vvn(0:nxmaxn-1, 0:nymaxn-1, 1:nzmax, cm_index, 1:maxnests))
1326            call handle_nf90_err(ncret)
1327
1328            ncret = nf90_def_var(ncid, 'wwn', NF90_FLOAT, &
1329&                                       dim4dids, ncvarid)
1330            call handle_nf90_err(ncret)
1331            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1332&                                        shuffle=0,     &
1333&                                        deflate=1,     &
1334&                                        deflate_level=DEF_LEVEL)
1335            call handle_nf90_err(ncret)
1336            ncret = nf90_put_var(ncid, ncvarid, &
1337&                                wwn(0:nxmaxn-1, 0:nymaxn-1, 1:nzmax, cm_index, 1:maxnests))
1338            call handle_nf90_err(ncret)
1339
1340            ncret = nf90_def_var(ncid, 'ttn', NF90_FLOAT, &
1341&                                       dim4dids, ncvarid)
1342            call handle_nf90_err(ncret)
1343            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1344&                                        shuffle=0,     &
1345&                                        deflate=1,     &
1346&                                        deflate_level=DEF_LEVEL)
1347            call handle_nf90_err(ncret)
1348            ncret = nf90_put_var(ncid, ncvarid, &
1349&                                ttn(0:nxmaxn-1, 0:nymaxn-1, 1:nzmax, cm_index, 1:maxnests))
1350            call handle_nf90_err(ncret)
1351
1352            ncret = nf90_def_var(ncid, 'qvn', NF90_FLOAT, &
1353&                                       dim4dids, ncvarid)
1354            call handle_nf90_err(ncret)
1355            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1356&                                        shuffle=0,     &
1357&                                        deflate=1,     &
1358&                                        deflate_level=DEF_LEVEL)
1359            call handle_nf90_err(ncret)
1360            ncret = nf90_put_var(ncid, ncvarid, &
1361&                                qvn(0:nxmaxn-1, 0:nymaxn-1, 1:nzmax, cm_index, 1:maxnests))
1362            call handle_nf90_err(ncret)
1363
1364            ncret = nf90_def_var(ncid, 'pvn', NF90_FLOAT, &
1365&                                       dim4dids, ncvarid)
1366            call handle_nf90_err(ncret)
1367            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1368&                                        shuffle=0,     &
1369&                                        deflate=1,     &
1370&                                        deflate_level=DEF_LEVEL)
1371            call handle_nf90_err(ncret)
1372            ncret = nf90_put_var(ncid, ncvarid, &
1373&                                pvn(0:nxmaxn-1, 0:nymaxn-1, 1:nzmax, cm_index, 1:maxnests))
1374            call handle_nf90_err(ncret)
1375
1376            ncret = nf90_def_var(ncid, 'rhon', NF90_FLOAT, &
1377&                                       dim4dids, ncvarid)
1378            call handle_nf90_err(ncret)
1379            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1380&                                        shuffle=0,     &
1381&                                        deflate=1,     &
1382&                                        deflate_level=DEF_LEVEL)
1383            call handle_nf90_err(ncret)
1384            ncret = nf90_put_var(ncid, ncvarid, &
1385&                                rhon(0:nxmaxn-1, 0:nymaxn-1, 1:nzmax, cm_index, 1:maxnests))
1386            call handle_nf90_err(ncret)
1387
1388            ncret = nf90_def_var(ncid, 'drhodzn', NF90_FLOAT, &
1389&                                       dim4dids, ncvarid)
1390            call handle_nf90_err(ncret)
1391            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1392&                                        shuffle=0,     &
1393&                                        deflate=1,     &
1394&                                        deflate_level=DEF_LEVEL)
1395            call handle_nf90_err(ncret)
1396            ncret = nf90_put_var(ncid, ncvarid, &
1397&                                drhodzn(0:nxmaxn-1, 0:nymaxn-1, 1:nzmax, cm_index, 1:maxnests))
1398            call handle_nf90_err(ncret)
1399
1400
1401            ! Note the new dimensions
1402            dim4dids = (/nxmaxn_dimid, nymaxn_dimid, nuvzmax_dimid, maxnests_dimid/)
1403
1404            ncret = nf90_def_var(ncid, 'tthn', NF90_FLOAT, &
1405&                                       dim4dids, ncvarid)
1406            call handle_nf90_err(ncret)
1407            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1408&                                        shuffle=0,     &
1409&                                        deflate=1,     &
1410&                                        deflate_level=DEF_LEVEL)
1411            call handle_nf90_err(ncret)
1412            ncret = nf90_put_var(ncid, ncvarid, &
1413&                                tthn(0:nxmaxn-1, 0:nymaxn-1, 1:nuvzmax, cm_index, 1:maxnests))
1414            call handle_nf90_err(ncret)
1415
1416            ncret = nf90_def_var(ncid, 'qvhn', NF90_FLOAT, &
1417&                                       dim4dids, ncvarid)
1418            call handle_nf90_err(ncret)
1419            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1420&                                        shuffle=0,     &
1421&                                        deflate=1,     &
1422&                                        deflate_level=DEF_LEVEL)
1423            call handle_nf90_err(ncret)
1424            ncret = nf90_put_var(ncid, ncvarid, &
1425&                                qvhn(0:nxmaxn-1, 0:nymaxn-1, 1:nuvzmax, cm_index, 1:maxnests))
1426            call handle_nf90_err(ncret)
1427
1428            ! Note the new dimensions
1429            dim4dids = (/nxmaxn_dimid, nymaxn_dimid, zero_to_nzmax_dimid, maxnests_dimid/)
1430
1431            ncret = nf90_def_var(ncid, 'cloudsn', NF90_INT, &
1432&                                       dim4dids, ncvarid)
1433            call handle_nf90_err(ncret)
1434            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1435&                                        shuffle=0,     &
1436&                                        deflate=1,     &
1437&                                        deflate_level=DEF_LEVEL)
1438            call handle_nf90_err(ncret)
1439            ncret = nf90_put_var(ncid, ncvarid, &
1440&                                cloudsn(0:nxmaxn-1, 0:nymaxn-1, 0:nzmax, cm_index, 1:maxnests))
1441            call handle_nf90_err(ncret)
1442
1443            ! Note the new dimensions
1444            dim3dids = (/nxmaxn_dimid, nymaxn_dimid, maxnests_dimid/)
1445
1446            ncret = nf90_def_var(ncid, 'cloudsnh', NF90_INT, &
1447&                                       dim3dids, ncvarid)
1448            call handle_nf90_err(ncret)
1449            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1450&                                        shuffle=0,     &
1451&                                        deflate=1,     &
1452&                                        deflate_level=DEF_LEVEL)
1453            call handle_nf90_err(ncret)
1454            ncret = nf90_put_var(ncid, ncvarid, &
1455&                                cloudsnh(0:nxmaxn-1, 0:nymaxn-1, cm_index, 1:maxnests))
1456            call handle_nf90_err(ncret)
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466            PRINT *, 'SUM(uun): ', SUM(uun(:,:,:,cm_index,:))
1467            PRINT *, 'SUM(qvhn): ', SUM(qvhn(:,:,:,cm_index,:))
1468            PRINT *, 'SUM(cloudsn): ', SUM(cloudsn(:,:,:,cm_index,:))
1469
1470
1471
1472            ! 2d nested fields
1473            WRITE(iounit) psn(:,:,:,cm_index,:)
1474            WRITE(iounit) sdn(:,:,:,cm_index,:)
1475            WRITE(iounit) msln(:,:,:,cm_index,:)
1476            WRITE(iounit) tccn(:,:,:,cm_index,:)
1477            WRITE(iounit) u10n(:,:,:,cm_index,:)
1478            WRITE(iounit) v10n(:,:,:,cm_index,:)
1479            WRITE(iounit) tt2n(:,:,:,cm_index,:)
1480            WRITE(iounit) td2n(:,:,:,cm_index,:)
1481            WRITE(iounit) lsprecn(:,:,:,cm_index,:)
1482            WRITE(iounit) convprecn(:,:,:,cm_index,:)
1483            WRITE(iounit) sshfn(:,:,:,cm_index,:)
1484            WRITE(iounit) ssrn(:,:,:,cm_index,:)
1485            WRITE(iounit) surfstrn(:,:,:,cm_index,:)
1486            WRITE(iounit) ustarn(:,:,:,cm_index,:)
1487            WRITE(iounit) wstarn(:,:,:,cm_index,:)
1488            WRITE(iounit) hmixn(:,:,:,cm_index,:)
1489            WRITE(iounit) tropopausen(:,:,:,cm_index,:)
1490            WRITE(iounit) olin(:,:,:,cm_index,:)
1491            WRITE(iounit) diffkn(:,:,:,cm_index,:)
1492            WRITE(iounit) vdepn(:,:,:,cm_index,:)
1493
1494            dim3dids = (/nxmax_dimid, nymax_dimid, maxnests_dimid/)
1495
1496            ncret = nf90_def_var(ncid, 'psn', NF90_FLOAT, &
1497&                                       dim3dids, ncvarid)
1498            call handle_nf90_err(ncret)
1499            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
1500&                                        shuffle=0,     &
1501&                                        deflate=1,     &
1502&                                        deflate_level=DEF_LEVEL)
1503            call handle_nf90_err(ncret)
1504            ncret = nf90_put_var(ncid, ncvarid, &
1505&                                psn(0:nxmaxn-1, 0:nymaxn-1, 1, cm_index, 1:maxnests))
1506            call handle_nf90_err(ncret)
1507
1508
1509
1510
1511
1512
1513            PRINT *, 'SUM(psn): ', SUM(psn(:,:,:,cm_index,:))
1514
1515
1516
1517
1518
1519
1520
1521            ! Auxiliary variables for nests
1522            WRITE(iounit) xresoln(:)
1523            WRITE(iounit) yresoln(:)
1524            WRITE(iounit) xln(:)
1525            WRITE(iounit) yln(:)
1526            WRITE(iounit) xrn(:)
1527            WRITE(iounit) yrn(:)
1528
1529            ! Variables for polar stereographic projection
1530            WRITE(iounit) xglobal, sglobal, nglobal
1531            WRITE(iounit) switchnorthg, switchsouthg
1532            WRITE(iounit) southpolemap(:)
1533            WRITE(iounit) northpolemap(:)
1534
1535            ! Variables declared in conv_mod (convection)
1536            WRITE(iounit) pconv(:)
1537            WRITE(iounit) phconv(:)
1538            WRITE(iounit) dpr(:)
1539            WRITE(iounit) pconv_hpa(:)
1540            WRITE(iounit) phconv_hpa(:)
1541            WRITE(iounit) ft(:)
1542            WRITE(iounit) fq(:)
1543            WRITE(iounit) fmass(:,:)
1544            WRITE(iounit) sub(:)
1545            WRITE(iounit) fmassfrac(:,:)
1546            WRITE(iounit) cbaseflux(:,:)
1547            WRITE(iounit) cbasefluxn(:,:,:)
1548            WRITE(iounit) tconv(:)
1549            WRITE(iounit) qconv(:)
1550            WRITE(iounit) qsconv(:)
1551            WRITE(iounit) psconv, tt2conv, td2conv
1552            WRITE(iounit) nconvlev, nconvtop
1553
1554        ELSE IF (op == 'LOAD') THEN
1555
1556            ! Read the preprocessed format version string and insure it
1557            ! matches this version
1558            READ (iounit) temp_preproc_format_version_str
1559            PRINT *, 'Reading preprocessed file format version: ', &
1560&                    temp_preproc_format_version_str
1561
1562            IF (TRIM(temp_preproc_format_version_str) == &
1563&                                        TRIM(PREPROC_FORMAT_VERSION_STR)) THEN
1564                CONTINUE
1565            ELSE
1566                ! PRINT *, ''  GK: causes relocation truncated to fit: R_X86_64_32
1567                PRINT *, 'Inconsistent preprocessing format version'
1568                PRINT *, 'Expected Version: ', PREPROC_FORMAT_VERSION_STR
1569                PRINT *, 'Detected Version: ', temp_preproc_format_version_str
1570                ! PRINT *, ''
1571                STOP
1572            END IF
1573
1574            ! Read the compiled max dimensions that were dumped from par_mod
1575            ! when creating the fp file, so that we can compare against
1576            ! current FLEXPART dimensions - they need to be the same, or else
1577            ! we abort.
1578            READ (iounit) temp_nxmax, temp_nymax, temp_nzmax, &
1579&                         temp_nuvzmax, temp_nwzmax
1580
1581            ! Get dimensions
1582            ncret = nf90_inq_dimid(ncid, 'nxmax', nxmax_dimid)
1583            call handle_nf90_err(ncret)
1584            ncret = nf90_inquire_dimension(ncid, nxmax_dimid, nxmax_dimname, &
1585&                                                temp_nxmax)
1586            call handle_nf90_err(ncret)
1587            PRINT *, 'temp_nxmax: ', temp_nxmax
1588
1589            ncret = nf90_inq_dimid(ncid, 'nymax', nymax_dimid)
1590            call handle_nf90_err(ncret)
1591            ncret = nf90_inquire_dimension(ncid, nymax_dimid, nymax_dimname, &
1592&                                                temp_nymax)
1593            call handle_nf90_err(ncret)
1594            PRINT *, 'temp_nymax: ', temp_nymax
1595
1596            ncret = nf90_inq_dimid(ncid, 'nzmax', nzmax_dimid)
1597            call handle_nf90_err(ncret)
1598            ncret = nf90_inquire_dimension(ncid, nzmax_dimid, nzmax_dimname, &
1599&                                                temp_nzmax)
1600            call handle_nf90_err(ncret)
1601            PRINT *, 'temp_nzmax: ', temp_nzmax
1602
1603            ncret = nf90_inq_dimid(ncid, 'nuvzmax', nuvzmax_dimid)
1604            call handle_nf90_err(ncret)
1605            ncret = nf90_inquire_dimension(ncid, nuvzmax_dimid, nuvzmax_dimname, &
1606&                                                temp_nuvzmax)
1607            call handle_nf90_err(ncret)
1608            PRINT *, 'temp_nuvzmax: ', temp_nuvzmax
1609
1610            ncret = nf90_inq_dimid(ncid, 'nwzmax', nwzmax_dimid)
1611            call handle_nf90_err(ncret)
1612            ncret = nf90_inquire_dimension(ncid, nwzmax_dimid, nwzmax_dimname, &
1613&                                                temp_nwzmax)
1614            call handle_nf90_err(ncret)
1615            PRINT *, 'temp_nwzmax: ', temp_nwzmax
1616
1617            ncret = nf90_inq_dimid(ncid, 'numclass', numclass_dimid)
1618            call handle_nf90_err(ncret)
1619            ncret = nf90_inquire_dimension(ncid, numclass_dimid, numclass_dimname, &
1620&                                                temp_numclass)
1621            call handle_nf90_err(ncret)
1622            PRINT *, 'temp_numclass: ', temp_numclass
1623
1624            ncret = nf90_inq_dimid(ncid, 'maxspec', maxspec_dimid)
1625            call handle_nf90_err(ncret)
1626            ncret = nf90_inquire_dimension(ncid, maxspec_dimid, maxspec_dimname, &
1627&                                                temp_maxspec)
1628            call handle_nf90_err(ncret)
1629            PRINT *, 'temp_maxspec: ', temp_maxspec
1630
1631
1632
1633            IF ( (temp_nxmax == nxmax) .AND. (temp_nymax == nymax) .AND. &
1634&                   (temp_nzmax == nzmax) .AND. &
1635&                   (temp_nuvzmax == nuvzmax) .AND. &
1636&                   (temp_nwzmax == nwzmax) .AND. &
1637&                   (temp_numclass == numclass) .AND. &
1638&                   (temp_maxspec == maxspec) ) THEN
1639                CONTINUE
1640            ELSE
1641                PRINT *, 'Incompatible dimensions between fp file and current FLEXPART!'
1642                ! PRINT *, ''
1643                PRINT *, '                  FP file     Compiled FP'
1644                PRINT *, 'nxmax:     ', temp_nxmax, '    ', nxmax
1645                PRINT *, 'nymax:     ', temp_nymax, '    ', nymax
1646                PRINT *, 'nzmax:     ', temp_nzmax, '    ', nzmax
1647                PRINT *, 'nuvzmax:     ', temp_nuvzmax, '    ', nuvzmax
1648                PRINT *, 'nwzmax:     ', temp_nwzmax, '    ', nwzmax
1649                PRINT *, 'numclass:     ', temp_numclass, '    ', numclass
1650                PRINT *, 'maxspec:     ', temp_maxspec, '    ', maxspec
1651                ! PRINT *, ''
1652                STOP
1653            END IF
1654
1655
1656
1657
1658            ! Scalar values
1659            READ(iounit) nx, ny, nxmin1, nymin1, nxfield
1660            READ(iounit) nuvz, nwz, nz, nmixz, nlev_ec
1661            READ(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst
1662
1663
1664
1665            ! Get the varid , then read into scalar variable
1666            ncret = nf90_inq_varid(ncid, 'nx', ncvarid)
1667            call handle_nf90_err(ncret)
1668            ncret = nf90_get_var(ncid, ncvarid, nx)
1669            call handle_nf90_err(ncret)
1670
1671            ncret = nf90_inq_varid(ncid, 'ny', ncvarid)
1672            call handle_nf90_err(ncret)
1673            ncret = nf90_get_var(ncid, ncvarid, ny)
1674            call handle_nf90_err(ncret)
1675
1676            ncret = nf90_inq_varid(ncid, 'nxmin1', ncvarid)
1677            call handle_nf90_err(ncret)
1678            ncret = nf90_get_var(ncid, ncvarid, nxmin1)
1679            call handle_nf90_err(ncret)
1680
1681            ncret = nf90_inq_varid(ncid, 'nymin1', ncvarid)
1682            call handle_nf90_err(ncret)
1683            ncret = nf90_get_var(ncid, ncvarid, nymin1)
1684            call handle_nf90_err(ncret)
1685
1686            ncret = nf90_inq_varid(ncid, 'nxfield', ncvarid)
1687            call handle_nf90_err(ncret)
1688            ncret = nf90_get_var(ncid, ncvarid, nxfield)
1689            call handle_nf90_err(ncret)
1690
1691            ncret = nf90_inq_varid(ncid, 'nuvz', ncvarid)
1692            call handle_nf90_err(ncret)
1693            ncret = nf90_get_var(ncid, ncvarid, nuvz)
1694            call handle_nf90_err(ncret)
1695
1696            ncret = nf90_inq_varid(ncid, 'nwz', ncvarid)
1697            call handle_nf90_err(ncret)
1698            ncret = nf90_get_var(ncid, ncvarid, nwz)
1699            call handle_nf90_err(ncret)
1700
1701            ncret = nf90_inq_varid(ncid, 'nz', ncvarid)
1702            call handle_nf90_err(ncret)
1703            ncret = nf90_get_var(ncid, ncvarid, nz)
1704            call handle_nf90_err(ncret)
1705
1706            ncret = nf90_inq_varid(ncid, 'nmixz', ncvarid)
1707            call handle_nf90_err(ncret)
1708            ncret = nf90_get_var(ncid, ncvarid, nmixz)
1709            call handle_nf90_err(ncret)
1710
1711            ncret = nf90_inq_varid(ncid, 'nlev_ec', ncvarid)
1712            call handle_nf90_err(ncret)
1713            ncret = nf90_get_var(ncid, ncvarid, nlev_ec)
1714            call handle_nf90_err(ncret)
1715
1716            ncret = nf90_inq_varid(ncid, 'dx', ncvarid)
1717            call handle_nf90_err(ncret)
1718            ncret = nf90_get_var(ncid, ncvarid, dx)
1719            call handle_nf90_err(ncret)
1720
1721            ncret = nf90_inq_varid(ncid, 'dy', ncvarid)
1722            call handle_nf90_err(ncret)
1723            ncret = nf90_get_var(ncid, ncvarid, dy)
1724            call handle_nf90_err(ncret)
1725
1726            ncret = nf90_inq_varid(ncid, 'xlon0', ncvarid)
1727            call handle_nf90_err(ncret)
1728            ncret = nf90_get_var(ncid, ncvarid, xlon0)
1729            call handle_nf90_err(ncret)
1730
1731            ncret = nf90_inq_varid(ncid, 'ylat0', ncvarid)
1732            call handle_nf90_err(ncret)
1733            ncret = nf90_get_var(ncid, ncvarid, ylat0)
1734            call handle_nf90_err(ncret)
1735
1736            ncret = nf90_inq_varid(ncid, 'dxconst', ncvarid)
1737            call handle_nf90_err(ncret)
1738            ncret = nf90_get_var(ncid, ncvarid, dxconst)
1739            call handle_nf90_err(ncret)
1740
1741            ncret = nf90_inq_varid(ncid, 'dyconst', ncvarid)
1742            call handle_nf90_err(ncret)
1743            ncret = nf90_get_var(ncid, ncvarid, dyconst)
1744            call handle_nf90_err(ncret)
1745
1746
1747
1748
1749
1750
1751            ! Fixed fields, static in time
1752            READ(iounit) oro, excessoro, lsm, xlanduse, height
1753
1754            ncret = nf90_inq_varid(ncid, 'oro', ncvarid)
1755            call handle_nf90_err(ncret)
1756            ncret = nf90_get_var(ncid, ncvarid, oro(0:nxmax-1,0:nymax-1))
1757            call handle_nf90_err(ncret)
1758
1759            ncret = nf90_inq_varid(ncid, 'excessoro', ncvarid)
1760            call handle_nf90_err(ncret)
1761            ncret = nf90_get_var(ncid, ncvarid, excessoro(0:nxmax-1,0:nymax-1))
1762            call handle_nf90_err(ncret)
1763
1764            ncret = nf90_inq_varid(ncid, 'lsm', ncvarid)
1765            call handle_nf90_err(ncret)
1766            ncret = nf90_get_var(ncid, ncvarid, lsm(0:nxmax-1,0:nymax-1))
1767            call handle_nf90_err(ncret)
1768
1769            ncret = nf90_inq_varid(ncid, 'xlanduse', ncvarid)
1770            call handle_nf90_err(ncret)
1771            ncret = nf90_get_var(ncid, ncvarid, xlanduse(0:nxmax-1,0:nymax-1, 1:numclass))
1772            call handle_nf90_err(ncret)
1773
1774            ncret = nf90_inq_varid(ncid, 'height', ncvarid)
1775            call handle_nf90_err(ncret)
1776            ncret = nf90_get_var(ncid, ncvarid, height(1:nzmax))
1777            call handle_nf90_err(ncret)
1778
1779
1780
1781
1782            ! 3d fields
1783            READ(iounit) uu(:,:,:,cm_index)
1784            READ(iounit) vv(:,:,:,cm_index)
1785            READ(iounit) uupol(:,:,:,cm_index)
1786            READ(iounit) vvpol(:,:,:,cm_index)
1787            READ(iounit) ww(:,:,:,cm_index)
1788            READ(iounit) tt(:,:,:,cm_index)
1789            READ(iounit) qv(:,:,:,cm_index)
1790            READ(iounit) pv(:,:,:,cm_index)
1791            READ(iounit) rho(:,:,:,cm_index)
1792            READ(iounit) drhodz(:,:,:,cm_index)
1793            READ(iounit) tth(:,:,:,cm_index)
1794            READ(iounit) qvh(:,:,:,cm_index)
1795            READ(iounit) pplev(:,:,:,cm_index)
1796            READ(iounit) clouds(:,:,:,cm_index)
1797            READ(iounit) cloudsh(:,:,cm_index)
1798
1799
1800
1801
1802            ! Get the varid and read the variable into the array
1803            ncret = nf90_inq_varid(ncid, 'uu', ncvarid)
1804            call handle_nf90_err(ncret)
1805            ncret = nf90_get_var(ncid, ncvarid, uu(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1806            call handle_nf90_err(ncret)
1807
1808            ncret = nf90_inq_varid(ncid, 'vv', ncvarid)
1809            call handle_nf90_err(ncret)
1810            ncret = nf90_get_var(ncid, ncvarid, vv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1811            call handle_nf90_err(ncret)
1812
1813            ncret = nf90_inq_varid(ncid, 'uupol', ncvarid)
1814            call handle_nf90_err(ncret)
1815            ncret = nf90_get_var(ncid, ncvarid, uupol(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1816            call handle_nf90_err(ncret)
1817
1818            ncret = nf90_inq_varid(ncid, 'vvpol', ncvarid)
1819            call handle_nf90_err(ncret)
1820            ncret = nf90_get_var(ncid, ncvarid, vvpol(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1821            call handle_nf90_err(ncret)
1822
1823            ncret = nf90_inq_varid(ncid, 'ww', ncvarid)
1824            call handle_nf90_err(ncret)
1825            ncret = nf90_get_var(ncid, ncvarid, ww(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1826            call handle_nf90_err(ncret)
1827
1828            ncret = nf90_inq_varid(ncid, 'tt', ncvarid)
1829            call handle_nf90_err(ncret)
1830            ncret = nf90_get_var(ncid, ncvarid, tt(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1831            call handle_nf90_err(ncret)
1832
1833            ncret = nf90_inq_varid(ncid, 'qv', ncvarid)
1834            call handle_nf90_err(ncret)
1835            ncret = nf90_get_var(ncid, ncvarid, qv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1836            call handle_nf90_err(ncret)
1837
1838            ncret = nf90_inq_varid(ncid, 'pv', ncvarid)
1839            call handle_nf90_err(ncret)
1840            ncret = nf90_get_var(ncid, ncvarid, pv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1841            call handle_nf90_err(ncret)
1842
1843            ncret = nf90_inq_varid(ncid, 'rho', ncvarid)
1844            call handle_nf90_err(ncret)
1845            ncret = nf90_get_var(ncid, ncvarid, rho(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1846            call handle_nf90_err(ncret)
1847
1848            ncret = nf90_inq_varid(ncid, 'drhodz', ncvarid)
1849            call handle_nf90_err(ncret)
1850            ncret = nf90_get_var(ncid, ncvarid, drhodz(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1851            call handle_nf90_err(ncret)
1852
1853            ncret = nf90_inq_varid(ncid, 'clouds', ncvarid)
1854            call handle_nf90_err(ncret)
1855            ncret = nf90_get_var(ncid, ncvarid, clouds(0:nxmax-1,0:nzmax-1,1:nzmax,cm_index))
1856            call handle_nf90_err(ncret)
1857
1858            ncret = nf90_inq_varid(ncid, 'tth', ncvarid)
1859            call handle_nf90_err(ncret)
1860            ncret = nf90_get_var(ncid, ncvarid, tth(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
1861            call handle_nf90_err(ncret)
1862
1863            ncret = nf90_inq_varid(ncid, 'qvh', ncvarid)
1864            call handle_nf90_err(ncret)
1865            ncret = nf90_get_var(ncid, ncvarid, qvh(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
1866            call handle_nf90_err(ncret)
1867
1868            ncret = nf90_inq_varid(ncid, 'pplev', ncvarid)
1869            call handle_nf90_err(ncret)
1870            ncret = nf90_get_var(ncid, ncvarid, pplev(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
1871            call handle_nf90_err(ncret)
1872
1873            ncret = nf90_inq_varid(ncid, 'cloudsh', ncvarid)
1874            call handle_nf90_err(ncret)
1875            ncret = nf90_get_var(ncid, ncvarid, cloudsh(0:nxmax-1,0:nymax-1,cm_index))
1876            call handle_nf90_err(ncret)
1877
1878
1879
1880
1881            PRINT *, 'SUM(tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
1882&                                        SUM(tt(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
1883
1884
1885            PRINT *, 'SUM(clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
1886&                                        SUM(clouds(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
1887
1888
1889
1890            ! 2d fields
1891            READ(iounit) ps(:,:,:,cm_index)
1892            READ(iounit) sd(:,:,:,cm_index)
1893            READ(iounit) msl(:,:,:,cm_index)
1894            READ(iounit) tcc(:,:,:,cm_index)
1895            READ(iounit) u10(:,:,:,cm_index)
1896            READ(iounit) v10(:,:,:,cm_index)
1897            READ(iounit) tt2(:,:,:,cm_index)
1898            READ(iounit) td2(:,:,:,cm_index)
1899            READ(iounit) lsprec(:,:,:,cm_index)
1900            READ(iounit) convprec(:,:,:,cm_index)
1901            READ(iounit) sshf(:,:,:,cm_index)
1902            READ(iounit) ssr(:,:,:,cm_index)
1903            READ(iounit) surfstr(:,:,:,cm_index)
1904            READ(iounit) ustar(:,:,:,cm_index)
1905            READ(iounit) wstar(:,:,:,cm_index)
1906            READ(iounit) hmix(:,:,:,cm_index)
1907            READ(iounit) tropopause(:,:,:,cm_index)
1908            READ(iounit) oli(:,:,:,cm_index)
1909            READ(iounit) diffk(:,:,:,cm_index)
1910            READ(iounit) vdep(:,:,:,cm_index)
1911
1912            ! Get the varid and read the variable into the array
1913            ncret = nf90_inq_varid(ncid, 'ps', ncvarid)
1914            call handle_nf90_err(ncret)
1915            ncret = nf90_get_var(ncid, ncvarid, ps(0:nxmax-1,0:nymax-1,1,cm_index))
1916            call handle_nf90_err(ncret)
1917
1918            ncret = nf90_inq_varid(ncid, 'sd', ncvarid)
1919            call handle_nf90_err(ncret)
1920            ncret = nf90_get_var(ncid, ncvarid, sd(0:nxmax-1,0:nymax-1,1,cm_index))
1921            call handle_nf90_err(ncret)
1922
1923            ncret = nf90_inq_varid(ncid, 'msl', ncvarid)
1924            call handle_nf90_err(ncret)
1925            ncret = nf90_get_var(ncid, ncvarid, msl(0:nxmax-1,0:nymax-1,1,cm_index))
1926            call handle_nf90_err(ncret)
1927
1928            ncret = nf90_inq_varid(ncid, 'tcc', ncvarid)
1929            call handle_nf90_err(ncret)
1930            ncret = nf90_get_var(ncid, ncvarid, tcc(0:nxmax-1,0:nymax-1,1,cm_index))
1931            call handle_nf90_err(ncret)
1932
1933            ncret = nf90_inq_varid(ncid, 'u10', ncvarid)
1934            call handle_nf90_err(ncret)
1935            ncret = nf90_get_var(ncid, ncvarid, u10(0:nxmax-1,0:nymax-1,1,cm_index))
1936            call handle_nf90_err(ncret)
1937
1938            ncret = nf90_inq_varid(ncid, 'v10', ncvarid)
1939            call handle_nf90_err(ncret)
1940            ncret = nf90_get_var(ncid, ncvarid, v10(0:nxmax-1,0:nymax-1,1,cm_index))
1941            call handle_nf90_err(ncret)
1942
1943            ncret = nf90_inq_varid(ncid, 'tt2', ncvarid)
1944            call handle_nf90_err(ncret)
1945            ncret = nf90_get_var(ncid, ncvarid, tt2(0:nxmax-1,0:nymax-1,1,cm_index))
1946            call handle_nf90_err(ncret)
1947
1948            ncret = nf90_inq_varid(ncid, 'td2', ncvarid)
1949            call handle_nf90_err(ncret)
1950            ncret = nf90_get_var(ncid, ncvarid, td2(0:nxmax-1,0:nymax-1,1,cm_index))
1951            call handle_nf90_err(ncret)
1952
1953            ncret = nf90_inq_varid(ncid, 'lsprec', ncvarid)
1954            call handle_nf90_err(ncret)
1955            ncret = nf90_get_var(ncid, ncvarid, lsprec(0:nxmax-1,0:nymax-1,1,cm_index))
1956            call handle_nf90_err(ncret)
1957
1958            ncret = nf90_inq_varid(ncid, 'convprec', ncvarid)
1959            call handle_nf90_err(ncret)
1960            ncret = nf90_get_var(ncid, ncvarid, convprec(0:nxmax-1,0:nymax-1,1,cm_index))
1961            call handle_nf90_err(ncret)
1962
1963            ncret = nf90_inq_varid(ncid, 'sshf', ncvarid)
1964            call handle_nf90_err(ncret)
1965            ncret = nf90_get_var(ncid, ncvarid, sshf(0:nxmax-1,0:nymax-1,1,cm_index))
1966            call handle_nf90_err(ncret)
1967
1968            ncret = nf90_inq_varid(ncid, 'ssr', ncvarid)
1969            call handle_nf90_err(ncret)
1970            ncret = nf90_get_var(ncid, ncvarid, ssr(0:nxmax-1,0:nymax-1,1,cm_index))
1971            call handle_nf90_err(ncret)
1972
1973            ncret = nf90_inq_varid(ncid, 'surfstr', ncvarid)
1974            call handle_nf90_err(ncret)
1975            ncret = nf90_get_var(ncid, ncvarid, surfstr(0:nxmax-1,0:nymax-1,1,cm_index))
1976            call handle_nf90_err(ncret)
1977
1978            ncret = nf90_inq_varid(ncid, 'ustar', ncvarid)
1979            call handle_nf90_err(ncret)
1980            ncret = nf90_get_var(ncid, ncvarid, ustar(0:nxmax-1,0:nymax-1,1,cm_index))
1981            call handle_nf90_err(ncret)
1982
1983            ncret = nf90_inq_varid(ncid, 'wstar', ncvarid)
1984            call handle_nf90_err(ncret)
1985            ncret = nf90_get_var(ncid, ncvarid, wstar(0:nxmax-1,0:nymax-1,1,cm_index))
1986            call handle_nf90_err(ncret)
1987
1988            ncret = nf90_inq_varid(ncid, 'hmix', ncvarid)
1989            call handle_nf90_err(ncret)
1990            ncret = nf90_get_var(ncid, ncvarid, hmix(0:nxmax-1,0:nymax-1,1,cm_index))
1991            call handle_nf90_err(ncret)
1992
1993            ncret = nf90_inq_varid(ncid, 'tropopause', ncvarid)
1994            call handle_nf90_err(ncret)
1995            ncret = nf90_get_var(ncid, ncvarid, tropopause(0:nxmax-1,0:nymax-1,1,cm_index))
1996            call handle_nf90_err(ncret)
1997
1998            ncret = nf90_inq_varid(ncid, 'oli', ncvarid)
1999            call handle_nf90_err(ncret)
2000            ncret = nf90_get_var(ncid, ncvarid, oli(0:nxmax-1,0:nymax-1,1,cm_index))
2001            call handle_nf90_err(ncret)
2002
2003            ncret = nf90_inq_varid(ncid, 'diffk', ncvarid)
2004            call handle_nf90_err(ncret)
2005            ncret = nf90_get_var(ncid, ncvarid, diffk(0:nxmax-1,0:nymax-1,1,cm_index))
2006            call handle_nf90_err(ncret)
2007
2008
2009
2010
2011            PRINT *, 'SUM(ps(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
2012&                                        SUM(ps(0:nxmax-1,0:nymax-1,1, cm_index))
2013
2014            PRINT *, 'SUM(wstar(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
2015&                                        SUM(wstar(0:nxmax-1,0:nymax-1,1, cm_index))
2016
2017
2018            ncret = nf90_inq_varid(ncid, 'vdep', ncvarid)
2019            call handle_nf90_err(ncret)
2020            ncret = nf90_get_var(ncid, ncvarid, vdep(0:nxmax-1,0:nymax-1,1:maxspec, cm_index))
2021            call handle_nf90_err(ncret)
2022
2023
2024
2025
2026
2027            ! 1d fields
2028            READ(iounit) z0(:)
2029            READ(iounit) akm(:)
2030            READ(iounit) bkm(:)
2031            READ(iounit) akz(:)
2032            READ(iounit) bkz(:)
2033            READ(iounit) aknew(:)
2034            READ(iounit) bknew(:)
2035
2036
2037            ncret = nf90_inq_varid(ncid, 'z0', ncvarid)
2038            call handle_nf90_err(ncret)
2039            ncret = nf90_get_var(ncid, ncvarid, z0(1:numclass))
2040            call handle_nf90_err(ncret)
2041
2042            ncret = nf90_inq_varid(ncid, 'akm', ncvarid)
2043            call handle_nf90_err(ncret)
2044            ncret = nf90_get_var(ncid, ncvarid, akm(1:nwzmax))
2045            call handle_nf90_err(ncret)
2046
2047            ncret = nf90_inq_varid(ncid, 'bkm', ncvarid)
2048            call handle_nf90_err(ncret)
2049            ncret = nf90_get_var(ncid, ncvarid, bkm(1:nwzmax))
2050            call handle_nf90_err(ncret)
2051
2052            ncret = nf90_inq_varid(ncid, 'akz', ncvarid)
2053            call handle_nf90_err(ncret)
2054            ncret = nf90_get_var(ncid, ncvarid, akz(1:nuvzmax))
2055            call handle_nf90_err(ncret)
2056
2057            ncret = nf90_inq_varid(ncid, 'bkz', ncvarid)
2058            call handle_nf90_err(ncret)
2059            ncret = nf90_get_var(ncid, ncvarid, bkz(1:nuvzmax))
2060            call handle_nf90_err(ncret)
2061
2062            ncret = nf90_inq_varid(ncid, 'aknew', ncvarid)
2063            call handle_nf90_err(ncret)
2064            ncret = nf90_get_var(ncid, ncvarid, aknew(1:nzmax))
2065            call handle_nf90_err(ncret)
2066
2067            ncret = nf90_inq_varid(ncid, 'bknew', ncvarid)
2068            call handle_nf90_err(ncret)
2069            ncret = nf90_get_var(ncid, ncvarid, bknew(1:nzmax))
2070            call handle_nf90_err(ncret)
2071
2072
2073
2074            PRINT *, 'SUM(bknew(1:nzmax)): ', &
2075&                                        SUM(bknew(1:nzmax))
2076
2077
2078
2079
2080            ! Now the nested input grid variables
2081            ! Get the compiled values that were written into the FP file, and
2082            ! make sure they are equal to the current compiled values, to make
2083            ! sure we are working with consistent arrays
2084            ncret = nf90_inq_dimid(ncid, 'maxnests', maxnests_dimid)
2085            call handle_nf90_err(ncret)
2086            ncret = nf90_inquire_dimension(ncid, maxnests_dimid, maxnests_dimname, &
2087&                                                temp_maxnests)
2088            call handle_nf90_err(ncret)
2089            PRINT *, 'temp_maxnests: ', temp_maxnests
2090
2091            ncret = nf90_inq_dimid(ncid, 'nxmaxn', nxmaxn_dimid)
2092            call handle_nf90_err(ncret)
2093            ncret = nf90_inquire_dimension(ncid, nxmaxn_dimid, nxmaxn_dimname, &
2094&                                                temp_nxmaxn)
2095            call handle_nf90_err(ncret)
2096            PRINT *, 'temp_nxmaxn: ', temp_nxmaxn
2097
2098            ncret = nf90_inq_dimid(ncid, 'nymaxn', nymaxn_dimid)
2099            call handle_nf90_err(ncret)
2100            ncret = nf90_inquire_dimension(ncid, nymaxn_dimid, nymaxn_dimname, &
2101&                                                temp_nymaxn)
2102            call handle_nf90_err(ncret)
2103            PRINT *, 'temp_nymaxn: ', temp_nymaxn
2104
2105            ! Note that maxspec_dimid and numclass_dimid were checked above
2106
2107
2108
2109
2110            IF ( (temp_nxmaxn == nxmaxn) .AND. (temp_nymaxn == nymaxn) .AND. &
2111&                   (temp_maxnests == maxnests) ) THEN
2112                CONTINUE
2113            ELSE
2114                PRINT *, 'Incompatible dimensions between fp file and current FLEXPART!'
2115                ! PRINT *, ''
2116                PRINT *, '                  FP file     Compiled FP'
2117                PRINT *, 'nxmaxn:     ', temp_nxmaxn, '    ', nxmaxn
2118                PRINT *, 'nymaxn:     ', temp_nymaxn, '    ', nymaxn
2119                PRINT *, 'maxnests:     ', temp_maxnests, '    ', maxnests
2120                STOP
2121            END IF
2122
2123
2124
2125            ! Nested, scalar values (for each nest)
2126            READ(iounit) nxn(:)
2127            READ(iounit) nyn(:)
2128            READ(iounit) dxn(:)
2129            READ(iounit) dyn(:)
2130            READ(iounit) xlon0n(:)
2131            READ(iounit) ylat0n(:)
2132
2133
2134            ncret = nf90_inq_varid(ncid, 'nxn', ncvarid)
2135            call handle_nf90_err(ncret)
2136            ncret = nf90_get_var(ncid, ncvarid, nxn(1:maxnests))
2137            call handle_nf90_err(ncret)
2138
2139            ncret = nf90_inq_varid(ncid, 'nyn', ncvarid)
2140            call handle_nf90_err(ncret)
2141            ncret = nf90_get_var(ncid, ncvarid, nyn(1:maxnests))
2142            call handle_nf90_err(ncret)
2143
2144            ncret = nf90_inq_varid(ncid, 'dxn', ncvarid)
2145            call handle_nf90_err(ncret)
2146            ncret = nf90_get_var(ncid, ncvarid, dxn(1:maxnests))
2147            call handle_nf90_err(ncret)
2148
2149            ncret = nf90_inq_varid(ncid, 'dyn', ncvarid)
2150            call handle_nf90_err(ncret)
2151            ncret = nf90_get_var(ncid, ncvarid, dyn(1:maxnests))
2152            call handle_nf90_err(ncret)
2153
2154            ncret = nf90_inq_varid(ncid, 'xlon0n', ncvarid)
2155            call handle_nf90_err(ncret)
2156            ncret = nf90_get_var(ncid, ncvarid, xlon0n(1:maxnests))
2157            call handle_nf90_err(ncret)
2158
2159            ncret = nf90_inq_varid(ncid, 'ylat0n', ncvarid)
2160            call handle_nf90_err(ncret)
2161            ncret = nf90_get_var(ncid, ncvarid, ylat0n(1:maxnests))
2162            call handle_nf90_err(ncret)
2163
2164
2165
2166            ! Nested fields, static over time
2167            READ(iounit) oron, excessoron, lsmn, xlandusen
2168
2169            ncret = nf90_inq_varid(ncid, 'oron', ncvarid)
2170            call handle_nf90_err(ncret)
2171            ncret = nf90_get_var(ncid, ncvarid, oron(0:nxmaxn-1,0:nymaxn-1,1:maxnests))
2172            call handle_nf90_err(ncret)
2173
2174            ncret = nf90_inq_varid(ncid, 'excessoron', ncvarid)
2175            call handle_nf90_err(ncret)
2176            ncret = nf90_get_var(ncid, ncvarid, excessoron(0:nxmaxn-1,0:nymaxn-1,1:maxnests))
2177            call handle_nf90_err(ncret)
2178
2179            ncret = nf90_inq_varid(ncid, 'lsmn', ncvarid)
2180            call handle_nf90_err(ncret)
2181            ncret = nf90_get_var(ncid, ncvarid, lsmn(0:nxmaxn-1,0:nymaxn-1,1:maxnests))
2182            call handle_nf90_err(ncret)
2183
2184            ncret = nf90_inq_varid(ncid, 'xlandusen', ncvarid)
2185            call handle_nf90_err(ncret)
2186            ncret = nf90_get_var(ncid, ncvarid, xlandusen(0:nxmaxn-1,0:nymaxn-1,1:numclass,1:maxnests))
2187            call handle_nf90_err(ncret)
2188
2189
2190            PRINT *, 'SUM(oron): ', SUM(oron)
2191
2192
2193
2194
2195            ! 3d nested fields
2196            READ(iounit) uun(:,:,:,cm_index,:)
2197            READ(iounit) vvn(:,:,:,cm_index,:)
2198            READ(iounit) wwn(:,:,:,cm_index,:)
2199            READ(iounit) ttn(:,:,:,cm_index,:)
2200            READ(iounit) qvn(:,:,:,cm_index,:)
2201            READ(iounit) pvn(:,:,:,cm_index,:)
2202            READ(iounit) cloudsn(:,:,:,cm_index,:)
2203            READ(iounit) cloudsnh(:,:,cm_index,:)
2204            READ(iounit) rhon(:,:,:,cm_index,:)
2205            READ(iounit) drhodzn(:,:,:,cm_index,:)
2206            READ(iounit) tthn(:,:,:,cm_index,:)
2207            READ(iounit) qvhn(:,:,:,cm_index,:)
2208
2209
2210            ncret = nf90_inq_varid(ncid, 'uun', ncvarid)
2211            call handle_nf90_err(ncret)
2212            ncret = nf90_get_var(ncid, ncvarid, uun(0:nxmaxn-1,0:nymaxn-1,1:nzmax,cm_index,1:maxnests))
2213            call handle_nf90_err(ncret)
2214
2215            ncret = nf90_inq_varid(ncid, 'vvn', ncvarid)
2216            call handle_nf90_err(ncret)
2217            ncret = nf90_get_var(ncid, ncvarid, vvn(0:nxmaxn-1,0:nymaxn-1,1:nzmax,cm_index,1:maxnests))
2218            call handle_nf90_err(ncret)
2219
2220            ncret = nf90_inq_varid(ncid, 'wwn', ncvarid)
2221            call handle_nf90_err(ncret)
2222            ncret = nf90_get_var(ncid, ncvarid, wwn(0:nxmaxn-1,0:nymaxn-1,1:nzmax,cm_index,1:maxnests))
2223            call handle_nf90_err(ncret)
2224
2225            ncret = nf90_inq_varid(ncid, 'ttn', ncvarid)
2226            call handle_nf90_err(ncret)
2227            ncret = nf90_get_var(ncid, ncvarid, ttn(0:nxmaxn-1,0:nymaxn-1,1:nzmax,cm_index,1:maxnests))
2228            call handle_nf90_err(ncret)
2229
2230            ncret = nf90_inq_varid(ncid, 'qvn', ncvarid)
2231            call handle_nf90_err(ncret)
2232            ncret = nf90_get_var(ncid, ncvarid, qvn(0:nxmaxn-1,0:nymaxn-1,1:nzmax,cm_index,1:maxnests))
2233            call handle_nf90_err(ncret)
2234
2235            ncret = nf90_inq_varid(ncid, 'pvn', ncvarid)
2236            call handle_nf90_err(ncret)
2237            ncret = nf90_get_var(ncid, ncvarid, pvn(0:nxmaxn-1,0:nymaxn-1,1:nzmax,cm_index,1:maxnests))
2238            call handle_nf90_err(ncret)
2239
2240            ncret = nf90_inq_varid(ncid, 'rhon', ncvarid)
2241            call handle_nf90_err(ncret)
2242            ncret = nf90_get_var(ncid, ncvarid, rhon(0:nxmaxn-1,0:nymaxn-1,1:nzmax,cm_index,1:maxnests))
2243            call handle_nf90_err(ncret)
2244
2245            ncret = nf90_inq_varid(ncid, 'drhodzn', ncvarid)
2246            call handle_nf90_err(ncret)
2247            ncret = nf90_get_var(ncid, ncvarid, drhodzn(0:nxmaxn-1,0:nymaxn-1,1:nzmax,cm_index,1:maxnests))
2248            call handle_nf90_err(ncret)
2249
2250            ncret = nf90_inq_varid(ncid, 'tthn', ncvarid)
2251            call handle_nf90_err(ncret)
2252            ncret = nf90_get_var(ncid, ncvarid, tthn(0:nxmaxn-1,0:nymaxn-1,1:nuvzmax,cm_index,1:maxnests))
2253            call handle_nf90_err(ncret)
2254
2255            ncret = nf90_inq_varid(ncid, 'qvhn', ncvarid)
2256            call handle_nf90_err(ncret)
2257            ncret = nf90_get_var(ncid, ncvarid, qvhn(0:nxmaxn-1,0:nymaxn-1,1:nuvzmax,cm_index,1:maxnests))
2258            call handle_nf90_err(ncret)
2259
2260            ncret = nf90_inq_varid(ncid, 'cloudsn', ncvarid)
2261            call handle_nf90_err(ncret)
2262            ncret = nf90_get_var(ncid, ncvarid, cloudsn(0:nxmaxn-1,0:nymaxn-1,0:nzmax,cm_index,1:maxnests))
2263            call handle_nf90_err(ncret)
2264
2265            ncret = nf90_inq_varid(ncid, 'cloudsnh', ncvarid)
2266            call handle_nf90_err(ncret)
2267            ncret = nf90_get_var(ncid, ncvarid, cloudsnh(0:nxmaxn-1,0:nymaxn-1,cm_index,1:maxnests))
2268            call handle_nf90_err(ncret)
2269
2270
2271
2272
2273            PRINT *, 'SUM(uun): ', SUM(uun(:,:,:,cm_index,:))
2274            PRINT *, 'SUM(qvhn): ', SUM(qvhn(:,:,:,cm_index,:))
2275            PRINT *, 'SUM(cloudsn): ', SUM(cloudsn(:,:,:,cm_index,:))
2276
2277
2278
2279
2280            ! 2d nested fields
2281            READ(iounit) psn(:,:,:,cm_index,:)
2282            READ(iounit) sdn(:,:,:,cm_index,:)
2283            READ(iounit) msln(:,:,:,cm_index,:)
2284            READ(iounit) tccn(:,:,:,cm_index,:)
2285            READ(iounit) u10n(:,:,:,cm_index,:)
2286            READ(iounit) v10n(:,:,:,cm_index,:)
2287            READ(iounit) tt2n(:,:,:,cm_index,:)
2288            READ(iounit) td2n(:,:,:,cm_index,:)
2289            READ(iounit) lsprecn(:,:,:,cm_index,:)
2290            READ(iounit) convprecn(:,:,:,cm_index,:)
2291            READ(iounit) sshfn(:,:,:,cm_index,:)
2292            READ(iounit) ssrn(:,:,:,cm_index,:)
2293            READ(iounit) surfstrn(:,:,:,cm_index,:)
2294            READ(iounit) ustarn(:,:,:,cm_index,:)
2295            READ(iounit) wstarn(:,:,:,cm_index,:)
2296            READ(iounit) hmixn(:,:,:,cm_index,:)
2297            READ(iounit) tropopausen(:,:,:,cm_index,:)
2298            READ(iounit) olin(:,:,:,cm_index,:)
2299            READ(iounit) diffkn(:,:,:,cm_index,:)
2300            READ(iounit) vdepn(:,:,:,cm_index,:)
2301
2302            ncret = nf90_inq_varid(ncid, 'psn', ncvarid)
2303            call handle_nf90_err(ncret)
2304            ncret = nf90_get_var(ncid, ncvarid, psn(0:nxmaxn-1,0:nymaxn-1,1,cm_index,1:maxnests))
2305            call handle_nf90_err(ncret)
2306
2307
2308
2309            PRINT *, 'SUM(psn): ', SUM(psn(:,:,:,cm_index,:))
2310
2311
2312
2313            ! Auxiliary variables for nests
2314            READ(iounit) xresoln(:)
2315            READ(iounit) yresoln(:)
2316            READ(iounit) xln(:)
2317            READ(iounit) yln(:)
2318            READ(iounit) xrn(:)
2319            READ(iounit) yrn(:)
2320
2321            ! Variables for polar stereographic projection
2322            READ(iounit) xglobal, sglobal, nglobal
2323            READ(iounit) switchnorthg, switchsouthg
2324            READ(iounit) southpolemap(:)
2325            READ(iounit) northpolemap(:)
2326
2327            ! Variables declared in conv_mod (convection)
2328            READ(iounit) pconv(:)
2329            READ(iounit) phconv(:)
2330            READ(iounit) dpr(:)
2331            READ(iounit) pconv_hpa(:)
2332            READ(iounit) phconv_hpa(:)
2333            READ(iounit) ft(:)
2334            READ(iounit) fq(:)
2335            READ(iounit) fmass(:,:)
2336            READ(iounit) sub(:)
2337            READ(iounit) fmassfrac(:,:)
2338            READ(iounit) cbaseflux(:,:)
2339            READ(iounit) cbasefluxn(:,:,:)
2340            READ(iounit) tconv(:)
2341            READ(iounit) qconv(:)
2342            READ(iounit) qsconv(:)
2343            READ(iounit) psconv, tt2conv, td2conv
2344            READ(iounit) nconvlev, nconvtop
2345
2346        ELSE
2347            STOP 'fpio(): Illegal operation'
2348           
2349        ENDIF
2350    END SUBROUTINE fpio
2351
2352    SUBROUTINE fpmetbinary_filetext(filename, cm_index)
2353
2354        ! This is a utility subroutine meant to be used for testing purposes.
2355        ! It facilitates the text output of variables read in from the
2356        ! specified .fp file.  This routine will easily cause the program
2357        ! to crash due memory allocation issues, particularly when you are
2358        ! trying to text print 3d arrays in a single formatted statetment.
2359       
2360        CHARACTER(LEN=*), INTENT(IN) :: filename
2361        INTEGER, INTENT(IN) :: cm_index           ! Index of last dimension in
2362                                                  ! most com_mod variables.
2363                                                  ! Should be 1 or 2
2364
2365        !OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', status='REPLACE', &
2366        !    form="formatted", access="stream")
2367        OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', &
2368             form="formatted", access="APPEND")
2369
2370        WRITE(IOUNIT_TEXTOUT, *) 'oro: ', oro
2371        WRITE(IOUNIT_TEXTOUT, *) 'excessoro: ', excessoro
2372        WRITE(IOUNIT_TEXTOUT, *) 'lsm: ', lsm
2373        WRITE(IOUNIT_TEXTOUT, *) 'xlanduse: ', xlanduse
2374        WRITE(IOUNIT_TEXTOUT, *) 'height: ', height
2375
2376        WRITE(IOUNIT_TEXTOUT, *) 'uu: ', uu(:,:,:,cm_index)
2377        WRITE(IOUNIT_TEXTOUT, *) 'vv: ', vv(:,:,:,cm_index)
2378        WRITE(IOUNIT_TEXTOUT, *) 'uupol: ', uupol(:,:,:,cm_index)
2379        WRITE(IOUNIT_TEXTOUT, *) 'vvpol: ', vvpol(:,:,:,cm_index)
2380        WRITE(IOUNIT_TEXTOUT, *) 'ww: ', ww(:,:,:,cm_index)
2381        WRITE(IOUNIT_TEXTOUT, *) 'tt: ', tt(:,:,:,cm_index)
2382        WRITE(IOUNIT_TEXTOUT, *) 'qv: ', qv(:,:,:,cm_index)
2383        WRITE(IOUNIT_TEXTOUT, *) 'pv: ', pv(:,:,:,cm_index)
2384        WRITE(IOUNIT_TEXTOUT, *) 'rho: ', rho(:,:,:,cm_index)
2385        WRITE(IOUNIT_TEXTOUT, *) 'drhodz: ', drhodz(:,:,:,cm_index)
2386        WRITE(IOUNIT_TEXTOUT, *) 'tth: ', tth(:,:,:,cm_index)
2387        WRITE(IOUNIT_TEXTOUT, *) 'qvh: ', qvh(:,:,:,cm_index)
2388        WRITE(IOUNIT_TEXTOUT, *) 'pplev: ', pplev(:,:,:,cm_index)
2389        WRITE(IOUNIT_TEXTOUT, *) 'clouds: ', clouds(:,:,:,cm_index)
2390        WRITE(IOUNIT_TEXTOUT, *) 'cloudsh: ', cloudsh(:,:,cm_index)
2391
2392
2393
2394
2395        CLOSE(IOUNIT_TEXTOUT)
2396    END SUBROUTINE fpmetbinary_filetext
2397
2398    subroutine handle_nf90_err(status)
2399
2400        ! Custom routine for checking NF90 error status
2401        ! and aborting if necessary
2402        use netcdf
2403        implicit none
2404        integer, intent (in) :: status
2405
2406        if (status /= nf90_noerr) then
2407            print *, trim(nf90_strerror(status))
2408            stop "Stopped..."
2409        endif
2410    end subroutine handle_nf90_err
2411
2412
2413
2414
2415END MODULE fpmetbinary_mod
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG