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

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

Intermediate work testing NC4 formatting of intermediate FP data

  • Property mode set to 100644
File size: 48.0 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
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        ! The lines here correspond to READ and WRITE in the dump/load routines
155
156        ! Scalar values
157        nx=0.0; ny=0.0; nxmin1=0.0; nymin1=0.0; nxfield=0.0
158        nuvz=0.0; nwz=0.0; nz=0.0; nmixz=0.0; nlev_ec=0.0
159        dx=0.0; dy=0.0; xlon0=0.0; ylat0=0.0; dxconst=0.0; dyconst=0.0
160
161        ! Fixed fields, static in time
162        oro=0.0; excessoro=0.0; lsm=0.0; xlanduse=0.0; height=0.0
163
164        ! 3d fields
165        uu(:,:,:,cm_index) = 0.0
166        vv(:,:,:,cm_index) = 0.0
167        uupol(:,:,:,cm_index) = 0.0
168        vvpol(:,:,:,cm_index) = 0.0
169        ww(:,:,:,cm_index) = 0.0
170        tt(:,:,:,cm_index) = 0.0
171        qv(:,:,:,cm_index) = 0.0
172        pv(:,:,:,cm_index) = 0.0
173        rho(:,:,:,cm_index) = 0.0
174        drhodz(:,:,:,cm_index) = 0.0
175        tth(:,:,:,cm_index) = 0.0
176        qvh(:,:,:,cm_index) = 0.0
177        pplev(:,:,:,cm_index) = 0.0
178        clouds(:,:,:,cm_index) = 0.0
179        cloudsh(:,:,cm_index) = 0.0
180
181        ! 2d fields
182        ps(:,:,:,cm_index) = 0.0
183        sd(:,:,:,cm_index) = 0.0
184        msl(:,:,:,cm_index) = 0.0
185        tcc(:,:,:,cm_index) = 0.0
186        u10(:,:,:,cm_index) = 0.0
187        v10(:,:,:,cm_index) = 0.0
188        tt2(:,:,:,cm_index) = 0.0
189        td2(:,:,:,cm_index) = 0.0
190        lsprec(:,:,:,cm_index) = 0.0
191        convprec(:,:,:,cm_index) = 0.0
192        sshf(:,:,:,cm_index) = 0.0
193        ssr(:,:,:,cm_index) = 0.0
194        surfstr(:,:,:,cm_index) = 0.0
195        ustar(:,:,:,cm_index) = 0.0
196        wstar(:,:,:,cm_index) = 0.0
197        hmix(:,:,:,cm_index) = 0.0
198        tropopause(:,:,:,cm_index) = 0.0
199        oli(:,:,:,cm_index) = 0.0
200        diffk(:,:,:,cm_index) = 0.0
201        vdep(:,:,:,cm_index) = 0.0
202
203        ! 1d fields
204        z0(:) = 0.0
205        akm(:) = 0.0
206        bkm(:) = 0.0
207        akz(:) = 0.0
208        bkz(:) = 0.0
209        aknew(:) = 0.0
210        bknew(:) = 0.0
211
212        ! Nested, scalar values (for each nest)
213        nxn(:) = 0.0
214        nyn(:) = 0.0
215        dxn(:) = 0.0
216        dyn(:) = 0.0
217        xlon0n(:) = 0.0
218        ylat0n(:) = 0.0
219
220        ! Nested fields, static in time
221        oron=0.0; excessoron=0.0; lsmn=0.0; xlandusen=0.0
222
223        ! 3d nested fields
224        uun(:,:,:,cm_index,:) = 0.0
225        wwn(:,:,:,cm_index,:) = 0.0
226        ttn(:,:,:,cm_index,:) = 0.0
227        qvn(:,:,:,cm_index,:) = 0.0
228        pvn(:,:,:,cm_index,:) = 0.0
229        cloudsn(:,:,:,cm_index,:) = 0.0
230        cloudsnh(:,:,cm_index,:) = 0.0
231        rhon(:,:,:,cm_index,:) = 0.0
232        drhodzn(:,:,:,cm_index,:) = 0.0
233        tthn(:,:,:,cm_index,:) = 0.0
234        qvhn(:,:,:,cm_index,:) = 0.0
235
236        ! 2d nested fields
237        psn(:,:,:,cm_index,:) = 0.0
238        sdn(:,:,:,cm_index,:) = 0.0
239        msln(:,:,:,cm_index,:) = 0.0
240        tccn(:,:,:,cm_index,:) = 0.0
241        u10n(:,:,:,cm_index,:) = 0.0
242        v10n(:,:,:,cm_index,:) = 0.0
243        tt2n(:,:,:,cm_index,:) = 0.0
244        td2n(:,:,:,cm_index,:) = 0.0
245        lsprecn(:,:,:,cm_index,:) = 0.0
246        convprecn(:,:,:,cm_index,:) = 0.0
247        sshfn(:,:,:,cm_index,:) = 0.0
248        ssrn(:,:,:,cm_index,:) = 0.0
249        surfstrn(:,:,:,cm_index,:) = 0.0
250        ustarn(:,:,:,cm_index,:) = 0.0
251        wstarn(:,:,:,cm_index,:) = 0.0
252        hmixn(:,:,:,cm_index,:) = 0.0
253        tropopausen(:,:,:,cm_index,:) = 0.0
254        olin(:,:,:,cm_index,:) = 0.0
255        diffkn(:,:,:,cm_index,:) = 0.0
256        vdepn(:,:,:,cm_index,:) = 0.0
257
258        ! Auxiliary variables for nests
259        xresoln(:) = 0.0
260        yresoln(:) = 0.0
261        xln(:) = 0.0
262        yln(:) = 0.0
263        xrn(:) = 0.0
264        yrn(:) = 0.0
265
266        ! Variables for polar stereographic projection
267        xglobal=.FALSE.; sglobal=.FALSE.; nglobal=.FALSE.
268        switchnorthg=0.0; switchsouthg=0.0
269        southpolemap(:) = 0.0
270        northpolemap(:) = 0.0
271
272        ! Variables declared in conv_mod (convection)
273        pconv(:) = 0.0
274        phconv(:) = 0.0
275        dpr(:) = 0.0
276        pconv_hpa(:) = 0.0
277        phconv_hpa(:) = 0.0
278        ft(:) = 0.0
279        fq(:) = 0.0
280        fmass(:,:) = 0.0
281        sub(:) = 0.0
282        fmassfrac(:,:) = 0.0
283        cbaseflux(:,:) = 0.0
284        cbasefluxn(:,:,:) = 0.0
285        tconv(:) = 0.0
286        qconv(:) = 0.0
287        qsconv(:) = 0.0
288        psconv=0.0; tt2conv=0.0; td2conv=0.0
289        nconvlev=0.0; nconvtop=0.0
290
291    END SUBROUTINE fpmetbinary_zero
292
293    SUBROUTINE fpio(iounit, ncid, op, cm_index)
294        IMPLICIT NONE
295        INTEGER, INTENT(IN) :: ncid               ! NetCDF file id
296        INTEGER, INTENT(IN) :: iounit
297        CHARACTER(LEN=4), INTENT(IN) :: op        ! Operation - DUMP or LOAD
298        INTEGER, INTENT(IN) :: cm_index           ! Index of last dimension in
299                                                  ! most com_mod variables.
300                                                  ! Should be 1 or 2
301
302
303        INTEGER :: ncret          ! Return value from NetCDF calls
304        INTEGER :: ncvarid          ! NetCDF variable ID
305
306        INTEGER :: nxmax_dimid, nymax_dimid, nzmax_dimid, nuvzmax_dimid, nwzmax_dimid
307
308        INTEGER, DIMENSION(1) :: dim1dids    ! Dimension IDs for 1D arrays
309        INTEGER, DIMENSION(2) :: dim2dids    ! Dimension IDs for 2D arrays
310        INTEGER, DIMENSION(3) :: dim3dids    ! Dimension IDs for 3D arrays
311
312        ! These are used when loading in dimensions from NC file
313        CHARACTER(LEN=NF90_MAX_NAME) :: nxmax_dimname, nymax_dimname, nzmax_dimname, &
314&                                       nuvzmax_dimname, nwzmax_dimname
315
316        ! These are temporary variables, used in the LOAD option, for
317        ! comparing against the current values in FLEXPART of nxmax, nymax, ...
318        INTEGER :: temp_nxmax, temp_nymax, temp_nzmax, &
319&                  temp_nuvzmax, temp_nwzmax
320
321        CHARACTER(LEN=12) :: temp_preproc_format_version_str
322
323        CHARACTER(LEN=128) :: errmesg
324
325        INTEGER, PARAMETER :: DEF_LEVEL = 9
326
327        if (op == 'DUMP') THEN
328
329
330            ! Write the preprocessing format version string
331            WRITE (iounit) PREPROC_FORMAT_VERSION_STR
332
333            ! Write the compiled max dimensions from par_mod - these are
334            ! not meant to be reassigned during a LOAD, but used as "header"
335            ! information to provide the structure of arrays
336            WRITE (iounit) nxmax, nymax, nzmax, nuvzmax, nwzmax
337
338            ncret = nf90_def_dim(ncid, 'nxmax', nxmax, nxmax_dimid)
339            ncret = nf90_def_dim(ncid, 'nymax', nymax, nymax_dimid)
340            ncret = nf90_def_dim(ncid, 'nzmax', nzmax, nzmax_dimid)
341            ncret = nf90_def_dim(ncid, 'nuvzmax', nuvzmax, nuvzmax_dimid)
342            ncret = nf90_def_dim(ncid, 'nwzmax', nwzmax, nwzmax_dimid)
343
344
345
346            ! Scalar values
347            WRITE(iounit) nx, ny, nxmin1, nymin1, nxfield
348            WRITE(iounit) nuvz, nwz, nz, nmixz, nlev_ec
349            WRITE(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst
350
351            ncret = nf90_def_var(ncid, 'nx', NF90_INT, ncvarid)
352            ncret = nf90_put_var(ncid, ncvarid, nx)
353
354            ncret = nf90_def_var(ncid, 'ny', NF90_INT, ncvarid)
355            ncret = nf90_put_var(ncid, ncvarid, ny)
356
357            ncret = nf90_def_var(ncid, 'nxmin1', NF90_INT, ncvarid)
358            ncret = nf90_put_var(ncid, ncvarid, nxmin1)
359
360            ncret = nf90_def_var(ncid, 'nymin1', NF90_INT, ncvarid)
361            ncret = nf90_put_var(ncid, ncvarid, nymin1)
362
363            ncret = nf90_def_var(ncid, 'nxfield', NF90_INT, ncvarid)
364            ncret = nf90_put_var(ncid, ncvarid, nxfield)
365
366            ncret = nf90_def_var(ncid, 'nuvz', NF90_INT, ncvarid)
367            ncret = nf90_put_var(ncid, ncvarid, nuvz)
368
369            ncret = nf90_def_var(ncid, 'nwz', NF90_INT, ncvarid)
370            ncret = nf90_put_var(ncid, ncvarid, nwz)
371
372            ncret = nf90_def_var(ncid, 'nz', NF90_INT, ncvarid)
373            ncret = nf90_put_var(ncid, ncvarid, nz)
374
375            ncret = nf90_def_var(ncid, 'nmixz', NF90_INT, ncvarid)
376            ncret = nf90_put_var(ncid, ncvarid, nmixz)
377
378            ncret = nf90_def_var(ncid, 'nlev_', NF90_INT, ncvarid)
379            ncret = nf90_put_var(ncid, ncvarid, nlev_ec)
380
381            ncret = nf90_def_var(ncid, 'dx', NF90_FLOAT, ncvarid)
382            ncret = nf90_put_var(ncid, ncvarid, dx)
383
384            ncret = nf90_def_var(ncid, 'dy', NF90_FLOAT, ncvarid)
385            ncret = nf90_put_var(ncid, ncvarid, dy)
386
387            ncret = nf90_def_var(ncid, 'xlon0', NF90_FLOAT, ncvarid)
388            ncret = nf90_put_var(ncid, ncvarid, xlon0)
389
390            ncret = nf90_def_var(ncid, 'ylat0', NF90_FLOAT, ncvarid)
391            ncret = nf90_put_var(ncid, ncvarid, ylat0)
392
393            ncret = nf90_def_var(ncid, 'dxconst', NF90_FLOAT, ncvarid)
394            ncret = nf90_put_var(ncid, ncvarid, dxconst)
395
396            ncret = nf90_def_var(ncid, 'dyconst', NF90_FLOAT, ncvarid)
397            ncret = nf90_put_var(ncid, ncvarid, dyconst)
398
399
400
401            ! Fixed fields, static in time
402            WRITE(iounit) oro, excessoro, lsm, xlanduse, height
403
404            dim2dids = (/nxmax_dimid, nymax_dimid/)
405
406            ncret = nf90_def_var(ncid, 'oro', NF90_FLOAT, &
407&                                       dim2dids, ncvarid)
408            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
409&                                        shuffle=0,     &
410&                                        deflate=1,     &
411&                                        deflate_level=DEF_LEVEL)
412            ncret = nf90_put_var(ncid, ncvarid, &
413&                                oro(0:nxmax-1, 0:nymax-1))
414
415            ncret = nf90_def_var(ncid, 'excessoro', NF90_FLOAT, &
416&                                       dim2dids, ncvarid)
417            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
418&                                        shuffle=0,     &
419&                                        deflate=1,     &
420&                                        deflate_level=DEF_LEVEL)
421            ncret = nf90_put_var(ncid, ncvarid, &
422&                                excessoro(0:nxmax-1, 0:nymax-1))
423
424            ncret = nf90_def_var(ncid, 'lsm', NF90_FLOAT, &
425&                                       dim2dids, ncvarid)
426            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
427&                                        shuffle=0,     &
428&                                        deflate=1,     &
429&                                        deflate_level=DEF_LEVEL)
430            ncret = nf90_put_var(ncid, ncvarid, &
431&                                lsm(0:nxmax-1, 0:nymax-1))
432
433            dim3dids = (/nxmax_dimid, nymax_dimid, numclass/)
434            ! numclass comes from par_mod - number of land use classes
435            ncret = nf90_def_var(ncid, 'xlanduse', NF90_FLOAT, &
436&                                       dim3dids, ncvarid)
437            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
438&                                        shuffle=0,     &
439&                                        deflate=1,     &
440&                                        deflate_level=DEF_LEVEL)
441            ncret = nf90_put_var(ncid, ncvarid, &
442&                                xlanduse(0:nxmax-1, 0:nymax-1, 1:numclass))
443
444            dim1dids = (/nzmax_dimid/)
445            ncret = nf90_def_var(ncid, 'height', NF90_FLOAT, &
446&                                       dim1dids, ncvarid)
447            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
448&                                        shuffle=0,     &
449&                                        deflate=1,     &
450&                                        deflate_level=DEF_LEVEL)
451            ncret = nf90_put_var(ncid, ncvarid, &
452&                                height(1:nzmax))
453
454
455
456
457            ! 3d fields
458            WRITE(iounit) uu(:,:,:,cm_index)
459            WRITE(iounit) vv(:,:,:,cm_index)
460            WRITE(iounit) uupol(:,:,:,cm_index)
461            WRITE(iounit) vvpol(:,:,:,cm_index)
462            WRITE(iounit) ww(:,:,:,cm_index)
463            WRITE(iounit) tt(:,:,:,cm_index)
464            WRITE(iounit) qv(:,:,:,cm_index)
465            WRITE(iounit) pv(:,:,:,cm_index)
466            WRITE(iounit) rho(:,:,:,cm_index)
467            WRITE(iounit) drhodz(:,:,:,cm_index)
468            WRITE(iounit) tth(:,:,:,cm_index)
469            WRITE(iounit) qvh(:,:,:,cm_index)
470            WRITE(iounit) pplev(:,:,:,cm_index)
471            WRITE(iounit) clouds(:,:,:,cm_index)
472            WRITE(iounit) cloudsh(:,:,cm_index)
473
474            dim3dids = (/nxmax_dimid, nymax_dimid, nzmax_dimid/)
475
476            ncret = nf90_def_var(ncid, 'uu', NF90_FLOAT, &
477&                                       dim3dids, ncvarid)
478            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
479&                                        shuffle=0,     &
480&                                        deflate=1,     &
481&                                        deflate_level=DEF_LEVEL)
482            ncret = nf90_put_var(ncid, ncvarid, &
483&                                uu(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
484
485            ncret = nf90_def_var(ncid, 'vv', NF90_FLOAT, &
486&                                       dim3dids, ncvarid)
487            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
488&                                        shuffle=0,     &
489&                                        deflate=1,     &
490&                                        deflate_level=DEF_LEVEL)
491            ncret = nf90_put_var(ncid, ncvarid, &
492&                                vv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
493
494            ncret = nf90_def_var(ncid, 'uupol', NF90_FLOAT, &
495&                                       dim3dids, ncvarid)
496            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
497&                                        shuffle=0,     &
498&                                        deflate=1,     &
499&                                        deflate_level=DEF_LEVEL)
500            ncret = nf90_put_var(ncid, ncvarid, &
501&                                uupol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
502
503            ncret = nf90_def_var(ncid, 'vvpol', NF90_FLOAT, &
504&                                       dim3dids, ncvarid)
505            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
506&                                        shuffle=0,     &
507&                                        deflate=1,     &
508&                                        deflate_level=DEF_LEVEL)
509            ncret = nf90_put_var(ncid, ncvarid, &
510&                                vvpol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
511
512            ncret = nf90_def_var(ncid, 'ww', NF90_FLOAT, &
513&                                       dim3dids, ncvarid)
514            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
515&                                        shuffle=0,     &
516&                                        deflate=1,     &
517&                                        deflate_level=DEF_LEVEL)
518            ncret = nf90_put_var(ncid, ncvarid, &
519&                                ww(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
520
521            ncret = nf90_def_var(ncid, 'tt', NF90_FLOAT, &
522&                                       dim3dids, ncvarid)
523            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
524&                                        shuffle=0,     &
525&                                        deflate=1,     &
526&                                        deflate_level=DEF_LEVEL)
527            ncret = nf90_put_var(ncid, ncvarid, &
528&                                tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
529
530            ncret = nf90_def_var(ncid, 'qv', NF90_FLOAT, &
531&                                       dim3dids, ncvarid)
532            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
533&                                        shuffle=0,     &
534&                                        deflate=1,     &
535&                                        deflate_level=DEF_LEVEL)
536            ncret = nf90_put_var(ncid, ncvarid, &
537&                                qv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
538
539            ncret = nf90_def_var(ncid, 'pv', NF90_FLOAT, &
540&                                       dim3dids, ncvarid)
541            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
542&                                        shuffle=0,     &
543&                                        deflate=1,     &
544&                                        deflate_level=DEF_LEVEL)
545            ncret = nf90_put_var(ncid, ncvarid, &
546&                                pv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
547
548            ncret = nf90_def_var(ncid, 'rho', NF90_FLOAT, &
549&                                       dim3dids, ncvarid)
550            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
551&                                        shuffle=0,     &
552&                                        deflate=1,     &
553&                                        deflate_level=DEF_LEVEL)
554            ncret = nf90_put_var(ncid, ncvarid, &
555&                                rho(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
556
557            ncret = nf90_def_var(ncid, 'drhodz', NF90_FLOAT, &
558&                                       dim3dids, ncvarid)
559            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
560&                                        shuffle=0,     &
561&                                        deflate=1,     &
562&                                        deflate_level=DEF_LEVEL)
563            ncret = nf90_put_var(ncid, ncvarid, &
564&                                drhodz(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
565
566            ncret = nf90_def_var(ncid, 'clouds', NF90_BYTE, &
567&                                       dim3dids, ncvarid)
568            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
569&                                        shuffle=0,     &
570&                                        deflate=1,     &
571&                                        deflate_level=DEF_LEVEL)
572            ncret = nf90_put_var(ncid, ncvarid, &
573&                                clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
574
575
576
577            ! Note the change in z dimension for the following
578            dim3dids = (/nxmax_dimid, nymax_dimid, nuvzmax_dimid/)
579
580            ncret = nf90_def_var(ncid, 'tth', NF90_FLOAT, &
581&                                       dim3dids, ncvarid)
582            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
583&                                        shuffle=0,     &
584&                                        deflate=1,     &
585&                                        deflate_level=DEF_LEVEL)
586            ncret = nf90_put_var(ncid, ncvarid, &
587&                                tth(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
588
589            ncret = nf90_def_var(ncid, 'qvh', NF90_FLOAT, &
590&                                       dim3dids, ncvarid)
591            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
592&                                        shuffle=0,     &
593&                                        deflate=1,     &
594&                                        deflate_level=DEF_LEVEL)
595            ncret = nf90_put_var(ncid, ncvarid, &
596&                                qvh(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
597
598            ncret = nf90_def_var(ncid, 'pplev', NF90_FLOAT, &
599&                                       dim3dids, ncvarid)
600            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
601&                                        shuffle=0,     &
602&                                        deflate=1,     &
603&                                        deflate_level=DEF_LEVEL)
604            ncret = nf90_put_var(ncid, ncvarid, &
605&                                pplev(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
606
607
608            dim2dids = (/nxmax_dimid, nymax_dimid/)
609            ncret = nf90_def_var(ncid, 'cloudsh', NF90_INT, &
610&                                       dim2dids, ncvarid)
611            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
612&                                        shuffle=0,     &
613&                                        deflate=1,     &
614&                                        deflate_level=DEF_LEVEL)
615            ncret = nf90_put_var(ncid, ncvarid, &
616&                                cloudsh(0:nxmax-1, 0:nymax-1, cm_index))
617
618
619
620            PRINT *, 'SUM(tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
621&                                        SUM(tt(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
622
623            PRINT *, 'SUM(clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
624&                                        SUM(clouds(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
625
626
627
628            ! 2d fields
629            WRITE(iounit) ps(:,:,:,cm_index)
630            WRITE(iounit) sd(:,:,:,cm_index)
631            WRITE(iounit) msl(:,:,:,cm_index)
632            WRITE(iounit) tcc(:,:,:,cm_index)
633            WRITE(iounit) u10(:,:,:,cm_index)
634            WRITE(iounit) v10(:,:,:,cm_index)
635            WRITE(iounit) tt2(:,:,:,cm_index)
636            WRITE(iounit) td2(:,:,:,cm_index)
637            WRITE(iounit) lsprec(:,:,:,cm_index)
638            WRITE(iounit) convprec(:,:,:,cm_index)
639            WRITE(iounit) sshf(:,:,:,cm_index)
640            WRITE(iounit) ssr(:,:,:,cm_index)
641            WRITE(iounit) surfstr(:,:,:,cm_index)
642            WRITE(iounit) ustar(:,:,:,cm_index)
643            WRITE(iounit) wstar(:,:,:,cm_index)
644            WRITE(iounit) hmix(:,:,:,cm_index)
645            WRITE(iounit) tropopause(:,:,:,cm_index)
646            WRITE(iounit) oli(:,:,:,cm_index)
647            WRITE(iounit) diffk(:,:,:,cm_index)
648            WRITE(iounit) vdep(:,:,:,cm_index)
649
650            ! 1d fields
651            WRITE(iounit) z0(:)
652            WRITE(iounit) akm(:)
653            WRITE(iounit) bkm(:)
654            WRITE(iounit) akz(:)
655            WRITE(iounit) bkz(:)
656            WRITE(iounit) aknew(:)
657            WRITE(iounit) bknew(:)
658
659            ! Nested, scalar values (for each nest)
660            WRITE(iounit) nxn(:)
661            WRITE(iounit) nyn(:)
662            WRITE(iounit) dxn(:)
663            WRITE(iounit) dyn(:)
664            WRITE(iounit) xlon0n(:)
665            WRITE(iounit) ylat0n(:)
666
667            ! Nested fields, static over time
668            WRITE(iounit) oron, excessoron, lsmn, xlandusen
669
670            ! 3d nested fields
671            WRITE(iounit) uun(:,:,:,cm_index,:)
672            WRITE(iounit) vvn(:,:,:,cm_index,:)
673            WRITE(iounit) wwn(:,:,:,cm_index,:)
674            WRITE(iounit) ttn(:,:,:,cm_index,:)
675            WRITE(iounit) qvn(:,:,:,cm_index,:)
676            WRITE(iounit) pvn(:,:,:,cm_index,:)
677            WRITE(iounit) cloudsn(:,:,:,cm_index,:)
678            WRITE(iounit) cloudsnh(:,:,cm_index,:)
679            WRITE(iounit) rhon(:,:,:,cm_index,:)
680            WRITE(iounit) drhodzn(:,:,:,cm_index,:)
681            WRITE(iounit) tthn(:,:,:,cm_index,:)
682            WRITE(iounit) qvhn(:,:,:,cm_index,:)
683
684            ! 2d nested fields
685            WRITE(iounit) psn(:,:,:,cm_index,:)
686            WRITE(iounit) sdn(:,:,:,cm_index,:)
687            WRITE(iounit) msln(:,:,:,cm_index,:)
688            WRITE(iounit) tccn(:,:,:,cm_index,:)
689            WRITE(iounit) u10n(:,:,:,cm_index,:)
690            WRITE(iounit) v10n(:,:,:,cm_index,:)
691            WRITE(iounit) tt2n(:,:,:,cm_index,:)
692            WRITE(iounit) td2n(:,:,:,cm_index,:)
693            WRITE(iounit) lsprecn(:,:,:,cm_index,:)
694            WRITE(iounit) convprecn(:,:,:,cm_index,:)
695            WRITE(iounit) sshfn(:,:,:,cm_index,:)
696            WRITE(iounit) ssrn(:,:,:,cm_index,:)
697            WRITE(iounit) surfstrn(:,:,:,cm_index,:)
698            WRITE(iounit) ustarn(:,:,:,cm_index,:)
699            WRITE(iounit) wstarn(:,:,:,cm_index,:)
700            WRITE(iounit) hmixn(:,:,:,cm_index,:)
701            WRITE(iounit) tropopausen(:,:,:,cm_index,:)
702            WRITE(iounit) olin(:,:,:,cm_index,:)
703            WRITE(iounit) diffkn(:,:,:,cm_index,:)
704            WRITE(iounit) vdepn(:,:,:,cm_index,:)
705
706            ! Auxiliary variables for nests
707            WRITE(iounit) xresoln(:)
708            WRITE(iounit) yresoln(:)
709            WRITE(iounit) xln(:)
710            WRITE(iounit) yln(:)
711            WRITE(iounit) xrn(:)
712            WRITE(iounit) yrn(:)
713
714            ! Variables for polar stereographic projection
715            WRITE(iounit) xglobal, sglobal, nglobal
716            WRITE(iounit) switchnorthg, switchsouthg
717            WRITE(iounit) southpolemap(:)
718            WRITE(iounit) northpolemap(:)
719
720            ! Variables declared in conv_mod (convection)
721            WRITE(iounit) pconv(:)
722            WRITE(iounit) phconv(:)
723            WRITE(iounit) dpr(:)
724            WRITE(iounit) pconv_hpa(:)
725            WRITE(iounit) phconv_hpa(:)
726            WRITE(iounit) ft(:)
727            WRITE(iounit) fq(:)
728            WRITE(iounit) fmass(:,:)
729            WRITE(iounit) sub(:)
730            WRITE(iounit) fmassfrac(:,:)
731            WRITE(iounit) cbaseflux(:,:)
732            WRITE(iounit) cbasefluxn(:,:,:)
733            WRITE(iounit) tconv(:)
734            WRITE(iounit) qconv(:)
735            WRITE(iounit) qsconv(:)
736            WRITE(iounit) psconv, tt2conv, td2conv
737            WRITE(iounit) nconvlev, nconvtop
738
739        ELSE IF (op == 'LOAD') THEN
740
741            ! Read the preprocessed format version string and insure it
742            ! matches this version
743            READ (iounit) temp_preproc_format_version_str
744            PRINT *, 'Reading preprocessed file format version: ', &
745&                    temp_preproc_format_version_str
746
747            IF (TRIM(temp_preproc_format_version_str) == &
748&                                        TRIM(PREPROC_FORMAT_VERSION_STR)) THEN
749                CONTINUE
750            ELSE
751                ! PRINT *, ''  GK: causes relocation truncated to fit: R_X86_64_32
752                PRINT *, 'Inconsistent preprocessing format version'
753                PRINT *, 'Expected Version: ', PREPROC_FORMAT_VERSION_STR
754                PRINT *, 'Detected Version: ', temp_preproc_format_version_str
755                ! PRINT *, ''
756                STOP
757            END IF
758
759            ! Read the compiled max dimensions that were dumped from par_mod
760            ! when creating the fp file, so that we can compare against
761            ! current FLEXPART dimensions - they need to be the same, or else
762            ! we abort.
763            READ (iounit) temp_nxmax, temp_nymax, temp_nzmax, &
764&                         temp_nuvzmax, temp_nwzmax
765
766            ! Get dimensions
767            ncret = nf90_inq_dimid(ncid, 'nxmax', nxmax_dimid)
768            ncret = nf90_inquire_dimension(ncid, nxmax_dimid, nxmax_dimname, &
769&                                                temp_nxmax)
770            PRINT *, 'temp_nxmax: ', temp_nxmax
771
772            ncret = nf90_inq_dimid(ncid, 'nymax', nymax_dimid)
773            ncret = nf90_inquire_dimension(ncid, nymax_dimid, nymax_dimname, &
774&                                                temp_nymax)
775            PRINT *, 'temp_nymax: ', temp_nymax
776
777            ncret = nf90_inq_dimid(ncid, 'nzmax', nzmax_dimid)
778            ncret = nf90_inquire_dimension(ncid, nzmax_dimid, nzmax_dimname, &
779&                                                temp_nzmax)
780            PRINT *, 'temp_nzmax: ', temp_nzmax
781
782            ncret = nf90_inq_dimid(ncid, 'nuvzmax', nuvzmax_dimid)
783            ncret = nf90_inquire_dimension(ncid, nuvzmax_dimid, nuvzmax_dimname, &
784&                                                temp_nuvzmax)
785            PRINT *, 'temp_nuvzmax: ', temp_nuvzmax
786
787            ncret = nf90_inq_dimid(ncid, 'nwzmax', nwzmax_dimid)
788            ncret = nf90_inquire_dimension(ncid, nwzmax_dimid, nwzmax_dimname, &
789&                                                temp_nwzmax)
790            PRINT *, 'temp_nwzmax: ', temp_nwzmax
791
792
793
794            IF ( (temp_nxmax == nxmax) .AND. (temp_nymax == nymax) .AND. &
795&                   (temp_nzmax == nzmax) .AND. &
796&                   (temp_nuvzmax == nuvzmax) .AND. &
797&                   (temp_nwzmax == nwzmax) ) THEN
798                CONTINUE
799            ELSE
800                PRINT *, 'Incompatible dimensions between fp file and current FLEXPART!'
801                ! PRINT *, ''
802                PRINT *, '                  FP file     Compiled FP'
803                PRINT *, 'nxmax:     ', temp_nxmax, '    ', nxmax
804                PRINT *, 'nymax:     ', temp_nymax, '    ', nymax
805                PRINT *, 'nzmax:     ', temp_nzmax, '    ', nzmax
806                PRINT *, 'nuvzmax:     ', temp_nuvzmax, '    ', nuvzmax
807                PRINT *, 'nwzmax:     ', temp_nwzmax, '    ', nwzmax
808                ! PRINT *, ''
809                STOP
810            END IF
811
812
813            ! Scalar values
814            READ(iounit) nx, ny, nxmin1, nymin1, nxfield
815            READ(iounit) nuvz, nwz, nz, nmixz, nlev_ec
816            READ(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst
817
818
819
820            ! Get the varid , then read into scalar variable
821            ncret = nf90_inq_varid(ncid, 'nx', ncvarid)
822            ncret = nf90_get_var(ncid, ncvarid, nx)
823
824            ncret = nf90_inq_varid(ncid, 'ny', ncvarid)
825            ncret = nf90_get_var(ncid, ncvarid, ny)
826
827            ncret = nf90_inq_varid(ncid, 'nxmin1', ncvarid)
828            ncret = nf90_get_var(ncid, ncvarid, nxmin1)
829
830            ncret = nf90_inq_varid(ncid, 'nymin1', ncvarid)
831            ncret = nf90_get_var(ncid, ncvarid, nymin1)
832
833            ncret = nf90_inq_varid(ncid, 'nxfield', ncvarid)
834            ncret = nf90_get_var(ncid, ncvarid, nxfield)
835
836            ncret = nf90_inq_varid(ncid, 'nuvz', ncvarid)
837            ncret = nf90_get_var(ncid, ncvarid, nuvz)
838
839            ncret = nf90_inq_varid(ncid, 'nwz', ncvarid)
840            ncret = nf90_get_var(ncid, ncvarid, nwz)
841
842            ncret = nf90_inq_varid(ncid, 'nz', ncvarid)
843            ncret = nf90_get_var(ncid, ncvarid, nz)
844
845            ncret = nf90_inq_varid(ncid, 'nmixz', ncvarid)
846            ncret = nf90_get_var(ncid, ncvarid, nmixz)
847
848            ncret = nf90_inq_varid(ncid, 'nlev_ec', ncvarid)
849            ncret = nf90_get_var(ncid, ncvarid, nlev_ec)
850
851            ncret = nf90_inq_varid(ncid, 'dx', ncvarid)
852            ncret = nf90_get_var(ncid, ncvarid, dx)
853
854            ncret = nf90_inq_varid(ncid, 'dy', ncvarid)
855            ncret = nf90_get_var(ncid, ncvarid, dy)
856
857            ncret = nf90_inq_varid(ncid, 'xlon0', ncvarid)
858            ncret = nf90_get_var(ncid, ncvarid, xlon0)
859
860            ncret = nf90_inq_varid(ncid, 'ylat0', ncvarid)
861            ncret = nf90_get_var(ncid, ncvarid, ylat0)
862
863            ncret = nf90_inq_varid(ncid, 'dxconst', ncvarid)
864            ncret = nf90_get_var(ncid, ncvarid, dxconst)
865
866            ncret = nf90_inq_varid(ncid, 'dyconst', ncvarid)
867            ncret = nf90_get_var(ncid, ncvarid, dyconst)
868
869
870
871
872
873
874            ! Fixed fields, static in time
875            READ(iounit) oro, excessoro, lsm, xlanduse, height
876
877            ncret = nf90_inq_varid(ncid, 'oro', ncvarid)
878            ncret = nf90_get_var(ncid, ncvarid, oro(0:nxmax-1,0:nymax-1))
879
880            ncret = nf90_inq_varid(ncid, 'excessoro', ncvarid)
881            ncret = nf90_get_var(ncid, ncvarid, excessoro(0:nxmax-1,0:nymax-1))
882
883            ncret = nf90_inq_varid(ncid, 'lsm', ncvarid)
884            ncret = nf90_get_var(ncid, ncvarid, lsm(0:nxmax-1,0:nymax-1))
885
886            ncret = nf90_inq_varid(ncid, 'xlanduse', ncvarid)
887            ncret = nf90_get_var(ncid, ncvarid, xlanduse(0:nxmax-1,0:nymax-1, 1:numclass))
888
889            ncret = nf90_inq_varid(ncid, 'height', ncvarid)
890            ncret = nf90_get_var(ncid, ncvarid, height(1:nzmax))
891
892
893
894
895            ! 3d fields
896            READ(iounit) uu(:,:,:,cm_index)
897            READ(iounit) vv(:,:,:,cm_index)
898            READ(iounit) uupol(:,:,:,cm_index)
899            READ(iounit) vvpol(:,:,:,cm_index)
900            READ(iounit) ww(:,:,:,cm_index)
901            READ(iounit) tt(:,:,:,cm_index)
902            READ(iounit) qv(:,:,:,cm_index)
903            READ(iounit) pv(:,:,:,cm_index)
904            READ(iounit) rho(:,:,:,cm_index)
905            READ(iounit) drhodz(:,:,:,cm_index)
906            READ(iounit) tth(:,:,:,cm_index)
907            READ(iounit) qvh(:,:,:,cm_index)
908            READ(iounit) pplev(:,:,:,cm_index)
909            READ(iounit) clouds(:,:,:,cm_index)
910            READ(iounit) cloudsh(:,:,cm_index)
911
912
913
914
915            ! Get the varid and read the variable into the array
916            ncret = nf90_inq_varid(ncid, 'uu', ncvarid)
917            ncret = nf90_get_var(ncid, ncvarid, uu(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
918
919            ncret = nf90_inq_varid(ncid, 'vv', ncvarid)
920            ncret = nf90_get_var(ncid, ncvarid, vv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
921
922            ncret = nf90_inq_varid(ncid, 'uupol', ncvarid)
923            ncret = nf90_get_var(ncid, ncvarid, uupol(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
924
925            ncret = nf90_inq_varid(ncid, 'vvpol', ncvarid)
926            ncret = nf90_get_var(ncid, ncvarid, vvpol(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
927
928            ncret = nf90_inq_varid(ncid, 'ww', ncvarid)
929            ncret = nf90_get_var(ncid, ncvarid, ww(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
930
931            ncret = nf90_inq_varid(ncid, 'tt', ncvarid)
932            ncret = nf90_get_var(ncid, ncvarid, tt(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
933
934            ncret = nf90_inq_varid(ncid, 'qv', ncvarid)
935            ncret = nf90_get_var(ncid, ncvarid, qv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
936
937            ncret = nf90_inq_varid(ncid, 'pv', ncvarid)
938            ncret = nf90_get_var(ncid, ncvarid, pv(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
939
940            ncret = nf90_inq_varid(ncid, 'rho', ncvarid)
941            ncret = nf90_get_var(ncid, ncvarid, rho(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
942
943            ncret = nf90_inq_varid(ncid, 'drhodz', ncvarid)
944            ncret = nf90_get_var(ncid, ncvarid, drhodz(0:nxmax-1,0:nymax-1,1:nzmax,cm_index))
945
946            ncret = nf90_inq_varid(ncid, 'clouds', ncvarid)
947            ncret = nf90_get_var(ncid, ncvarid, clouds(0:nxmax-1,0:nzmax-1,1:nzmax,cm_index))
948
949            ncret = nf90_inq_varid(ncid, 'tth', ncvarid)
950            ncret = nf90_get_var(ncid, ncvarid, tth(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
951
952            ncret = nf90_inq_varid(ncid, 'qvh', ncvarid)
953            ncret = nf90_get_var(ncid, ncvarid, qvh(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
954
955            ncret = nf90_inq_varid(ncid, 'pplev', ncvarid)
956            ncret = nf90_get_var(ncid, ncvarid, pplev(0:nxmax-1,0:nymax-1,1:nuvzmax,cm_index))
957
958            ncret = nf90_inq_varid(ncid, 'cloudsh', ncvarid)
959            ncret = nf90_get_var(ncid, ncvarid, cloudsh(0:nxmax-1,0:nymax-1,cm_index))
960
961
962
963
964            PRINT *, 'SUM(tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
965&                                        SUM(tt(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
966
967
968            PRINT *, 'SUM(clouds(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index)): ', &
969&                                        SUM(clouds(0:nxmax-1,0:nymax-1,1:nzmax, cm_index))
970
971
972
973
974            ! 2d fields
975            READ(iounit) ps(:,:,:,cm_index)
976            READ(iounit) sd(:,:,:,cm_index)
977            READ(iounit) msl(:,:,:,cm_index)
978            READ(iounit) tcc(:,:,:,cm_index)
979            READ(iounit) u10(:,:,:,cm_index)
980            READ(iounit) v10(:,:,:,cm_index)
981            READ(iounit) tt2(:,:,:,cm_index)
982            READ(iounit) td2(:,:,:,cm_index)
983            READ(iounit) lsprec(:,:,:,cm_index)
984            READ(iounit) convprec(:,:,:,cm_index)
985            READ(iounit) sshf(:,:,:,cm_index)
986            READ(iounit) ssr(:,:,:,cm_index)
987            READ(iounit) surfstr(:,:,:,cm_index)
988            READ(iounit) ustar(:,:,:,cm_index)
989            READ(iounit) wstar(:,:,:,cm_index)
990            READ(iounit) hmix(:,:,:,cm_index)
991            READ(iounit) tropopause(:,:,:,cm_index)
992            READ(iounit) oli(:,:,:,cm_index)
993            READ(iounit) diffk(:,:,:,cm_index)
994            READ(iounit) vdep(:,:,:,cm_index)
995
996            ! 1d fields
997            READ(iounit) z0(:)
998            READ(iounit) akm(:)
999            READ(iounit) bkm(:)
1000            READ(iounit) akz(:)
1001            READ(iounit) bkz(:)
1002            READ(iounit) aknew(:)
1003            READ(iounit) bknew(:)
1004
1005
1006            ! Nested, scalar values (for each nest)
1007            READ(iounit) nxn(:)
1008            READ(iounit) nyn(:)
1009            READ(iounit) dxn(:)
1010            READ(iounit) dyn(:)
1011            READ(iounit) xlon0n(:)
1012            READ(iounit) ylat0n(:)
1013
1014
1015            ! Nested fields, static over time
1016            READ(iounit) oron, excessoron, lsmn, xlandusen
1017
1018            ! 3d nested fields
1019            READ(iounit) uun(:,:,:,cm_index,:)
1020            READ(iounit) vvn(:,:,:,cm_index,:)
1021            READ(iounit) wwn(:,:,:,cm_index,:)
1022            READ(iounit) ttn(:,:,:,cm_index,:)
1023            READ(iounit) qvn(:,:,:,cm_index,:)
1024            READ(iounit) pvn(:,:,:,cm_index,:)
1025            READ(iounit) cloudsn(:,:,:,cm_index,:)
1026            READ(iounit) cloudsnh(:,:,cm_index,:)
1027            READ(iounit) rhon(:,:,:,cm_index,:)
1028            READ(iounit) drhodzn(:,:,:,cm_index,:)
1029            READ(iounit) tthn(:,:,:,cm_index,:)
1030            READ(iounit) qvhn(:,:,:,cm_index,:)
1031
1032            ! 2d nested fields
1033            READ(iounit) psn(:,:,:,cm_index,:)
1034            READ(iounit) sdn(:,:,:,cm_index,:)
1035            READ(iounit) msln(:,:,:,cm_index,:)
1036            READ(iounit) tccn(:,:,:,cm_index,:)
1037            READ(iounit) u10n(:,:,:,cm_index,:)
1038            READ(iounit) v10n(:,:,:,cm_index,:)
1039            READ(iounit) tt2n(:,:,:,cm_index,:)
1040            READ(iounit) td2n(:,:,:,cm_index,:)
1041            READ(iounit) lsprecn(:,:,:,cm_index,:)
1042            READ(iounit) convprecn(:,:,:,cm_index,:)
1043            READ(iounit) sshfn(:,:,:,cm_index,:)
1044            READ(iounit) ssrn(:,:,:,cm_index,:)
1045            READ(iounit) surfstrn(:,:,:,cm_index,:)
1046            READ(iounit) ustarn(:,:,:,cm_index,:)
1047            READ(iounit) wstarn(:,:,:,cm_index,:)
1048            READ(iounit) hmixn(:,:,:,cm_index,:)
1049            READ(iounit) tropopausen(:,:,:,cm_index,:)
1050            READ(iounit) olin(:,:,:,cm_index,:)
1051            READ(iounit) diffkn(:,:,:,cm_index,:)
1052            READ(iounit) vdepn(:,:,:,cm_index,:)
1053
1054            ! Auxiliary variables for nests
1055            READ(iounit) xresoln(:)
1056            READ(iounit) yresoln(:)
1057            READ(iounit) xln(:)
1058            READ(iounit) yln(:)
1059            READ(iounit) xrn(:)
1060            READ(iounit) yrn(:)
1061
1062            ! Variables for polar stereographic projection
1063            READ(iounit) xglobal, sglobal, nglobal
1064            READ(iounit) switchnorthg, switchsouthg
1065            READ(iounit) southpolemap(:)
1066            READ(iounit) northpolemap(:)
1067
1068            ! Variables declared in conv_mod (convection)
1069            READ(iounit) pconv(:)
1070            READ(iounit) phconv(:)
1071            READ(iounit) dpr(:)
1072            READ(iounit) pconv_hpa(:)
1073            READ(iounit) phconv_hpa(:)
1074            READ(iounit) ft(:)
1075            READ(iounit) fq(:)
1076            READ(iounit) fmass(:,:)
1077            READ(iounit) sub(:)
1078            READ(iounit) fmassfrac(:,:)
1079            READ(iounit) cbaseflux(:,:)
1080            READ(iounit) cbasefluxn(:,:,:)
1081            READ(iounit) tconv(:)
1082            READ(iounit) qconv(:)
1083            READ(iounit) qsconv(:)
1084            READ(iounit) psconv, tt2conv, td2conv
1085            READ(iounit) nconvlev, nconvtop
1086
1087        ELSE
1088            STOP 'fpio(): Illegal operation'
1089           
1090        ENDIF
1091    END SUBROUTINE fpio
1092
1093    SUBROUTINE fpmetbinary_filetext(filename, cm_index)
1094
1095        ! This is a utility subroutine meant to be used for testing purposes.
1096        ! It facilitates the text output of variables read in from the
1097        ! specified .fp file.  This routine will easily cause the program
1098        ! to crash due memory allocation issues, particularly when you are
1099        ! trying to text print 3d arrays in a single formatted statetment.
1100       
1101        CHARACTER(LEN=*), INTENT(IN) :: filename
1102        INTEGER, INTENT(IN) :: cm_index           ! Index of last dimension in
1103                                                  ! most com_mod variables.
1104                                                  ! Should be 1 or 2
1105
1106        !OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', status='REPLACE', &
1107        !    form="formatted", access="stream")
1108        OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', &
1109             form="formatted", access="APPEND")
1110
1111        WRITE(IOUNIT_TEXTOUT, *) 'oro: ', oro
1112        WRITE(IOUNIT_TEXTOUT, *) 'excessoro: ', excessoro
1113        WRITE(IOUNIT_TEXTOUT, *) 'lsm: ', lsm
1114        WRITE(IOUNIT_TEXTOUT, *) 'xlanduse: ', xlanduse
1115        WRITE(IOUNIT_TEXTOUT, *) 'height: ', height
1116
1117        WRITE(IOUNIT_TEXTOUT, *) 'uu: ', uu(:,:,:,cm_index)
1118        WRITE(IOUNIT_TEXTOUT, *) 'vv: ', vv(:,:,:,cm_index)
1119        WRITE(IOUNIT_TEXTOUT, *) 'uupol: ', uupol(:,:,:,cm_index)
1120        WRITE(IOUNIT_TEXTOUT, *) 'vvpol: ', vvpol(:,:,:,cm_index)
1121        WRITE(IOUNIT_TEXTOUT, *) 'ww: ', ww(:,:,:,cm_index)
1122        WRITE(IOUNIT_TEXTOUT, *) 'tt: ', tt(:,:,:,cm_index)
1123        WRITE(IOUNIT_TEXTOUT, *) 'qv: ', qv(:,:,:,cm_index)
1124        WRITE(IOUNIT_TEXTOUT, *) 'pv: ', pv(:,:,:,cm_index)
1125        WRITE(IOUNIT_TEXTOUT, *) 'rho: ', rho(:,:,:,cm_index)
1126        WRITE(IOUNIT_TEXTOUT, *) 'drhodz: ', drhodz(:,:,:,cm_index)
1127        WRITE(IOUNIT_TEXTOUT, *) 'tth: ', tth(:,:,:,cm_index)
1128        WRITE(IOUNIT_TEXTOUT, *) 'qvh: ', qvh(:,:,:,cm_index)
1129        WRITE(IOUNIT_TEXTOUT, *) 'pplev: ', pplev(:,:,:,cm_index)
1130        WRITE(IOUNIT_TEXTOUT, *) 'clouds: ', clouds(:,:,:,cm_index)
1131        WRITE(IOUNIT_TEXTOUT, *) 'cloudsh: ', cloudsh(:,:,cm_index)
1132
1133
1134
1135
1136        CLOSE(IOUNIT_TEXTOUT)
1137    END SUBROUTINE fpmetbinary_filetext
1138
1139
1140END MODULE fpmetbinary_mod
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG