source: flexpart.git/flexpart_code/fpmetbinary_mod.F90 @ 41a2981

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

Incremental save in branch fp9.3.1-20161214-nc4

  • Property mode set to 100644
File size: 68.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&                  maxnests_dimid, nxmaxn_dimid, nymaxn_dimid
311
312        INTEGER, DIMENSION(1) :: dim1dids    ! Dimension IDs for 1D arrays
313        INTEGER, DIMENSION(2) :: dim2dids    ! Dimension IDs for 2D arrays
314        INTEGER, DIMENSION(3) :: dim3dids    ! Dimension IDs for 3D arrays
315
316        ! These are used when loading in dimensions from NC file
317        CHARACTER(LEN=NF90_MAX_NAME) :: nxmax_dimname, nymax_dimname, nzmax_dimname, &
318&                                       nuvzmax_dimname, nwzmax_dimname,&
319&                                       maxnests_dimname, nxmaxn_dimname, nymaxn_dimname
320
321        ! These are temporary variables, used in the LOAD option, for
322        ! comparing against the current values in FLEXPART of nxmax, nymax, ...
323        INTEGER :: temp_nxmax, temp_nymax, temp_nzmax, &
324&                  temp_nuvzmax, temp_nwzmax, &
325&                   temp_maxnests, temp_nxmaxn, temp_nymaxn
326
327        CHARACTER(LEN=12) :: temp_preproc_format_version_str
328
329        CHARACTER(LEN=128) :: errmesg
330
331        INTEGER, PARAMETER :: DEF_LEVEL = 3
332
333        if (op == 'DUMP') THEN
334
335
336            ! Write the preprocessing format version string
337            WRITE (iounit) PREPROC_FORMAT_VERSION_STR
338
339            ! Write the compiled max dimensions from par_mod - these are
340            ! not meant to be reassigned during a LOAD, but used as "header"
341            ! information to provide the structure of arrays
342            WRITE (iounit) nxmax, nymax, nzmax, nuvzmax, nwzmax
343
344            ncret = nf90_def_dim(ncid, 'nxmax', nxmax, nxmax_dimid)
345            ncret = nf90_def_dim(ncid, 'nymax', nymax, nymax_dimid)
346            ncret = nf90_def_dim(ncid, 'nzmax', nzmax, nzmax_dimid)
347            ncret = nf90_def_dim(ncid, 'nuvzmax', nuvzmax, nuvzmax_dimid)
348            ncret = nf90_def_dim(ncid, 'nwzmax', nwzmax, nwzmax_dimid)
349
350
351
352            ! Scalar values
353            WRITE(iounit) nx, ny, nxmin1, nymin1, nxfield
354            WRITE(iounit) nuvz, nwz, nz, nmixz, nlev_ec
355            WRITE(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst
356
357            ncret = nf90_def_var(ncid, 'nx', NF90_INT, ncvarid)
358            ncret = nf90_put_var(ncid, ncvarid, nx)
359
360            ncret = nf90_def_var(ncid, 'ny', NF90_INT, ncvarid)
361            ncret = nf90_put_var(ncid, ncvarid, ny)
362
363            ncret = nf90_def_var(ncid, 'nxmin1', NF90_INT, ncvarid)
364            ncret = nf90_put_var(ncid, ncvarid, nxmin1)
365
366            ncret = nf90_def_var(ncid, 'nymin1', NF90_INT, ncvarid)
367            ncret = nf90_put_var(ncid, ncvarid, nymin1)
368
369            ncret = nf90_def_var(ncid, 'nxfield', NF90_INT, ncvarid)
370            ncret = nf90_put_var(ncid, ncvarid, nxfield)
371
372            ncret = nf90_def_var(ncid, 'nuvz', NF90_INT, ncvarid)
373            ncret = nf90_put_var(ncid, ncvarid, nuvz)
374
375            ncret = nf90_def_var(ncid, 'nwz', NF90_INT, ncvarid)
376            ncret = nf90_put_var(ncid, ncvarid, nwz)
377
378            ncret = nf90_def_var(ncid, 'nz', NF90_INT, ncvarid)
379            ncret = nf90_put_var(ncid, ncvarid, nz)
380
381            ncret = nf90_def_var(ncid, 'nmixz', NF90_INT, ncvarid)
382            ncret = nf90_put_var(ncid, ncvarid, nmixz)
383
384            ncret = nf90_def_var(ncid, 'nlev_', NF90_INT, ncvarid)
385            ncret = nf90_put_var(ncid, ncvarid, nlev_ec)
386
387            ncret = nf90_def_var(ncid, 'dx', NF90_FLOAT, ncvarid)
388            ncret = nf90_put_var(ncid, ncvarid, dx)
389
390            ncret = nf90_def_var(ncid, 'dy', NF90_FLOAT, ncvarid)
391            ncret = nf90_put_var(ncid, ncvarid, dy)
392
393            ncret = nf90_def_var(ncid, 'xlon0', NF90_FLOAT, ncvarid)
394            ncret = nf90_put_var(ncid, ncvarid, xlon0)
395
396            ncret = nf90_def_var(ncid, 'ylat0', NF90_FLOAT, ncvarid)
397            ncret = nf90_put_var(ncid, ncvarid, ylat0)
398
399            ncret = nf90_def_var(ncid, 'dxconst', NF90_FLOAT, ncvarid)
400            ncret = nf90_put_var(ncid, ncvarid, dxconst)
401
402            ncret = nf90_def_var(ncid, 'dyconst', NF90_FLOAT, ncvarid)
403            ncret = nf90_put_var(ncid, ncvarid, dyconst)
404
405
406
407            ! Fixed fields, static in time
408            WRITE(iounit) oro, excessoro, lsm, xlanduse, height
409
410            dim2dids = (/nxmax_dimid, nymax_dimid/)
411
412            ncret = nf90_def_var(ncid, 'oro', NF90_FLOAT, &
413&                                       dim2dids, ncvarid)
414            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
415&                                        shuffle=0,     &
416&                                        deflate=1,     &
417&                                        deflate_level=DEF_LEVEL)
418            ncret = nf90_put_var(ncid, ncvarid, &
419&                                oro(0:nxmax-1, 0:nymax-1))
420
421            ncret = nf90_def_var(ncid, 'excessoro', NF90_FLOAT, &
422&                                       dim2dids, ncvarid)
423            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
424&                                        shuffle=0,     &
425&                                        deflate=1,     &
426&                                        deflate_level=DEF_LEVEL)
427            ncret = nf90_put_var(ncid, ncvarid, &
428&                                excessoro(0:nxmax-1, 0:nymax-1))
429
430            ncret = nf90_def_var(ncid, 'lsm', NF90_FLOAT, &
431&                                       dim2dids, ncvarid)
432            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
433&                                        shuffle=0,     &
434&                                        deflate=1,     &
435&                                        deflate_level=DEF_LEVEL)
436            ncret = nf90_put_var(ncid, ncvarid, &
437&                                lsm(0:nxmax-1, 0:nymax-1))
438
439            dim3dids = (/nxmax_dimid, nymax_dimid, numclass/)
440            ! numclass comes from par_mod - number of land use classes
441            ncret = nf90_def_var(ncid, 'xlanduse', NF90_FLOAT, &
442&                                       dim3dids, ncvarid)
443            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
444&                                        shuffle=0,     &
445&                                        deflate=1,     &
446&                                        deflate_level=DEF_LEVEL)
447            ncret = nf90_put_var(ncid, ncvarid, &
448&                                xlanduse(0:nxmax-1, 0:nymax-1, 1:numclass))
449
450            dim1dids = (/nzmax_dimid/)
451            ncret = nf90_def_var(ncid, 'height', NF90_FLOAT, &
452&                                       dim1dids, ncvarid)
453            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
454&                                        shuffle=0,     &
455&                                        deflate=1,     &
456&                                        deflate_level=DEF_LEVEL)
457            ncret = nf90_put_var(ncid, ncvarid, &
458&                                height(1:nzmax))
459
460
461
462
463            ! 3d fields
464            WRITE(iounit) uu(:,:,:,cm_index)
465            WRITE(iounit) vv(:,:,:,cm_index)
466            WRITE(iounit) uupol(:,:,:,cm_index)
467            WRITE(iounit) vvpol(:,:,:,cm_index)
468            WRITE(iounit) ww(:,:,:,cm_index)
469            WRITE(iounit) tt(:,:,:,cm_index)
470            WRITE(iounit) qv(:,:,:,cm_index)
471            WRITE(iounit) pv(:,:,:,cm_index)
472            WRITE(iounit) rho(:,:,:,cm_index)
473            WRITE(iounit) drhodz(:,:,:,cm_index)
474            WRITE(iounit) tth(:,:,:,cm_index)
475            WRITE(iounit) qvh(:,:,:,cm_index)
476            WRITE(iounit) pplev(:,:,:,cm_index)
477            WRITE(iounit) clouds(:,:,:,cm_index)
478            WRITE(iounit) cloudsh(:,:,cm_index)
479
480            dim3dids = (/nxmax_dimid, nymax_dimid, nzmax_dimid/)
481
482            ncret = nf90_def_var(ncid, 'uu', NF90_FLOAT, &
483&                                       dim3dids, ncvarid)
484            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
485&                                        shuffle=0,     &
486&                                        deflate=1,     &
487&                                        deflate_level=DEF_LEVEL)
488            ncret = nf90_put_var(ncid, ncvarid, &
489&                                uu(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
490
491            ncret = nf90_def_var(ncid, 'vv', NF90_FLOAT, &
492&                                       dim3dids, ncvarid)
493            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
494&                                        shuffle=0,     &
495&                                        deflate=1,     &
496&                                        deflate_level=DEF_LEVEL)
497            ncret = nf90_put_var(ncid, ncvarid, &
498&                                vv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
499
500            ncret = nf90_def_var(ncid, 'uupol', NF90_FLOAT, &
501&                                       dim3dids, ncvarid)
502            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
503&                                        shuffle=0,     &
504&                                        deflate=1,     &
505&                                        deflate_level=DEF_LEVEL)
506            ncret = nf90_put_var(ncid, ncvarid, &
507&                                uupol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
508
509            ncret = nf90_def_var(ncid, 'vvpol', NF90_FLOAT, &
510&                                       dim3dids, ncvarid)
511            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
512&                                        shuffle=0,     &
513&                                        deflate=1,     &
514&                                        deflate_level=DEF_LEVEL)
515            ncret = nf90_put_var(ncid, ncvarid, &
516&                                vvpol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
517
518            ncret = nf90_def_var(ncid, 'ww', NF90_FLOAT, &
519&                                       dim3dids, ncvarid)
520            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
521&                                        shuffle=0,     &
522&                                        deflate=1,     &
523&                                        deflate_level=DEF_LEVEL)
524            ncret = nf90_put_var(ncid, ncvarid, &
525&                                ww(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
526
527            ncret = nf90_def_var(ncid, 'tt', NF90_FLOAT, &
528&                                       dim3dids, ncvarid)
529            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
530&                                        shuffle=0,     &
531&                                        deflate=1,     &
532&                                        deflate_level=DEF_LEVEL)
533            ncret = nf90_put_var(ncid, ncvarid, &
534&                                tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
535
536            ncret = nf90_def_var(ncid, 'qv', NF90_FLOAT, &
537&                                       dim3dids, ncvarid)
538            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
539&                                        shuffle=0,     &
540&                                        deflate=1,     &
541&                                        deflate_level=DEF_LEVEL)
542            ncret = nf90_put_var(ncid, ncvarid, &
543&                                qv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
544
545            ncret = nf90_def_var(ncid, 'pv', NF90_FLOAT, &
546&                                       dim3dids, ncvarid)
547            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
548&                                        shuffle=0,     &
549&                                        deflate=1,     &
550&                                        deflate_level=DEF_LEVEL)
551            ncret = nf90_put_var(ncid, ncvarid, &
552&                                pv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
553
554            ncret = nf90_def_var(ncid, 'rho', NF90_FLOAT, &
555&                                       dim3dids, ncvarid)
556            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
557&                                        shuffle=0,     &
558&                                        deflate=1,     &
559&                                        deflate_level=DEF_LEVEL)
560            ncret = nf90_put_var(ncid, ncvarid, &
561&                                rho(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
562
563            ncret = nf90_def_var(ncid, 'drhodz', NF90_FLOAT, &
564&                                       dim3dids, ncvarid)
565            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
566&                                        shuffle=0,     &
567&                                        deflate=1,     &
568&                                        deflate_level=DEF_LEVEL)
569            ncret = nf90_put_var(ncid, ncvarid, &
570&                                drhodz(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
571
572            ncret = nf90_def_var(ncid, 'clouds', NF90_BYTE, &
573&                                       dim3dids, ncvarid)
574            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
575&                                        shuffle=0,     &
576&                                        deflate=1,     &
577&                                        deflate_level=DEF_LEVEL)
578            ncret = nf90_put_var(ncid, ncvarid, &
579&                                clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
580
581
582
583            ! Note the change in z dimension for the following
584            dim3dids = (/nxmax_dimid, nymax_dimid, nuvzmax_dimid/)
585
586            ncret = nf90_def_var(ncid, 'tth', NF90_FLOAT, &
587&                                       dim3dids, ncvarid)
588            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
589&                                        shuffle=0,     &
590&                                        deflate=1,     &
591&                                        deflate_level=DEF_LEVEL)
592            ncret = nf90_put_var(ncid, ncvarid, &
593&                                tth(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
594
595            ncret = nf90_def_var(ncid, 'qvh', NF90_FLOAT, &
596&                                       dim3dids, ncvarid)
597            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
598&                                        shuffle=0,     &
599&                                        deflate=1,     &
600&                                        deflate_level=DEF_LEVEL)
601            ncret = nf90_put_var(ncid, ncvarid, &
602&                                qvh(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
603
604            ncret = nf90_def_var(ncid, 'pplev', NF90_FLOAT, &
605&                                       dim3dids, ncvarid)
606            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
607&                                        shuffle=0,     &
608&                                        deflate=1,     &
609&                                        deflate_level=DEF_LEVEL)
610            ncret = nf90_put_var(ncid, ncvarid, &
611&                                pplev(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
612
613
614            dim2dids = (/nxmax_dimid, nymax_dimid/)
615            ncret = nf90_def_var(ncid, 'cloudsh', NF90_INT, &
616&                                       dim2dids, ncvarid)
617            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
618&                                        shuffle=0,     &
619&                                        deflate=1,     &
620&                                        deflate_level=DEF_LEVEL)
621            ncret = nf90_put_var(ncid, ncvarid, &
622&                                cloudsh(0:nxmax-1, 0:nymax-1, cm_index))
623
624
625
626            PRINT *, 'SUM(tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
627&                                        SUM(tt(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
628
629            PRINT *, 'SUM(clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
630&                                        SUM(clouds(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
631
632
633
634            ! 2d fields
635            WRITE(iounit) ps(:,:,:,cm_index)
636            WRITE(iounit) sd(:,:,:,cm_index)
637            WRITE(iounit) msl(:,:,:,cm_index)
638            WRITE(iounit) tcc(:,:,:,cm_index)
639            WRITE(iounit) u10(:,:,:,cm_index)
640            WRITE(iounit) v10(:,:,:,cm_index)
641            WRITE(iounit) tt2(:,:,:,cm_index)
642            WRITE(iounit) td2(:,:,:,cm_index)
643            WRITE(iounit) lsprec(:,:,:,cm_index)
644            WRITE(iounit) convprec(:,:,:,cm_index)
645            WRITE(iounit) sshf(:,:,:,cm_index)
646            WRITE(iounit) ssr(:,:,:,cm_index)
647            WRITE(iounit) surfstr(:,:,:,cm_index)
648            WRITE(iounit) ustar(:,:,:,cm_index)
649            WRITE(iounit) wstar(:,:,:,cm_index)
650            WRITE(iounit) hmix(:,:,:,cm_index)
651            WRITE(iounit) tropopause(:,:,:,cm_index)
652            WRITE(iounit) oli(:,:,:,cm_index)
653            WRITE(iounit) diffk(:,:,:,cm_index)
654            WRITE(iounit) vdep(:,:,:,cm_index)
655
656            dim2dids = (/nxmax_dimid, nymax_dimid/)
657
658            ncret = nf90_def_var(ncid, 'ps', NF90_FLOAT, &
659&                                       dim2dids, ncvarid)
660            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
661&                                        shuffle=0,     &
662&                                        deflate=1,     &
663&                                        deflate_level=DEF_LEVEL)
664            ncret = nf90_put_var(ncid, ncvarid, &
665&                                ps(0:nxmax-1, 0:nymax-1, 1, cm_index))
666
667            ncret = nf90_def_var(ncid, 'sd', NF90_FLOAT, &
668&                                       dim2dids, ncvarid)
669            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
670&                                        shuffle=0,     &
671&                                        deflate=1,     &
672&                                        deflate_level=DEF_LEVEL)
673            ncret = nf90_put_var(ncid, ncvarid, &
674&                                sd(0:nxmax-1, 0:nymax-1, 1, cm_index))
675
676            ncret = nf90_def_var(ncid, 'msl', NF90_FLOAT, &
677&                                       dim2dids, ncvarid)
678            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
679&                                        shuffle=0,     &
680&                                        deflate=1,     &
681&                                        deflate_level=DEF_LEVEL)
682            ncret = nf90_put_var(ncid, ncvarid, &
683&                                msl(0:nxmax-1, 0:nymax-1, 1, cm_index))
684
685            ncret = nf90_def_var(ncid, 'tcc', NF90_FLOAT, &
686&                                       dim2dids, ncvarid)
687            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
688&                                        shuffle=0,     &
689&                                        deflate=1,     &
690&                                        deflate_level=DEF_LEVEL)
691            ncret = nf90_put_var(ncid, ncvarid, &
692&                                tcc(0:nxmax-1, 0:nymax-1, 1, cm_index))
693
694            ncret = nf90_def_var(ncid, 'u10', NF90_FLOAT, &
695&                                       dim2dids, ncvarid)
696            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
697&                                        shuffle=0,     &
698&                                        deflate=1,     &
699&                                        deflate_level=DEF_LEVEL)
700            ncret = nf90_put_var(ncid, ncvarid, &
701&                                u10(0:nxmax-1, 0:nymax-1, 1, cm_index))
702
703            ncret = nf90_def_var(ncid, 'v10', NF90_FLOAT, &
704&                                       dim2dids, ncvarid)
705            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
706&                                        shuffle=0,     &
707&                                        deflate=1,     &
708&                                        deflate_level=DEF_LEVEL)
709            ncret = nf90_put_var(ncid, ncvarid, &
710&                                v10(0:nxmax-1, 0:nymax-1, 1, cm_index))
711
712            ncret = nf90_def_var(ncid, 'tt2', NF90_FLOAT, &
713&                                       dim2dids, ncvarid)
714            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
715&                                        shuffle=0,     &
716&                                        deflate=1,     &
717&                                        deflate_level=DEF_LEVEL)
718            ncret = nf90_put_var(ncid, ncvarid, &
719&                                tt2(0:nxmax-1, 0:nymax-1, 1, cm_index))
720
721            ncret = nf90_def_var(ncid, 'td2', NF90_FLOAT, &
722&                                       dim2dids, ncvarid)
723            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
724&                                        shuffle=0,     &
725&                                        deflate=1,     &
726&                                        deflate_level=DEF_LEVEL)
727            ncret = nf90_put_var(ncid, ncvarid, &
728&                                td2(0:nxmax-1, 0:nymax-1, 1, cm_index))
729
730            ncret = nf90_def_var(ncid, 'lsprec', NF90_FLOAT, &
731&                                       dim2dids, ncvarid)
732            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
733&                                        shuffle=0,     &
734&                                        deflate=1,     &
735&                                        deflate_level=DEF_LEVEL)
736            ncret = nf90_put_var(ncid, ncvarid, &
737&                                lsprec(0:nxmax-1, 0:nymax-1, 1, cm_index))
738
739            ncret = nf90_def_var(ncid, 'convprec', NF90_FLOAT, &
740&                                       dim2dids, ncvarid)
741            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
742&                                        shuffle=0,     &
743&                                        deflate=1,     &
744&                                        deflate_level=DEF_LEVEL)
745            ncret = nf90_put_var(ncid, ncvarid, &
746&                                convprec(0:nxmax-1, 0:nymax-1, 1, cm_index))
747
748            ncret = nf90_def_var(ncid, 'sshf', NF90_FLOAT, &
749&                                       dim2dids, ncvarid)
750            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
751&                                        shuffle=0,     &
752&                                        deflate=1,     &
753&                                        deflate_level=DEF_LEVEL)
754            ncret = nf90_put_var(ncid, ncvarid, &
755&                                sshf(0:nxmax-1, 0:nymax-1, 1, cm_index))
756
757            ncret = nf90_def_var(ncid, 'ssr', NF90_FLOAT, &
758&                                       dim2dids, ncvarid)
759            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
760&                                        shuffle=0,     &
761&                                        deflate=1,     &
762&                                        deflate_level=DEF_LEVEL)
763            ncret = nf90_put_var(ncid, ncvarid, &
764&                                ssr(0:nxmax-1, 0:nymax-1, 1, cm_index))
765
766            ncret = nf90_def_var(ncid, 'surfstr', NF90_FLOAT, &
767&                                       dim2dids, ncvarid)
768            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
769&                                        shuffle=0,     &
770&                                        deflate=1,     &
771&                                        deflate_level=DEF_LEVEL)
772            ncret = nf90_put_var(ncid, ncvarid, &
773&                                surfstr(0:nxmax-1, 0:nymax-1, 1, cm_index))
774
775            ncret = nf90_def_var(ncid, 'ustar', NF90_FLOAT, &
776&                                       dim2dids, ncvarid)
777            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
778&                                        shuffle=0,     &
779&                                        deflate=1,     &
780&                                        deflate_level=DEF_LEVEL)
781            ncret = nf90_put_var(ncid, ncvarid, &
782&                                ustar(0:nxmax-1, 0:nymax-1, 1, cm_index))
783
784            ncret = nf90_def_var(ncid, 'wstar', NF90_FLOAT, &
785&                                       dim2dids, ncvarid)
786            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
787&                                        shuffle=0,     &
788&                                        deflate=1,     &
789&                                        deflate_level=DEF_LEVEL)
790            ncret = nf90_put_var(ncid, ncvarid, &
791&                                wstar(0:nxmax-1, 0:nymax-1, 1, cm_index))
792
793            ncret = nf90_def_var(ncid, 'hmix', NF90_FLOAT, &
794&                                       dim2dids, ncvarid)
795            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
796&                                        shuffle=0,     &
797&                                        deflate=1,     &
798&                                        deflate_level=DEF_LEVEL)
799            ncret = nf90_put_var(ncid, ncvarid, &
800&                                hmix(0:nxmax-1, 0:nymax-1, 1, cm_index))
801
802            ncret = nf90_def_var(ncid, 'tropopause', NF90_FLOAT, &
803&                                       dim2dids, ncvarid)
804            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
805&                                        shuffle=0,     &
806&                                        deflate=1,     &
807&                                        deflate_level=DEF_LEVEL)
808            ncret = nf90_put_var(ncid, ncvarid, &
809&                                tropopause(0:nxmax-1, 0:nymax-1, 1, cm_index))
810
811            ncret = nf90_def_var(ncid, 'oli', NF90_FLOAT, &
812&                                       dim2dids, ncvarid)
813            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
814&                                        shuffle=0,     &
815&                                        deflate=1,     &
816&                                        deflate_level=DEF_LEVEL)
817            ncret = nf90_put_var(ncid, ncvarid, &
818&                                oli(0:nxmax-1, 0:nymax-1, 1, cm_index))
819
820            ncret = nf90_def_var(ncid, 'diffk', NF90_FLOAT, &
821&                                       dim2dids, ncvarid)
822            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
823&                                        shuffle=0,     &
824&                                        deflate=1,     &
825&                                        deflate_level=DEF_LEVEL)
826            ncret = nf90_put_var(ncid, ncvarid, &
827&                                diffk(0:nxmax-1, 0:nymax-1, 1, cm_index))
828
829
830
831            PRINT *, 'SUM(ps(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
832&                                        SUM(ps(0:nxmax-1,0:nymax-1,1, cm_index))
833
834            PRINT *, 'SUM(wstar(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
835&                                        SUM(wstar(0:nxmax-1,0:nymax-1,1, cm_index))
836
837
838            dim3dids = (/nxmax_dimid, nymax_dimid, maxspec/)
839
840            ncret = nf90_def_var(ncid, 'vdep', NF90_FLOAT, &
841&                                       dim3dids, ncvarid)
842            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
843&                                        shuffle=0,     &
844&                                        deflate=1,     &
845&                                        deflate_level=DEF_LEVEL)
846            ncret = nf90_put_var(ncid, ncvarid, &
847&                                vdep(0:nxmax-1, 0:nymax-1, 1:maxspec, cm_index))
848
849
850
851            ! 1d fields
852            WRITE(iounit) z0(:)
853            WRITE(iounit) akm(:)
854            WRITE(iounit) bkm(:)
855            WRITE(iounit) akz(:)
856            WRITE(iounit) bkz(:)
857            WRITE(iounit) aknew(:)
858            WRITE(iounit) bknew(:)
859
860            dim1dids = (/numclass/)
861
862            ncret = nf90_def_var(ncid, 'z0', NF90_FLOAT, &
863&                                       dim1dids, ncvarid)
864            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
865&                                        shuffle=0,     &
866&                                        deflate=1,     &
867&                                        deflate_level=DEF_LEVEL)
868            ncret = nf90_put_var(ncid, ncvarid, &
869&                                z0(1:numclass))
870
871
872            dim1dids = (/nwzmax/)
873
874            ncret = nf90_def_var(ncid, 'akm', NF90_FLOAT, &
875&                                       dim1dids, ncvarid)
876            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
877&                                        shuffle=0,     &
878&                                        deflate=1,     &
879&                                        deflate_level=DEF_LEVEL)
880            ncret = nf90_put_var(ncid, ncvarid, &
881&                                akm(1:nwzmax))
882
883            ncret = nf90_def_var(ncid, 'bkm', NF90_FLOAT, &
884&                                       dim1dids, ncvarid)
885            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
886&                                        shuffle=0,     &
887&                                        deflate=1,     &
888&                                        deflate_level=DEF_LEVEL)
889            ncret = nf90_put_var(ncid, ncvarid, &
890&                                bkm(1:nwzmax))
891
892
893            dim1dids = (/nuvzmax/)
894
895            ncret = nf90_def_var(ncid, 'akz', NF90_FLOAT, &
896&                                       dim1dids, ncvarid)
897            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
898&                                        shuffle=0,     &
899&                                        deflate=1,     &
900&                                        deflate_level=DEF_LEVEL)
901            ncret = nf90_put_var(ncid, ncvarid, &
902&                                akz(1:nuvzmax))
903
904            ncret = nf90_def_var(ncid, 'bkz', NF90_FLOAT, &
905&                                       dim1dids, ncvarid)
906            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
907&                                        shuffle=0,     &
908&                                        deflate=1,     &
909&                                        deflate_level=DEF_LEVEL)
910            ncret = nf90_put_var(ncid, ncvarid, &
911&                                bkz(1:nuvzmax))
912
913
914            dim1dids = (/nzmax/)
915
916            ncret = nf90_def_var(ncid, 'aknew', NF90_FLOAT, &
917&                                       dim1dids, ncvarid)
918            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
919&                                        shuffle=0,     &
920&                                        deflate=1,     &
921&                                        deflate_level=DEF_LEVEL)
922            ncret = nf90_put_var(ncid, ncvarid, &
923&                                aknew(1:nzmax))
924
925            ncret = nf90_def_var(ncid, 'bknew', NF90_FLOAT, &
926&                                       dim1dids, ncvarid)
927            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
928&                                        shuffle=0,     &
929&                                        deflate=1,     &
930&                                        deflate_level=DEF_LEVEL)
931            ncret = nf90_put_var(ncid, ncvarid, &
932&                                bknew(1:nzmax))
933
934
935            PRINT *, 'SUM(bknew(1:nzmax)): ', &
936&                                        SUM(bknew(1:nzmax))
937
938
939
940            ! Getting ready to add in nested code
941
942            ! These are compiled max dimensions from par_mod - these are
943            ! not meant to be reassigned during a LOAD, but used as "header"
944            ! information to provide the structure of arrays
945            ncret = nf90_def_dim(ncid, 'maxnests', maxnests, maxnests_dimid)
946            ncret = nf90_def_dim(ncid, 'nxmaxn', nxmaxn, nxmaxn_dimid)
947            ncret = nf90_def_dim(ncid, 'nymaxn', nymaxn, nymaxn_dimid)
948
949
950
951            ! Nested, scalar values (for each nest)
952            WRITE(iounit) nxn(:)
953            WRITE(iounit) nyn(:)
954            WRITE(iounit) dxn(:)
955            WRITE(iounit) dyn(:)
956            WRITE(iounit) xlon0n(:)
957            WRITE(iounit) ylat0n(:)
958
959            ! Nested fields, static over time
960            WRITE(iounit) oron, excessoron, lsmn, xlandusen
961
962            ! 3d nested fields
963            WRITE(iounit) uun(:,:,:,cm_index,:)
964            WRITE(iounit) vvn(:,:,:,cm_index,:)
965            WRITE(iounit) wwn(:,:,:,cm_index,:)
966            WRITE(iounit) ttn(:,:,:,cm_index,:)
967            WRITE(iounit) qvn(:,:,:,cm_index,:)
968            WRITE(iounit) pvn(:,:,:,cm_index,:)
969            WRITE(iounit) cloudsn(:,:,:,cm_index,:)
970            WRITE(iounit) cloudsnh(:,:,cm_index,:)
971            WRITE(iounit) rhon(:,:,:,cm_index,:)
972            WRITE(iounit) drhodzn(:,:,:,cm_index,:)
973            WRITE(iounit) tthn(:,:,:,cm_index,:)
974            WRITE(iounit) qvhn(:,:,:,cm_index,:)
975
976            ! 2d nested fields
977            WRITE(iounit) psn(:,:,:,cm_index,:)
978            WRITE(iounit) sdn(:,:,:,cm_index,:)
979            WRITE(iounit) msln(:,:,:,cm_index,:)
980            WRITE(iounit) tccn(:,:,:,cm_index,:)
981            WRITE(iounit) u10n(:,:,:,cm_index,:)
982            WRITE(iounit) v10n(:,:,:,cm_index,:)
983            WRITE(iounit) tt2n(:,:,:,cm_index,:)
984            WRITE(iounit) td2n(:,:,:,cm_index,:)
985            WRITE(iounit) lsprecn(:,:,:,cm_index,:)
986            WRITE(iounit) convprecn(:,:,:,cm_index,:)
987            WRITE(iounit) sshfn(:,:,:,cm_index,:)
988            WRITE(iounit) ssrn(:,:,:,cm_index,:)
989            WRITE(iounit) surfstrn(:,:,:,cm_index,:)
990            WRITE(iounit) ustarn(:,:,:,cm_index,:)
991            WRITE(iounit) wstarn(:,:,:,cm_index,:)
992            WRITE(iounit) hmixn(:,:,:,cm_index,:)
993            WRITE(iounit) tropopausen(:,:,:,cm_index,:)
994            WRITE(iounit) olin(:,:,:,cm_index,:)
995            WRITE(iounit) diffkn(:,:,:,cm_index,:)
996            WRITE(iounit) vdepn(:,:,:,cm_index,:)
997
998            ! Auxiliary variables for nests
999            WRITE(iounit) xresoln(:)
1000            WRITE(iounit) yresoln(:)
1001            WRITE(iounit) xln(:)
1002            WRITE(iounit) yln(:)
1003            WRITE(iounit) xrn(:)
1004            WRITE(iounit) yrn(:)
1005
1006            ! Variables for polar stereographic projection
1007            WRITE(iounit) xglobal, sglobal, nglobal
1008            WRITE(iounit) switchnorthg, switchsouthg
1009            WRITE(iounit) southpolemap(:)
1010            WRITE(iounit) northpolemap(:)
1011
1012            ! Variables declared in conv_mod (convection)
1013            WRITE(iounit) pconv(:)
1014            WRITE(iounit) phconv(:)
1015            WRITE(iounit) dpr(:)
1016            WRITE(iounit) pconv_hpa(:)
1017            WRITE(iounit) phconv_hpa(:)
1018            WRITE(iounit) ft(:)
1019            WRITE(iounit) fq(:)
1020            WRITE(iounit) fmass(:,:)
1021            WRITE(iounit) sub(:)
1022            WRITE(iounit) fmassfrac(:,:)
1023            WRITE(iounit) cbaseflux(:,:)
1024            WRITE(iounit) cbasefluxn(:,:,:)
1025            WRITE(iounit) tconv(:)
1026            WRITE(iounit) qconv(:)
1027            WRITE(iounit) qsconv(:)
1028            WRITE(iounit) psconv, tt2conv, td2conv
1029            WRITE(iounit) nconvlev, nconvtop
1030
1031        ELSE IF (op == 'LOAD') THEN
1032
1033            ! Read the preprocessed format version string and insure it
1034            ! matches this version
1035            READ (iounit) temp_preproc_format_version_str
1036            PRINT *, 'Reading preprocessed file format version: ', &
1037&                    temp_preproc_format_version_str
1038
1039            IF (TRIM(temp_preproc_format_version_str) == &
1040&                                        TRIM(PREPROC_FORMAT_VERSION_STR)) THEN
1041                CONTINUE
1042            ELSE
1043                ! PRINT *, ''  GK: causes relocation truncated to fit: R_X86_64_32
1044                PRINT *, 'Inconsistent preprocessing format version'
1045                PRINT *, 'Expected Version: ', PREPROC_FORMAT_VERSION_STR
1046                PRINT *, 'Detected Version: ', temp_preproc_format_version_str
1047                ! PRINT *, ''
1048                STOP
1049            END IF
1050
1051            ! Read the compiled max dimensions that were dumped from par_mod
1052            ! when creating the fp file, so that we can compare against
1053            ! current FLEXPART dimensions - they need to be the same, or else
1054            ! we abort.
1055            READ (iounit) temp_nxmax, temp_nymax, temp_nzmax, &
1056&                         temp_nuvzmax, temp_nwzmax
1057
1058            ! Get dimensions
1059            ncret = nf90_inq_dimid(ncid, 'nxmax', nxmax_dimid)
1060            ncret = nf90_inquire_dimension(ncid, nxmax_dimid, nxmax_dimname, &
1061&                                                temp_nxmax)
1062            PRINT *, 'temp_nxmax: ', temp_nxmax
1063
1064            ncret = nf90_inq_dimid(ncid, 'nymax', nymax_dimid)
1065            ncret = nf90_inquire_dimension(ncid, nymax_dimid, nymax_dimname, &
1066&                                                temp_nymax)
1067            PRINT *, 'temp_nymax: ', temp_nymax
1068
1069            ncret = nf90_inq_dimid(ncid, 'nzmax', nzmax_dimid)
1070            ncret = nf90_inquire_dimension(ncid, nzmax_dimid, nzmax_dimname, &
1071&                                                temp_nzmax)
1072            PRINT *, 'temp_nzmax: ', temp_nzmax
1073
1074            ncret = nf90_inq_dimid(ncid, 'nuvzmax', nuvzmax_dimid)
1075            ncret = nf90_inquire_dimension(ncid, nuvzmax_dimid, nuvzmax_dimname, &
1076&                                                temp_nuvzmax)
1077            PRINT *, 'temp_nuvzmax: ', temp_nuvzmax
1078
1079            ncret = nf90_inq_dimid(ncid, 'nwzmax', nwzmax_dimid)
1080            ncret = nf90_inquire_dimension(ncid, nwzmax_dimid, nwzmax_dimname, &
1081&                                                temp_nwzmax)
1082            PRINT *, 'temp_nwzmax: ', temp_nwzmax
1083
1084
1085
1086            IF ( (temp_nxmax == nxmax) .AND. (temp_nymax == nymax) .AND. &
1087&                   (temp_nzmax == nzmax) .AND. &
1088&                   (temp_nuvzmax == nuvzmax) .AND. &
1089&                   (temp_nwzmax == nwzmax) ) THEN
1090                CONTINUE
1091            ELSE
1092                PRINT *, 'Incompatible dimensions between fp file and current FLEXPART!'
1093                ! PRINT *, ''
1094                PRINT *, '                  FP file     Compiled FP'
1095                PRINT *, 'nxmax:     ', temp_nxmax, '    ', nxmax
1096                PRINT *, 'nymax:     ', temp_nymax, '    ', nymax
1097                PRINT *, 'nzmax:     ', temp_nzmax, '    ', nzmax
1098                PRINT *, 'nuvzmax:     ', temp_nuvzmax, '    ', nuvzmax
1099                PRINT *, 'nwzmax:     ', temp_nwzmax, '    ', nwzmax
1100                ! PRINT *, ''
1101                STOP
1102            END IF
1103
1104
1105
1106
1107            ! Scalar values
1108            READ(iounit) nx, ny, nxmin1, nymin1, nxfield
1109            READ(iounit) nuvz, nwz, nz, nmixz, nlev_ec
1110            READ(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst
1111
1112
1113
1114            ! Get the varid , then read into scalar variable
1115            ncret = nf90_inq_varid(ncid, 'nx', ncvarid)
1116            ncret = nf90_get_var(ncid, ncvarid, nx)
1117
1118            ncret = nf90_inq_varid(ncid, 'ny', ncvarid)
1119            ncret = nf90_get_var(ncid, ncvarid, ny)
1120
1121            ncret = nf90_inq_varid(ncid, 'nxmin1', ncvarid)
1122            ncret = nf90_get_var(ncid, ncvarid, nxmin1)
1123
1124            ncret = nf90_inq_varid(ncid, 'nymin1', ncvarid)
1125            ncret = nf90_get_var(ncid, ncvarid, nymin1)
1126
1127            ncret = nf90_inq_varid(ncid, 'nxfield', ncvarid)
1128            ncret = nf90_get_var(ncid, ncvarid, nxfield)
1129
1130            ncret = nf90_inq_varid(ncid, 'nuvz', ncvarid)
1131            ncret = nf90_get_var(ncid, ncvarid, nuvz)
1132
1133            ncret = nf90_inq_varid(ncid, 'nwz', ncvarid)
1134            ncret = nf90_get_var(ncid, ncvarid, nwz)
1135
1136            ncret = nf90_inq_varid(ncid, 'nz', ncvarid)
1137            ncret = nf90_get_var(ncid, ncvarid, nz)
1138
1139            ncret = nf90_inq_varid(ncid, 'nmixz', ncvarid)
1140            ncret = nf90_get_var(ncid, ncvarid, nmixz)
1141
1142            ncret = nf90_inq_varid(ncid, 'nlev_ec', ncvarid)
1143            ncret = nf90_get_var(ncid, ncvarid, nlev_ec)
1144
1145            ncret = nf90_inq_varid(ncid, 'dx', ncvarid)
1146            ncret = nf90_get_var(ncid, ncvarid, dx)
1147
1148            ncret = nf90_inq_varid(ncid, 'dy', ncvarid)
1149            ncret = nf90_get_var(ncid, ncvarid, dy)
1150
1151            ncret = nf90_inq_varid(ncid, 'xlon0', ncvarid)
1152            ncret = nf90_get_var(ncid, ncvarid, xlon0)
1153
1154            ncret = nf90_inq_varid(ncid, 'ylat0', ncvarid)
1155            ncret = nf90_get_var(ncid, ncvarid, ylat0)
1156
1157            ncret = nf90_inq_varid(ncid, 'dxconst', ncvarid)
1158            ncret = nf90_get_var(ncid, ncvarid, dxconst)
1159
1160            ncret = nf90_inq_varid(ncid, 'dyconst', ncvarid)
1161            ncret = nf90_get_var(ncid, ncvarid, dyconst)
1162
1163
1164
1165
1166
1167
1168            ! Fixed fields, static in time
1169            READ(iounit) oro, excessoro, lsm, xlanduse, height
1170
1171            ncret = nf90_inq_varid(ncid, 'oro', ncvarid)
1172            ncret = nf90_get_var(ncid, ncvarid, oro(0:nxmax-1,0:nymax-1))
1173
1174            ncret = nf90_inq_varid(ncid, 'excessoro', ncvarid)
1175            ncret = nf90_get_var(ncid, ncvarid, excessoro(0:nxmax-1,0:nymax-1))
1176
1177            ncret = nf90_inq_varid(ncid, 'lsm', ncvarid)
1178            ncret = nf90_get_var(ncid, ncvarid, lsm(0:nxmax-1,0:nymax-1))
1179
1180            ncret = nf90_inq_varid(ncid, 'xlanduse', ncvarid)
1181            ncret = nf90_get_var(ncid, ncvarid, xlanduse(0:nxmax-1,0:nymax-1, 1:numclass))
1182
1183            ncret = nf90_inq_varid(ncid, 'height', ncvarid)
1184            ncret = nf90_get_var(ncid, ncvarid, height(1:nzmax))
1185
1186
1187
1188
1189            ! 3d fields
1190            READ(iounit) uu(:,:,:,cm_index)
1191            READ(iounit) vv(:,:,:,cm_index)
1192            READ(iounit) uupol(:,:,:,cm_index)
1193            READ(iounit) vvpol(:,:,:,cm_index)
1194            READ(iounit) ww(:,:,:,cm_index)
1195            READ(iounit) tt(:,:,:,cm_index)
1196            READ(iounit) qv(:,:,:,cm_index)
1197            READ(iounit) pv(:,:,:,cm_index)
1198            READ(iounit) rho(:,:,:,cm_index)
1199            READ(iounit) drhodz(:,:,:,cm_index)
1200            READ(iounit) tth(:,:,:,cm_index)
1201            READ(iounit) qvh(:,:,:,cm_index)
1202            READ(iounit) pplev(:,:,:,cm_index)
1203            READ(iounit) clouds(:,:,:,cm_index)
1204            READ(iounit) cloudsh(:,:,cm_index)
1205
1206
1207
1208
1209            ! Get the varid and read the variable into the array
1210            ncret = nf90_inq_varid(ncid, 'uu', ncvarid)
1211            ncret = nf90_get_var(ncid, ncvarid, uu(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1212
1213            ncret = nf90_inq_varid(ncid, 'vv', ncvarid)
1214            ncret = nf90_get_var(ncid, ncvarid, vv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1215
1216            ncret = nf90_inq_varid(ncid, 'uupol', ncvarid)
1217            ncret = nf90_get_var(ncid, ncvarid, uupol(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1218
1219            ncret = nf90_inq_varid(ncid, 'vvpol', ncvarid)
1220            ncret = nf90_get_var(ncid, ncvarid, vvpol(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1221
1222            ncret = nf90_inq_varid(ncid, 'ww', ncvarid)
1223            ncret = nf90_get_var(ncid, ncvarid, ww(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1224
1225            ncret = nf90_inq_varid(ncid, 'tt', ncvarid)
1226            ncret = nf90_get_var(ncid, ncvarid, tt(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1227
1228            ncret = nf90_inq_varid(ncid, 'qv', ncvarid)
1229            ncret = nf90_get_var(ncid, ncvarid, qv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1230
1231            ncret = nf90_inq_varid(ncid, 'pv', ncvarid)
1232            ncret = nf90_get_var(ncid, ncvarid, pv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1233
1234            ncret = nf90_inq_varid(ncid, 'rho', ncvarid)
1235            ncret = nf90_get_var(ncid, ncvarid, rho(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1236
1237            ncret = nf90_inq_varid(ncid, 'drhodz', ncvarid)
1238            ncret = nf90_get_var(ncid, ncvarid, drhodz(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
1239
1240            ncret = nf90_inq_varid(ncid, 'clouds', ncvarid)
1241            ncret = nf90_get_var(ncid, ncvarid, clouds(0:nxmax-1,0:nzmax-1,1:nzmax,cm_index))
1242
1243            ncret = nf90_inq_varid(ncid, 'tth', ncvarid)
1244            ncret = nf90_get_var(ncid, ncvarid, tth(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
1245
1246            ncret = nf90_inq_varid(ncid, 'qvh', ncvarid)
1247            ncret = nf90_get_var(ncid, ncvarid, qvh(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
1248
1249            ncret = nf90_inq_varid(ncid, 'pplev', ncvarid)
1250            ncret = nf90_get_var(ncid, ncvarid, pplev(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
1251
1252            ncret = nf90_inq_varid(ncid, 'cloudsh', ncvarid)
1253            ncret = nf90_get_var(ncid, ncvarid, cloudsh(0:nxmax-1,0:nymax-1,cm_index))
1254
1255
1256
1257
1258            PRINT *, 'SUM(tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
1259&                                        SUM(tt(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
1260
1261
1262            PRINT *, 'SUM(clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
1263&                                        SUM(clouds(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
1264
1265
1266
1267            ! 2d fields
1268            READ(iounit) ps(:,:,:,cm_index)
1269            READ(iounit) sd(:,:,:,cm_index)
1270            READ(iounit) msl(:,:,:,cm_index)
1271            READ(iounit) tcc(:,:,:,cm_index)
1272            READ(iounit) u10(:,:,:,cm_index)
1273            READ(iounit) v10(:,:,:,cm_index)
1274            READ(iounit) tt2(:,:,:,cm_index)
1275            READ(iounit) td2(:,:,:,cm_index)
1276            READ(iounit) lsprec(:,:,:,cm_index)
1277            READ(iounit) convprec(:,:,:,cm_index)
1278            READ(iounit) sshf(:,:,:,cm_index)
1279            READ(iounit) ssr(:,:,:,cm_index)
1280            READ(iounit) surfstr(:,:,:,cm_index)
1281            READ(iounit) ustar(:,:,:,cm_index)
1282            READ(iounit) wstar(:,:,:,cm_index)
1283            READ(iounit) hmix(:,:,:,cm_index)
1284            READ(iounit) tropopause(:,:,:,cm_index)
1285            READ(iounit) oli(:,:,:,cm_index)
1286            READ(iounit) diffk(:,:,:,cm_index)
1287            READ(iounit) vdep(:,:,:,cm_index)
1288
1289            ! Get the varid and read the variable into the array
1290            ncret = nf90_inq_varid(ncid, 'ps', ncvarid)
1291            ncret = nf90_get_var(ncid, ncvarid, ps(0:nxmax-1,0:nymax-1,1,cm_index))
1292
1293            ncret = nf90_inq_varid(ncid, 'sd', ncvarid)
1294            ncret = nf90_get_var(ncid, ncvarid, sd(0:nxmax-1,0:nymax-1,1,cm_index))
1295
1296            ncret = nf90_inq_varid(ncid, 'msl', ncvarid)
1297            ncret = nf90_get_var(ncid, ncvarid, msl(0:nxmax-1,0:nymax-1,1,cm_index))
1298
1299            ncret = nf90_inq_varid(ncid, 'tcc', ncvarid)
1300            ncret = nf90_get_var(ncid, ncvarid, tcc(0:nxmax-1,0:nymax-1,1,cm_index))
1301
1302            ncret = nf90_inq_varid(ncid, 'u10', ncvarid)
1303            ncret = nf90_get_var(ncid, ncvarid, u10(0:nxmax-1,0:nymax-1,1,cm_index))
1304
1305            ncret = nf90_inq_varid(ncid, 'v10', ncvarid)
1306            ncret = nf90_get_var(ncid, ncvarid, v10(0:nxmax-1,0:nymax-1,1,cm_index))
1307
1308            ncret = nf90_inq_varid(ncid, 'tt2', ncvarid)
1309            ncret = nf90_get_var(ncid, ncvarid, tt2(0:nxmax-1,0:nymax-1,1,cm_index))
1310
1311            ncret = nf90_inq_varid(ncid, 'td2', ncvarid)
1312            ncret = nf90_get_var(ncid, ncvarid, td2(0:nxmax-1,0:nymax-1,1,cm_index))
1313
1314            ncret = nf90_inq_varid(ncid, 'lsprec', ncvarid)
1315            ncret = nf90_get_var(ncid, ncvarid, lsprec(0:nxmax-1,0:nymax-1,1,cm_index))
1316
1317            ncret = nf90_inq_varid(ncid, 'convprec', ncvarid)
1318            ncret = nf90_get_var(ncid, ncvarid, convprec(0:nxmax-1,0:nymax-1,1,cm_index))
1319
1320            ncret = nf90_inq_varid(ncid, 'sshf', ncvarid)
1321            ncret = nf90_get_var(ncid, ncvarid, sshf(0:nxmax-1,0:nymax-1,1,cm_index))
1322
1323            ncret = nf90_inq_varid(ncid, 'ssr', ncvarid)
1324            ncret = nf90_get_var(ncid, ncvarid, ssr(0:nxmax-1,0:nymax-1,1,cm_index))
1325
1326            ncret = nf90_inq_varid(ncid, 'surfstr', ncvarid)
1327            ncret = nf90_get_var(ncid, ncvarid, surfstr(0:nxmax-1,0:nymax-1,1,cm_index))
1328
1329            ncret = nf90_inq_varid(ncid, 'ustar', ncvarid)
1330            ncret = nf90_get_var(ncid, ncvarid, ustar(0:nxmax-1,0:nymax-1,1,cm_index))
1331
1332            ncret = nf90_inq_varid(ncid, 'wstar', ncvarid)
1333            ncret = nf90_get_var(ncid, ncvarid, wstar(0:nxmax-1,0:nymax-1,1,cm_index))
1334
1335            ncret = nf90_inq_varid(ncid, 'hmix', ncvarid)
1336            ncret = nf90_get_var(ncid, ncvarid, hmix(0:nxmax-1,0:nymax-1,1,cm_index))
1337
1338            ncret = nf90_inq_varid(ncid, 'tropopause', ncvarid)
1339            ncret = nf90_get_var(ncid, ncvarid, tropopause(0:nxmax-1,0:nymax-1,1,cm_index))
1340
1341            ncret = nf90_inq_varid(ncid, 'oli', ncvarid)
1342            ncret = nf90_get_var(ncid, ncvarid, oli(0:nxmax-1,0:nymax-1,1,cm_index))
1343
1344            ncret = nf90_inq_varid(ncid, 'diffk', ncvarid)
1345            ncret = nf90_get_var(ncid, ncvarid, diffk(0:nxmax-1,0:nymax-1,1,cm_index))
1346
1347
1348
1349
1350            PRINT *, 'SUM(ps(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
1351&                                        SUM(ps(0:nxmax-1,0:nymax-1,1, cm_index))
1352
1353            PRINT *, 'SUM(wstar(0:nxmax-1, 0:nymax-1, 1, cm_index)): ', &
1354&                                        SUM(wstar(0:nxmax-1,0:nymax-1,1, cm_index))
1355
1356
1357            ncret = nf90_inq_varid(ncid, 'vdep', ncvarid)
1358            ncret = nf90_get_var(ncid, ncvarid, vdep(0:nxmax-1,0:nymax-1,1:maxspec, cm_index))
1359
1360
1361
1362
1363
1364            ! 1d fields
1365            READ(iounit) z0(:)
1366            READ(iounit) akm(:)
1367            READ(iounit) bkm(:)
1368            READ(iounit) akz(:)
1369            READ(iounit) bkz(:)
1370            READ(iounit) aknew(:)
1371            READ(iounit) bknew(:)
1372
1373
1374            ncret = nf90_inq_varid(ncid, 'z0', ncvarid)
1375            ncret = nf90_get_var(ncid, ncvarid, z0(1:numclass))
1376
1377            ncret = nf90_inq_varid(ncid, 'akm', ncvarid)
1378            ncret = nf90_get_var(ncid, ncvarid, akm(1:nwzmax))
1379
1380            ncret = nf90_inq_varid(ncid, 'bkm', ncvarid)
1381            ncret = nf90_get_var(ncid, ncvarid, bkm(1:nwzmax))
1382
1383            ncret = nf90_inq_varid(ncid, 'akz', ncvarid)
1384            ncret = nf90_get_var(ncid, ncvarid, akz(1:nuvzmax))
1385
1386            ncret = nf90_inq_varid(ncid, 'bkz', ncvarid)
1387            ncret = nf90_get_var(ncid, ncvarid, bkz(1:nuvzmax))
1388
1389            ncret = nf90_inq_varid(ncid, 'aknew', ncvarid)
1390            ncret = nf90_get_var(ncid, ncvarid, aknew(1:nzmax))
1391
1392            ncret = nf90_inq_varid(ncid, 'bknew', ncvarid)
1393            ncret = nf90_get_var(ncid, ncvarid, bknew(1:nzmax))
1394
1395
1396
1397            PRINT *, 'SUM(bknew(1:nzmax)): ', &
1398&                                        SUM(bknew(1:nzmax))
1399
1400
1401
1402
1403            ! Now the nested input grid variables
1404            ! Get the compiled values that were written into the FP file, and
1405            ! make sure they are equal to the current compiled values, to make
1406            ! sure we are working with consistent arrays
1407            ncret = nf90_inq_dimid(ncid, 'maxnests', maxnests_dimid)
1408            ncret = nf90_inquire_dimension(ncid, maxnests_dimid, maxnests_dimname, &
1409&                                                temp_maxnests)
1410            PRINT *, 'temp_maxnests: ', temp_maxnests
1411
1412            ncret = nf90_inq_dimid(ncid, 'nxmaxn', nxmaxn_dimid)
1413            ncret = nf90_inquire_dimension(ncid, nxmaxn_dimid, nxmaxn_dimname, &
1414&                                                temp_nxmaxn)
1415            PRINT *, 'temp_nxmaxn: ', temp_nxmaxn
1416
1417            ncret = nf90_inq_dimid(ncid, 'nymaxn', nymaxn_dimid)
1418            ncret = nf90_inquire_dimension(ncid, nymaxn_dimid, nymaxn_dimname, &
1419&                                                temp_nymaxn)
1420            PRINT *, 'temp_nymaxn: ', temp_nymaxn
1421
1422
1423            IF ( (temp_nxmaxn == nxmaxn) .AND. (temp_nymaxn == nymaxn) .AND. &
1424&                   (temp_maxnests == maxnests) ) THEN
1425                CONTINUE
1426            ELSE
1427                PRINT *, 'Incompatible dimensions between fp file and current FLEXPART!'
1428                ! PRINT *, ''
1429                PRINT *, '                  FP file     Compiled FP'
1430                PRINT *, 'nxmaxn:     ', temp_nxmaxn, '    ', nxmaxn
1431                PRINT *, 'nymaxn:     ', temp_nymaxn, '    ', nymaxn
1432                PRINT *, 'maxnests:     ', temp_maxnests, '    ', maxnests
1433                STOP
1434            END IF
1435
1436
1437
1438            ! Nested, scalar values (for each nest)
1439            READ(iounit) nxn(:)
1440            READ(iounit) nyn(:)
1441            READ(iounit) dxn(:)
1442            READ(iounit) dyn(:)
1443            READ(iounit) xlon0n(:)
1444            READ(iounit) ylat0n(:)
1445
1446
1447            ! Nested fields, static over time
1448            READ(iounit) oron, excessoron, lsmn, xlandusen
1449
1450            ! 3d nested fields
1451            READ(iounit) uun(:,:,:,cm_index,:)
1452            READ(iounit) vvn(:,:,:,cm_index,:)
1453            READ(iounit) wwn(:,:,:,cm_index,:)
1454            READ(iounit) ttn(:,:,:,cm_index,:)
1455            READ(iounit) qvn(:,:,:,cm_index,:)
1456            READ(iounit) pvn(:,:,:,cm_index,:)
1457            READ(iounit) cloudsn(:,:,:,cm_index,:)
1458            READ(iounit) cloudsnh(:,:,cm_index,:)
1459            READ(iounit) rhon(:,:,:,cm_index,:)
1460            READ(iounit) drhodzn(:,:,:,cm_index,:)
1461            READ(iounit) tthn(:,:,:,cm_index,:)
1462            READ(iounit) qvhn(:,:,:,cm_index,:)
1463
1464            ! 2d nested fields
1465            READ(iounit) psn(:,:,:,cm_index,:)
1466            READ(iounit) sdn(:,:,:,cm_index,:)
1467            READ(iounit) msln(:,:,:,cm_index,:)
1468            READ(iounit) tccn(:,:,:,cm_index,:)
1469            READ(iounit) u10n(:,:,:,cm_index,:)
1470            READ(iounit) v10n(:,:,:,cm_index,:)
1471            READ(iounit) tt2n(:,:,:,cm_index,:)
1472            READ(iounit) td2n(:,:,:,cm_index,:)
1473            READ(iounit) lsprecn(:,:,:,cm_index,:)
1474            READ(iounit) convprecn(:,:,:,cm_index,:)
1475            READ(iounit) sshfn(:,:,:,cm_index,:)
1476            READ(iounit) ssrn(:,:,:,cm_index,:)
1477            READ(iounit) surfstrn(:,:,:,cm_index,:)
1478            READ(iounit) ustarn(:,:,:,cm_index,:)
1479            READ(iounit) wstarn(:,:,:,cm_index,:)
1480            READ(iounit) hmixn(:,:,:,cm_index,:)
1481            READ(iounit) tropopausen(:,:,:,cm_index,:)
1482            READ(iounit) olin(:,:,:,cm_index,:)
1483            READ(iounit) diffkn(:,:,:,cm_index,:)
1484            READ(iounit) vdepn(:,:,:,cm_index,:)
1485
1486            ! Auxiliary variables for nests
1487            READ(iounit) xresoln(:)
1488            READ(iounit) yresoln(:)
1489            READ(iounit) xln(:)
1490            READ(iounit) yln(:)
1491            READ(iounit) xrn(:)
1492            READ(iounit) yrn(:)
1493
1494            ! Variables for polar stereographic projection
1495            READ(iounit) xglobal, sglobal, nglobal
1496            READ(iounit) switchnorthg, switchsouthg
1497            READ(iounit) southpolemap(:)
1498            READ(iounit) northpolemap(:)
1499
1500            ! Variables declared in conv_mod (convection)
1501            READ(iounit) pconv(:)
1502            READ(iounit) phconv(:)
1503            READ(iounit) dpr(:)
1504            READ(iounit) pconv_hpa(:)
1505            READ(iounit) phconv_hpa(:)
1506            READ(iounit) ft(:)
1507            READ(iounit) fq(:)
1508            READ(iounit) fmass(:,:)
1509            READ(iounit) sub(:)
1510            READ(iounit) fmassfrac(:,:)
1511            READ(iounit) cbaseflux(:,:)
1512            READ(iounit) cbasefluxn(:,:,:)
1513            READ(iounit) tconv(:)
1514            READ(iounit) qconv(:)
1515            READ(iounit) qsconv(:)
1516            READ(iounit) psconv, tt2conv, td2conv
1517            READ(iounit) nconvlev, nconvtop
1518
1519        ELSE
1520            STOP 'fpio(): Illegal operation'
1521           
1522        ENDIF
1523    END SUBROUTINE fpio
1524
1525    SUBROUTINE fpmetbinary_filetext(filename, cm_index)
1526
1527        ! This is a utility subroutine meant to be used for testing purposes.
1528        ! It facilitates the text output of variables read in from the
1529        ! specified .fp file.  This routine will easily cause the program
1530        ! to crash due memory allocation issues, particularly when you are
1531        ! trying to text print 3d arrays in a single formatted statetment.
1532       
1533        CHARACTER(LEN=*), INTENT(IN) :: filename
1534        INTEGER, INTENT(IN) :: cm_index           ! Index of last dimension in
1535                                                  ! most com_mod variables.
1536                                                  ! Should be 1 or 2
1537
1538        !OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', status='REPLACE', &
1539        !    form="formatted", access="stream")
1540        OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', &
1541             form="formatted", access="APPEND")
1542
1543        WRITE(IOUNIT_TEXTOUT, *) 'oro: ', oro
1544        WRITE(IOUNIT_TEXTOUT, *) 'excessoro: ', excessoro
1545        WRITE(IOUNIT_TEXTOUT, *) 'lsm: ', lsm
1546        WRITE(IOUNIT_TEXTOUT, *) 'xlanduse: ', xlanduse
1547        WRITE(IOUNIT_TEXTOUT, *) 'height: ', height
1548
1549        WRITE(IOUNIT_TEXTOUT, *) 'uu: ', uu(:,:,:,cm_index)
1550        WRITE(IOUNIT_TEXTOUT, *) 'vv: ', vv(:,:,:,cm_index)
1551        WRITE(IOUNIT_TEXTOUT, *) 'uupol: ', uupol(:,:,:,cm_index)
1552        WRITE(IOUNIT_TEXTOUT, *) 'vvpol: ', vvpol(:,:,:,cm_index)
1553        WRITE(IOUNIT_TEXTOUT, *) 'ww: ', ww(:,:,:,cm_index)
1554        WRITE(IOUNIT_TEXTOUT, *) 'tt: ', tt(:,:,:,cm_index)
1555        WRITE(IOUNIT_TEXTOUT, *) 'qv: ', qv(:,:,:,cm_index)
1556        WRITE(IOUNIT_TEXTOUT, *) 'pv: ', pv(:,:,:,cm_index)
1557        WRITE(IOUNIT_TEXTOUT, *) 'rho: ', rho(:,:,:,cm_index)
1558        WRITE(IOUNIT_TEXTOUT, *) 'drhodz: ', drhodz(:,:,:,cm_index)
1559        WRITE(IOUNIT_TEXTOUT, *) 'tth: ', tth(:,:,:,cm_index)
1560        WRITE(IOUNIT_TEXTOUT, *) 'qvh: ', qvh(:,:,:,cm_index)
1561        WRITE(IOUNIT_TEXTOUT, *) 'pplev: ', pplev(:,:,:,cm_index)
1562        WRITE(IOUNIT_TEXTOUT, *) 'clouds: ', clouds(:,:,:,cm_index)
1563        WRITE(IOUNIT_TEXTOUT, *) 'cloudsh: ', cloudsh(:,:,cm_index)
1564
1565
1566
1567
1568        CLOSE(IOUNIT_TEXTOUT)
1569    END SUBROUTINE fpmetbinary_filetext
1570
1571
1572END MODULE fpmetbinary_mod
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG