source: flexpart.git/flexpart_code/fpmetbinary_mod.F90 @ 929ae37

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

Continued incremental deployment of NC4 preproc files

  • Property mode set to 100644
File size: 42.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
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(3) :: dim3dids    ! Dimension IDs for 3D arrays
309
310        ! These are used when loading in dimensions from NC file
311        CHARACTER(LEN=NF90_MAX_NAME) :: nxmax_dimname, nymax_dimname, nzmax_dimname, &
312&                                       nuvzmax_dimname, nwzmax_dimname
313
314        ! These are temporary variables, used in the LOAD option, for
315        ! comparing against the current values in FLEXPART of nxmax, nymax, ...
316        INTEGER :: temp_nxmax, temp_nymax, temp_nzmax, &
317&                  temp_nuvzmax, temp_nwzmax
318
319        CHARACTER(LEN=12) :: temp_preproc_format_version_str
320
321        CHARACTER(LEN=128) :: errmesg
322
323        INTEGER, PARAMETER :: DEF_LEVEL = 3
324
325        if (op == 'DUMP') THEN
326
327
328            ! Write the preprocessing format version string
329            WRITE (iounit) PREPROC_FORMAT_VERSION_STR
330
331            ! Write the compiled max dimensions from par_mod - these are
332            ! not meant to be reassigned during a LOAD, but used as "header"
333            ! information to provide the structure of arrays
334            WRITE (iounit) nxmax, nymax, nzmax, nuvzmax, nwzmax
335
336            ncret = nf90_def_dim(ncid, 'nxmax', nxmax, nxmax_dimid)
337            ncret = nf90_def_dim(ncid, 'nymax', nymax, nymax_dimid)
338            ncret = nf90_def_dim(ncid, 'nzmax', nzmax, nzmax_dimid)
339            ncret = nf90_def_dim(ncid, 'nuvzmax', nuvzmax, nuvzmax_dimid)
340            ncret = nf90_def_dim(ncid, 'nwzmax', nwzmax, nwzmax_dimid)
341
342
343
344            ! Scalar values
345            WRITE(iounit) nx, ny, nxmin1, nymin1, nxfield
346            WRITE(iounit) nuvz, nwz, nz, nmixz, nlev_ec
347            WRITE(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst
348
349            ncret = nf90_def_var(ncid, 'nx', NF90_INT, ncvarid)
350            ncret = nf90_put_var(ncid, ncvarid, nx)
351
352            ncret = nf90_def_var(ncid, 'ny', NF90_INT, ncvarid)
353            ncret = nf90_put_var(ncid, ncvarid, ny)
354
355            ncret = nf90_def_var(ncid, 'nxmin1', NF90_INT, ncvarid)
356            ncret = nf90_put_var(ncid, ncvarid, nxmin1)
357
358            ncret = nf90_def_var(ncid, 'nymin1', NF90_INT, ncvarid)
359            ncret = nf90_put_var(ncid, ncvarid, nymin1)
360
361            ncret = nf90_def_var(ncid, 'nxfield', NF90_INT, ncvarid)
362            ncret = nf90_put_var(ncid, ncvarid, nxfield)
363
364            ncret = nf90_def_var(ncid, 'nuvz', NF90_INT, ncvarid)
365            ncret = nf90_put_var(ncid, ncvarid, nuvz)
366
367            ncret = nf90_def_var(ncid, 'nwz', NF90_INT, ncvarid)
368            ncret = nf90_put_var(ncid, ncvarid, nwz)
369
370            ncret = nf90_def_var(ncid, 'nz', NF90_INT, ncvarid)
371            ncret = nf90_put_var(ncid, ncvarid, nz)
372
373            ncret = nf90_def_var(ncid, 'nmixz', NF90_INT, ncvarid)
374            ncret = nf90_put_var(ncid, ncvarid, nmixz)
375
376            ncret = nf90_def_var(ncid, 'nlev_', NF90_INT, ncvarid)
377            ncret = nf90_put_var(ncid, ncvarid, nlev_ec)
378
379            ncret = nf90_def_var(ncid, 'dx', NF90_FLOAT, ncvarid)
380            ncret = nf90_put_var(ncid, ncvarid, dx)
381
382            ncret = nf90_def_var(ncid, 'dy', NF90_FLOAT, ncvarid)
383            ncret = nf90_put_var(ncid, ncvarid, dy)
384
385            ncret = nf90_def_var(ncid, 'xlon0', NF90_FLOAT, ncvarid)
386            ncret = nf90_put_var(ncid, ncvarid, xlon0)
387
388            ncret = nf90_def_var(ncid, 'ylat0', NF90_FLOAT, ncvarid)
389            ncret = nf90_put_var(ncid, ncvarid, ylat0)
390
391            ncret = nf90_def_var(ncid, 'dxconst', NF90_FLOAT, ncvarid)
392            ncret = nf90_put_var(ncid, ncvarid, dxconst)
393
394            ncret = nf90_def_var(ncid, 'dyconst', NF90_FLOAT, ncvarid)
395            ncret = nf90_put_var(ncid, ncvarid, dyconst)
396
397
398
399            ! Fixed fields, static in time
400            WRITE(iounit) oro, excessoro, lsm, xlanduse, height
401
402            ! 3d fields
403            WRITE(iounit) uu(:,:,:,cm_index)
404            WRITE(iounit) vv(:,:,:,cm_index)
405            WRITE(iounit) uupol(:,:,:,cm_index)
406            WRITE(iounit) vvpol(:,:,:,cm_index)
407            WRITE(iounit) ww(:,:,:,cm_index)
408            WRITE(iounit) tt(:,:,:,cm_index)
409            WRITE(iounit) qv(:,:,:,cm_index)
410            WRITE(iounit) pv(:,:,:,cm_index)
411            WRITE(iounit) rho(:,:,:,cm_index)
412            WRITE(iounit) drhodz(:,:,:,cm_index)
413            WRITE(iounit) tth(:,:,:,cm_index)
414            WRITE(iounit) qvh(:,:,:,cm_index)
415            WRITE(iounit) pplev(:,:,:,cm_index)
416            WRITE(iounit) clouds(:,:,:,cm_index)
417            WRITE(iounit) cloudsh(:,:,cm_index)
418
419            dim3dids = (/nxmax_dimid, nymax_dimid, nzmax_dimid/)
420
421            ncret = nf90_def_var(ncid, 'uu', NF90_FLOAT, &
422&                                       dim3dids, 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&                                uu(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
429
430            ncret = nf90_def_var(ncid, 'vv', NF90_FLOAT, &
431&                                       dim3dids, 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&                                vv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
438
439            ncret = nf90_def_var(ncid, 'uupol', NF90_FLOAT, &
440&                                       dim3dids, ncvarid)
441            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
442&                                        shuffle=0,     &
443&                                        deflate=1,     &
444&                                        deflate_level=DEF_LEVEL)
445            ncret = nf90_put_var(ncid, ncvarid, &
446&                                uupol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
447
448            ncret = nf90_def_var(ncid, 'vvpol', NF90_FLOAT, &
449&                                       dim3dids, ncvarid)
450            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
451&                                        shuffle=0,     &
452&                                        deflate=1,     &
453&                                        deflate_level=DEF_LEVEL)
454            ncret = nf90_put_var(ncid, ncvarid, &
455&                                vvpol(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
456
457            ncret = nf90_def_var(ncid, 'ww', NF90_FLOAT, &
458&                                       dim3dids, ncvarid)
459            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
460&                                        shuffle=0,     &
461&                                        deflate=1,     &
462&                                        deflate_level=DEF_LEVEL)
463            ncret = nf90_put_var(ncid, ncvarid, &
464&                                ww(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
465
466            ncret = nf90_def_var(ncid, 'tt', NF90_FLOAT, &
467&                                       dim3dids, ncvarid)
468            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
469&                                        shuffle=0,     &
470&                                        deflate=1,     &
471&                                        deflate_level=DEF_LEVEL)
472            ncret = nf90_put_var(ncid, ncvarid, &
473&                                tt(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
474
475            ncret = nf90_def_var(ncid, 'qv', NF90_FLOAT, &
476&                                       dim3dids, ncvarid)
477            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
478&                                        shuffle=0,     &
479&                                        deflate=1,     &
480&                                        deflate_level=DEF_LEVEL)
481            ncret = nf90_put_var(ncid, ncvarid, &
482&                                qv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
483
484            ncret = nf90_def_var(ncid, 'pv', NF90_FLOAT, &
485&                                       dim3dids, ncvarid)
486            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
487&                                        shuffle=0,     &
488&                                        deflate=1,     &
489&                                        deflate_level=DEF_LEVEL)
490            ncret = nf90_put_var(ncid, ncvarid, &
491&                                pv(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
492
493            ncret = nf90_def_var(ncid, 'rho', NF90_FLOAT, &
494&                                       dim3dids, ncvarid)
495            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
496&                                        shuffle=0,     &
497&                                        deflate=1,     &
498&                                        deflate_level=DEF_LEVEL)
499            ncret = nf90_put_var(ncid, ncvarid, &
500&                                rho(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
501
502            ncret = nf90_def_var(ncid, 'drhodz', NF90_FLOAT, &
503&                                       dim3dids, ncvarid)
504            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
505&                                        shuffle=0,     &
506&                                        deflate=1,     &
507&                                        deflate_level=DEF_LEVEL)
508            ncret = nf90_put_var(ncid, ncvarid, &
509&                                drhodz(0:nxmax-1, 0:nymax-1, 1:nzmax, cm_index))
510
511
512            dim3dids = (/nxmax_dimid, nymax_dimid, nuvzmax_dimid/)
513
514            ncret = nf90_def_var(ncid, 'tth', NF90_FLOAT, &
515&                                       dim3dids, ncvarid)
516            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
517&                                        shuffle=0,     &
518&                                        deflate=1,     &
519&                                        deflate_level=DEF_LEVEL)
520            ncret = nf90_put_var(ncid, ncvarid, &
521&                                tth(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
522
523            ncret = nf90_def_var(ncid, 'qvh', NF90_FLOAT, &
524&                                       dim3dids, ncvarid)
525            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
526&                                        shuffle=0,     &
527&                                        deflate=1,     &
528&                                        deflate_level=DEF_LEVEL)
529            ncret = nf90_put_var(ncid, ncvarid, &
530&                                qvh(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
531
532            ncret = nf90_def_var(ncid, 'pplev', NF90_FLOAT, &
533&                                       dim3dids, ncvarid)
534            ncret = nf90_def_var_deflate(ncid, ncvarid,   &
535&                                        shuffle=0,     &
536&                                        deflate=1,     &
537&                                        deflate_level=DEF_LEVEL)
538            ncret = nf90_put_var(ncid, ncvarid, &
539&                                pplev(0:nxmax-1, 0:nymax-1, 1:nuvzmax, cm_index))
540
541
542
543
544
545            ! 2d fields
546            WRITE(iounit) ps(:,:,:,cm_index)
547            WRITE(iounit) sd(:,:,:,cm_index)
548            WRITE(iounit) msl(:,:,:,cm_index)
549            WRITE(iounit) tcc(:,:,:,cm_index)
550            WRITE(iounit) u10(:,:,:,cm_index)
551            WRITE(iounit) v10(:,:,:,cm_index)
552            WRITE(iounit) tt2(:,:,:,cm_index)
553            WRITE(iounit) td2(:,:,:,cm_index)
554            WRITE(iounit) lsprec(:,:,:,cm_index)
555            WRITE(iounit) convprec(:,:,:,cm_index)
556            WRITE(iounit) sshf(:,:,:,cm_index)
557            WRITE(iounit) ssr(:,:,:,cm_index)
558            WRITE(iounit) surfstr(:,:,:,cm_index)
559            WRITE(iounit) ustar(:,:,:,cm_index)
560            WRITE(iounit) wstar(:,:,:,cm_index)
561            WRITE(iounit) hmix(:,:,:,cm_index)
562            WRITE(iounit) tropopause(:,:,:,cm_index)
563            WRITE(iounit) oli(:,:,:,cm_index)
564            WRITE(iounit) diffk(:,:,:,cm_index)
565            WRITE(iounit) vdep(:,:,:,cm_index)
566
567            ! 1d fields
568            WRITE(iounit) z0(:)
569            WRITE(iounit) akm(:)
570            WRITE(iounit) bkm(:)
571            WRITE(iounit) akz(:)
572            WRITE(iounit) bkz(:)
573            WRITE(iounit) aknew(:)
574            WRITE(iounit) bknew(:)
575
576            ! Nested, scalar values (for each nest)
577            WRITE(iounit) nxn(:)
578            WRITE(iounit) nyn(:)
579            WRITE(iounit) dxn(:)
580            WRITE(iounit) dyn(:)
581            WRITE(iounit) xlon0n(:)
582            WRITE(iounit) ylat0n(:)
583
584            ! Nested fields, static over time
585            WRITE(iounit) oron, excessoron, lsmn, xlandusen
586
587            ! 3d nested fields
588            WRITE(iounit) uun(:,:,:,cm_index,:)
589            WRITE(iounit) vvn(:,:,:,cm_index,:)
590            WRITE(iounit) wwn(:,:,:,cm_index,:)
591            WRITE(iounit) ttn(:,:,:,cm_index,:)
592            WRITE(iounit) qvn(:,:,:,cm_index,:)
593            WRITE(iounit) pvn(:,:,:,cm_index,:)
594            WRITE(iounit) cloudsn(:,:,:,cm_index,:)
595            WRITE(iounit) cloudsnh(:,:,cm_index,:)
596            WRITE(iounit) rhon(:,:,:,cm_index,:)
597            WRITE(iounit) drhodzn(:,:,:,cm_index,:)
598            WRITE(iounit) tthn(:,:,:,cm_index,:)
599            WRITE(iounit) qvhn(:,:,:,cm_index,:)
600
601            ! 2d nested fields
602            WRITE(iounit) psn(:,:,:,cm_index,:)
603            WRITE(iounit) sdn(:,:,:,cm_index,:)
604            WRITE(iounit) msln(:,:,:,cm_index,:)
605            WRITE(iounit) tccn(:,:,:,cm_index,:)
606            WRITE(iounit) u10n(:,:,:,cm_index,:)
607            WRITE(iounit) v10n(:,:,:,cm_index,:)
608            WRITE(iounit) tt2n(:,:,:,cm_index,:)
609            WRITE(iounit) td2n(:,:,:,cm_index,:)
610            WRITE(iounit) lsprecn(:,:,:,cm_index,:)
611            WRITE(iounit) convprecn(:,:,:,cm_index,:)
612            WRITE(iounit) sshfn(:,:,:,cm_index,:)
613            WRITE(iounit) ssrn(:,:,:,cm_index,:)
614            WRITE(iounit) surfstrn(:,:,:,cm_index,:)
615            WRITE(iounit) ustarn(:,:,:,cm_index,:)
616            WRITE(iounit) wstarn(:,:,:,cm_index,:)
617            WRITE(iounit) hmixn(:,:,:,cm_index,:)
618            WRITE(iounit) tropopausen(:,:,:,cm_index,:)
619            WRITE(iounit) olin(:,:,:,cm_index,:)
620            WRITE(iounit) diffkn(:,:,:,cm_index,:)
621            WRITE(iounit) vdepn(:,:,:,cm_index,:)
622
623            ! Auxiliary variables for nests
624            WRITE(iounit) xresoln(:)
625            WRITE(iounit) yresoln(:)
626            WRITE(iounit) xln(:)
627            WRITE(iounit) yln(:)
628            WRITE(iounit) xrn(:)
629            WRITE(iounit) yrn(:)
630
631            ! Variables for polar stereographic projection
632            WRITE(iounit) xglobal, sglobal, nglobal
633            WRITE(iounit) switchnorthg, switchsouthg
634            WRITE(iounit) southpolemap(:)
635            WRITE(iounit) northpolemap(:)
636
637            ! Variables declared in conv_mod (convection)
638            WRITE(iounit) pconv(:)
639            WRITE(iounit) phconv(:)
640            WRITE(iounit) dpr(:)
641            WRITE(iounit) pconv_hpa(:)
642            WRITE(iounit) phconv_hpa(:)
643            WRITE(iounit) ft(:)
644            WRITE(iounit) fq(:)
645            WRITE(iounit) fmass(:,:)
646            WRITE(iounit) sub(:)
647            WRITE(iounit) fmassfrac(:,:)
648            WRITE(iounit) cbaseflux(:,:)
649            WRITE(iounit) cbasefluxn(:,:,:)
650            WRITE(iounit) tconv(:)
651            WRITE(iounit) qconv(:)
652            WRITE(iounit) qsconv(:)
653            WRITE(iounit) psconv, tt2conv, td2conv
654            WRITE(iounit) nconvlev, nconvtop
655
656        ELSE IF (op == 'LOAD') THEN
657
658            ! Read the preprocessed format version string and insure it
659            ! matches this version
660            READ (iounit) temp_preproc_format_version_str
661            PRINT *, 'Reading preprocessed file format version: ', &
662&                    temp_preproc_format_version_str
663
664            IF (TRIM(temp_preproc_format_version_str) == &
665&                                        TRIM(PREPROC_FORMAT_VERSION_STR)) THEN
666                CONTINUE
667            ELSE
668                ! PRINT *, ''  GK: causes relocation truncated to fit: R_X86_64_32
669                PRINT *, 'Inconsistent preprocessing format version'
670                PRINT *, 'Expected Version: ', PREPROC_FORMAT_VERSION_STR
671                PRINT *, 'Detected Version: ', temp_preproc_format_version_str
672                ! PRINT *, ''
673                STOP
674            END IF
675
676            ! Read the compiled max dimensions that were dumped from par_mod
677            ! when creating the fp file, so that we can compare against
678            ! current FLEXPART dimensions - they need to be the same, or else
679            ! we abort.
680            READ (iounit) temp_nxmax, temp_nymax, temp_nzmax, &
681&                         temp_nuvzmax, temp_nwzmax
682
683            ! Get dimensions
684            ncret = nf90_inq_dimid(ncid, 'nxmax', nxmax_dimid)
685            ncret = nf90_inquire_dimension(ncid, nxmax_dimid, nxmax_dimname, &
686&                                                temp_nxmax)
687            PRINT *, 'temp_nxmax: ', temp_nxmax
688
689            ncret = nf90_inq_dimid(ncid, 'nymax', nymax_dimid)
690            ncret = nf90_inquire_dimension(ncid, nymax_dimid, nymax_dimname, &
691&                                                temp_nymax)
692            PRINT *, 'temp_nymax: ', temp_nymax
693
694            ncret = nf90_inq_dimid(ncid, 'nzmax', nzmax_dimid)
695            ncret = nf90_inquire_dimension(ncid, nzmax_dimid, nzmax_dimname, &
696&                                                temp_nzmax)
697            PRINT *, 'temp_nzmax: ', temp_nzmax
698
699            ncret = nf90_inq_dimid(ncid, 'nuvzmax', nuvzmax_dimid)
700            ncret = nf90_inquire_dimension(ncid, nuvzmax_dimid, nuvzmax_dimname, &
701&                                                temp_nuvzmax)
702            PRINT *, 'temp_nuvzmax: ', temp_nuvzmax
703
704            ncret = nf90_inq_dimid(ncid, 'nwzmax', nwzmax_dimid)
705            ncret = nf90_inquire_dimension(ncid, nwzmax_dimid, nwzmax_dimname, &
706&                                                temp_nwzmax)
707            PRINT *, 'temp_nwzmax: ', temp_nwzmax
708
709
710
711            IF ( (temp_nxmax == nxmax) .AND. (temp_nymax == nymax) .AND. &
712&                   (temp_nzmax == nzmax) .AND. &
713&                   (temp_nuvzmax == nuvzmax) .AND. &
714&                   (temp_nwzmax == nwzmax) ) THEN
715                CONTINUE
716            ELSE
717                PRINT *, 'Incompatible dimensions between fp file and current FLEXPART!'
718                ! PRINT *, ''
719                PRINT *, '                  FP file     Compiled FP'
720                PRINT *, 'nxmax:     ', temp_nxmax, '    ', nxmax
721                PRINT *, 'nymax:     ', temp_nymax, '    ', nymax
722                PRINT *, 'nzmax:     ', temp_nzmax, '    ', nzmax
723                PRINT *, 'nuvzmax:     ', temp_nuvzmax, '    ', nuvzmax
724                PRINT *, 'nwzmax:     ', temp_nwzmax, '    ', nwzmax
725                ! PRINT *, ''
726                STOP
727            END IF
728
729
730            ! Scalar values
731            READ(iounit) nx, ny, nxmin1, nymin1, nxfield
732            READ(iounit) nuvz, nwz, nz, nmixz, nlev_ec
733            READ(iounit) dx, dy, xlon0, ylat0, dxconst, dyconst
734
735
736
737            ! Get the varid , then read into scalar variable
738            ncret = nf90_inq_varid(ncid, 'nx', ncvarid)
739            ncret = nf90_get_var(ncid, ncvarid, nx)
740
741            ncret = nf90_inq_varid(ncid, 'ny', ncvarid)
742            ncret = nf90_get_var(ncid, ncvarid, ny)
743
744            ncret = nf90_inq_varid(ncid, 'nxmin1', ncvarid)
745            ncret = nf90_get_var(ncid, ncvarid, nxmin1)
746
747            ncret = nf90_inq_varid(ncid, 'nymin1', ncvarid)
748            ncret = nf90_get_var(ncid, ncvarid, nymin1)
749
750            ncret = nf90_inq_varid(ncid, 'nxfield', ncvarid)
751            ncret = nf90_get_var(ncid, ncvarid, nxfield)
752
753            ncret = nf90_inq_varid(ncid, 'nuvz', ncvarid)
754            ncret = nf90_get_var(ncid, ncvarid, nuvz)
755
756            ncret = nf90_inq_varid(ncid, 'nwz', ncvarid)
757            ncret = nf90_get_var(ncid, ncvarid, nwz)
758
759            ncret = nf90_inq_varid(ncid, 'nz', ncvarid)
760            ncret = nf90_get_var(ncid, ncvarid, nz)
761
762            ncret = nf90_inq_varid(ncid, 'nmixz', ncvarid)
763            ncret = nf90_get_var(ncid, ncvarid, nmixz)
764
765            ncret = nf90_inq_varid(ncid, 'nlev_ec', ncvarid)
766            ncret = nf90_get_var(ncid, ncvarid, nlev_ec)
767
768            ncret = nf90_inq_varid(ncid, 'dx', ncvarid)
769            ncret = nf90_get_var(ncid, ncvarid, dx)
770
771            ncret = nf90_inq_varid(ncid, 'dy', ncvarid)
772            ncret = nf90_get_var(ncid, ncvarid, dy)
773
774            ncret = nf90_inq_varid(ncid, 'xlon0', ncvarid)
775            ncret = nf90_get_var(ncid, ncvarid, xlon0)
776
777            ncret = nf90_inq_varid(ncid, 'ylat0', ncvarid)
778            ncret = nf90_get_var(ncid, ncvarid, ylat0)
779
780            ncret = nf90_inq_varid(ncid, 'dxconst', ncvarid)
781            ncret = nf90_get_var(ncid, ncvarid, dxconst)
782
783            ncret = nf90_inq_varid(ncid, 'dyconst', ncvarid)
784            ncret = nf90_get_var(ncid, ncvarid, dyconst)
785
786
787
788
789
790
791            ! Fixed fields, static in time
792            READ(iounit) oro, excessoro, lsm, xlanduse, height
793
794            ! 3d fields
795            READ(iounit) uu(:,:,:,cm_index)
796            READ(iounit) vv(:,:,:,cm_index)
797            READ(iounit) uupol(:,:,:,cm_index)
798            READ(iounit) vvpol(:,:,:,cm_index)
799            READ(iounit) ww(:,:,:,cm_index)
800            READ(iounit) tt(:,:,:,cm_index)
801            READ(iounit) qv(:,:,:,cm_index)
802            READ(iounit) pv(:,:,:,cm_index)
803            READ(iounit) rho(:,:,:,cm_index)
804            READ(iounit) drhodz(:,:,:,cm_index)
805            READ(iounit) tth(:,:,:,cm_index)
806            READ(iounit) qvh(:,:,:,cm_index)
807            READ(iounit) pplev(:,:,:,cm_index)
808            READ(iounit) clouds(:,:,:,cm_index)
809            READ(iounit) cloudsh(:,:,cm_index)
810
811
812
813
814            ! Get the varid and read the variable into the array
815            ncret = nf90_inq_varid(ncid, 'uu', ncvarid)
816            ncret = nf90_get_var(ncid, ncvarid, uu)
817
818            ncret = nf90_inq_varid(ncid, 'vv', ncvarid)
819            ncret = nf90_get_var(ncid, ncvarid, vv)
820
821            ncret = nf90_inq_varid(ncid, 'uupol', ncvarid)
822            ncret = nf90_get_var(ncid, ncvarid, uupol)
823
824            ncret = nf90_inq_varid(ncid, 'vvpol', ncvarid)
825            ncret = nf90_get_var(ncid, ncvarid, vvpol)
826
827            ncret = nf90_inq_varid(ncid, 'ww', ncvarid)
828            ncret = nf90_get_var(ncid, ncvarid, ww)
829
830            ncret = nf90_inq_varid(ncid, 'tt', ncvarid)
831            ncret = nf90_get_var(ncid, ncvarid, tt)
832
833            ncret = nf90_inq_varid(ncid, 'qv', ncvarid)
834            ncret = nf90_get_var(ncid, ncvarid, qv)
835
836            ncret = nf90_inq_varid(ncid, 'pv', ncvarid)
837            ncret = nf90_get_var(ncid, ncvarid, pv)
838
839            ncret = nf90_inq_varid(ncid, 'rho', ncvarid)
840            ncret = nf90_get_var(ncid, ncvarid, rho)
841
842            ncret = nf90_inq_varid(ncid, 'drhodz', ncvarid)
843            ncret = nf90_get_var(ncid, ncvarid, drhodz)
844
845            ncret = nf90_inq_varid(ncid, 'tth', ncvarid)
846            ncret = nf90_get_var(ncid, ncvarid, tth)
847
848            ncret = nf90_inq_varid(ncid, 'qvh', ncvarid)
849            ncret = nf90_get_var(ncid, ncvarid, qvh)
850
851            ncret = nf90_inq_varid(ncid, 'pplev', ncvarid)
852            ncret = nf90_get_var(ncid, ncvarid, pplev)
853
854
855
856
857
858
859
860
861            ! 2d fields
862            READ(iounit) ps(:,:,:,cm_index)
863            READ(iounit) sd(:,:,:,cm_index)
864            READ(iounit) msl(:,:,:,cm_index)
865            READ(iounit) tcc(:,:,:,cm_index)
866            READ(iounit) u10(:,:,:,cm_index)
867            READ(iounit) v10(:,:,:,cm_index)
868            READ(iounit) tt2(:,:,:,cm_index)
869            READ(iounit) td2(:,:,:,cm_index)
870            READ(iounit) lsprec(:,:,:,cm_index)
871            READ(iounit) convprec(:,:,:,cm_index)
872            READ(iounit) sshf(:,:,:,cm_index)
873            READ(iounit) ssr(:,:,:,cm_index)
874            READ(iounit) surfstr(:,:,:,cm_index)
875            READ(iounit) ustar(:,:,:,cm_index)
876            READ(iounit) wstar(:,:,:,cm_index)
877            READ(iounit) hmix(:,:,:,cm_index)
878            READ(iounit) tropopause(:,:,:,cm_index)
879            READ(iounit) oli(:,:,:,cm_index)
880            READ(iounit) diffk(:,:,:,cm_index)
881            READ(iounit) vdep(:,:,:,cm_index)
882
883            ! 1d fields
884            READ(iounit) z0(:)
885            READ(iounit) akm(:)
886            READ(iounit) bkm(:)
887            READ(iounit) akz(:)
888            READ(iounit) bkz(:)
889            READ(iounit) aknew(:)
890            READ(iounit) bknew(:)
891
892
893            ! Nested, scalar values (for each nest)
894            READ(iounit) nxn(:)
895            READ(iounit) nyn(:)
896            READ(iounit) dxn(:)
897            READ(iounit) dyn(:)
898            READ(iounit) xlon0n(:)
899            READ(iounit) ylat0n(:)
900
901
902            ! Nested fields, static over time
903            READ(iounit) oron, excessoron, lsmn, xlandusen
904
905            ! 3d nested fields
906            READ(iounit) uun(:,:,:,cm_index,:)
907            READ(iounit) vvn(:,:,:,cm_index,:)
908            READ(iounit) wwn(:,:,:,cm_index,:)
909            READ(iounit) ttn(:,:,:,cm_index,:)
910            READ(iounit) qvn(:,:,:,cm_index,:)
911            READ(iounit) pvn(:,:,:,cm_index,:)
912            READ(iounit) cloudsn(:,:,:,cm_index,:)
913            READ(iounit) cloudsnh(:,:,cm_index,:)
914            READ(iounit) rhon(:,:,:,cm_index,:)
915            READ(iounit) drhodzn(:,:,:,cm_index,:)
916            READ(iounit) tthn(:,:,:,cm_index,:)
917            READ(iounit) qvhn(:,:,:,cm_index,:)
918
919            ! 2d nested fields
920            READ(iounit) psn(:,:,:,cm_index,:)
921            READ(iounit) sdn(:,:,:,cm_index,:)
922            READ(iounit) msln(:,:,:,cm_index,:)
923            READ(iounit) tccn(:,:,:,cm_index,:)
924            READ(iounit) u10n(:,:,:,cm_index,:)
925            READ(iounit) v10n(:,:,:,cm_index,:)
926            READ(iounit) tt2n(:,:,:,cm_index,:)
927            READ(iounit) td2n(:,:,:,cm_index,:)
928            READ(iounit) lsprecn(:,:,:,cm_index,:)
929            READ(iounit) convprecn(:,:,:,cm_index,:)
930            READ(iounit) sshfn(:,:,:,cm_index,:)
931            READ(iounit) ssrn(:,:,:,cm_index,:)
932            READ(iounit) surfstrn(:,:,:,cm_index,:)
933            READ(iounit) ustarn(:,:,:,cm_index,:)
934            READ(iounit) wstarn(:,:,:,cm_index,:)
935            READ(iounit) hmixn(:,:,:,cm_index,:)
936            READ(iounit) tropopausen(:,:,:,cm_index,:)
937            READ(iounit) olin(:,:,:,cm_index,:)
938            READ(iounit) diffkn(:,:,:,cm_index,:)
939            READ(iounit) vdepn(:,:,:,cm_index,:)
940
941            ! Auxiliary variables for nests
942            READ(iounit) xresoln(:)
943            READ(iounit) yresoln(:)
944            READ(iounit) xln(:)
945            READ(iounit) yln(:)
946            READ(iounit) xrn(:)
947            READ(iounit) yrn(:)
948
949            ! Variables for polar stereographic projection
950            READ(iounit) xglobal, sglobal, nglobal
951            READ(iounit) switchnorthg, switchsouthg
952            READ(iounit) southpolemap(:)
953            READ(iounit) northpolemap(:)
954
955            ! Variables declared in conv_mod (convection)
956            READ(iounit) pconv(:)
957            READ(iounit) phconv(:)
958            READ(iounit) dpr(:)
959            READ(iounit) pconv_hpa(:)
960            READ(iounit) phconv_hpa(:)
961            READ(iounit) ft(:)
962            READ(iounit) fq(:)
963            READ(iounit) fmass(:,:)
964            READ(iounit) sub(:)
965            READ(iounit) fmassfrac(:,:)
966            READ(iounit) cbaseflux(:,:)
967            READ(iounit) cbasefluxn(:,:,:)
968            READ(iounit) tconv(:)
969            READ(iounit) qconv(:)
970            READ(iounit) qsconv(:)
971            READ(iounit) psconv, tt2conv, td2conv
972            READ(iounit) nconvlev, nconvtop
973
974        ELSE
975            STOP 'fpio(): Illegal operation'
976           
977        ENDIF
978    END SUBROUTINE fpio
979
980    SUBROUTINE fpmetbinary_filetext(filename, cm_index)
981
982        ! This is a utility subroutine meant to be used for testing purposes.
983        ! It facilitates the text output of variables read in from the
984        ! specified .fp file.  This routine will easily cause the program
985        ! to crash due memory allocation issues, particularly when you are
986        ! trying to text print 3d arrays in a single formatted statetment.
987       
988        CHARACTER(LEN=*), INTENT(IN) :: filename
989        INTEGER, INTENT(IN) :: cm_index           ! Index of last dimension in
990                                                  ! most com_mod variables.
991                                                  ! Should be 1 or 2
992
993        !OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', status='REPLACE', &
994        !    form="formatted", access="stream")
995        OPEN(IOUNIT_TEXTOUT, file=filename, action='WRITE', &
996             form="formatted", access="APPEND")
997
998        WRITE(IOUNIT_TEXTOUT, *) 'oro: ', oro
999        WRITE(IOUNIT_TEXTOUT, *) 'excessoro: ', excessoro
1000        WRITE(IOUNIT_TEXTOUT, *) 'lsm: ', lsm
1001        WRITE(IOUNIT_TEXTOUT, *) 'xlanduse: ', xlanduse
1002        WRITE(IOUNIT_TEXTOUT, *) 'height: ', height
1003
1004        WRITE(IOUNIT_TEXTOUT, *) 'uu: ', uu(:,:,:,cm_index)
1005        WRITE(IOUNIT_TEXTOUT, *) 'vv: ', vv(:,:,:,cm_index)
1006        WRITE(IOUNIT_TEXTOUT, *) 'uupol: ', uupol(:,:,:,cm_index)
1007        WRITE(IOUNIT_TEXTOUT, *) 'vvpol: ', vvpol(:,:,:,cm_index)
1008        WRITE(IOUNIT_TEXTOUT, *) 'ww: ', ww(:,:,:,cm_index)
1009        WRITE(IOUNIT_TEXTOUT, *) 'tt: ', tt(:,:,:,cm_index)
1010        WRITE(IOUNIT_TEXTOUT, *) 'qv: ', qv(:,:,:,cm_index)
1011        WRITE(IOUNIT_TEXTOUT, *) 'pv: ', pv(:,:,:,cm_index)
1012        WRITE(IOUNIT_TEXTOUT, *) 'rho: ', rho(:,:,:,cm_index)
1013        WRITE(IOUNIT_TEXTOUT, *) 'drhodz: ', drhodz(:,:,:,cm_index)
1014        WRITE(IOUNIT_TEXTOUT, *) 'tth: ', tth(:,:,:,cm_index)
1015        WRITE(IOUNIT_TEXTOUT, *) 'qvh: ', qvh(:,:,:,cm_index)
1016        WRITE(IOUNIT_TEXTOUT, *) 'pplev: ', pplev(:,:,:,cm_index)
1017        WRITE(IOUNIT_TEXTOUT, *) 'clouds: ', clouds(:,:,:,cm_index)
1018        WRITE(IOUNIT_TEXTOUT, *) 'cloudsh: ', cloudsh(:,:,cm_index)
1019
1020
1021
1022
1023        CLOSE(IOUNIT_TEXTOUT)
1024    END SUBROUTINE fpmetbinary_filetext
1025
1026
1027END MODULE fpmetbinary_mod
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG