source: flexpart.git/flexpart_code/fpmetbinary_mod.F90 @ 1d9681b

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

Continued NC4 for FP development in branch fp9.3.1-20161214-nc4

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