source: palm/tags/release-3.4/UTIL/combine_plot_fields.f90 @ 818

Last change on this file since 818 was 108, checked in by letzel, 17 years ago
  • Improved coupler: evaporation - salinity-flux coupling for humidity = .T.,

avoid MPI hangs when coupled runs terminate, add DOC/app/chapter_3.8;

  • Optional calculation of km and kh from initial TKE e_init;
  • Default initialization of km,kh = 0.00001 for ocean = .T.;
  • Allow data_output_pr= q, wq, w"q", w*q* for humidity = .T.;
  • Bugfix: Rayleigh damping for ocean fixed.
File size: 25.4 KB
Line 
1 PROGRAM combine_plot_fields
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
6! loop for processing of output by coupled runs, id_string does not contain
7! modus any longer
8!
9! Former revisions:
10! -----------------
11! 18/01/06  Output of time-averaged data
12!
13! 25/05/05  Errors removed
14!
15! 26/04/05  Output in NetCDF format, iso2d and avs output only if parameter
16!           file exists
17!
18! 31/10/01  All comments and messages translated into English
19!
20! 23/02/99  Keine Bearbeitung komprimierter 3D-Daten
21! Ursprungsversion vom 28/07/97
22!
23!
24! Description:
25! ------------
26! This routine combines data of the PALM-subdomains into one file. In PALM
27! every processor element opens its own file and writes 2D- or 3D-binary-data
28! into it (different files are opened for xy-, xz-, yz-cross-sections and
29! 3D-data). For plotting or analyzing these PE-data have to be collected and
30! to be put into single files, which is done by this routine.
31! Output format is NetCDF. Additionally, a data are output in a binary format
32! readable by ISO2D-software (cross-sections) and by AVS (3D-data).
33!
34! This routine must be compiled with:
35! decalpha:
36!    f95 -cpp -D__netcdf -fast -r8 -I/usr/local/netcdf-3.5.1/include
37!    -L/usr/local/netcdf-3.5.1/lib -lnetcdf
38! IBM-Regatta:
39!    xlf95 -qsuffix=cpp=f90 -WF,-D__netcdf -qrealsize=8 -q64 -qmaxmem=-1 -Q
40!    -I /aws/dataformats/netcdf-3.6.0-p1/64-32/include-64
41!    -L/aws/dataformats/netcdf-3.6.0-p1/64-32/lib -lnetcdf -O3
42! IBM-Regatta KISTI:
43!    xlf95 -qsuffix=cpp=f90 -WF,-D__netcdf -qrealsize=8 -q64 -qmaxmem=-1 -Q
44!    -I /applic/netcdf64/src/f90
45!    -L/applic/lib/NETCDF64 -lnetcdf -O3
46! IBM-Regatta Yonsei (gfdl5):
47!    xlf95 -qsuffix=cpp=f90 -WF,-D__netcdf -qrealsize=8 -q64 -qmaxmem=-1 -Q
48!    -I /usr1/users/raasch/pub/netcdf-3.6.0-p1/include
49!    -L/usr1/users/raasch/pub/netcdf-3.6.0-p1/lib -lnetcdf -O3
50! IMUK:
51!    ifort combine...f90 -o combine...x
52!    -cpp -D__netcdf -I /muksoft/packages/netcdf/linux/include -axW -r8 -nbs
53!    -Vaxlib -L /muksoft/packages/netcdf/linux/lib -lnetcdf
54! NEC-SX6:
55!    sxf90 combine...f90 -o combine...x
56!    -I /pool/SX-6/netcdf/netcdf-3.6.0-p1/include  -C hopt -Wf '-A idbl4'
57!    -D__netcdf -L/pool/SX-6/netcdf/netcdf-3.6.0-p1/lib -lnetcdf
58! Sun Fire X4600
59!    pgf95 combine...f90 -o combine...x
60!    -Mpreprocess -D__netcdf -I /home/usr5/mkanda/netcdf-3.6.1/src/f90 -r8
61!    -fast -L/home/usr5/mkanda/netcdf-3.6.1/src/libsrc -lnetcdf
62! FIMM:
63!    ifort combine...f90 -o combine...x
64!    -axW -cpp -openmp -r8 -nbs -convert little_endian -D__netcdf
65!    -I /local/netcdf/include -Vaxlib -L/local/netcdf/lib -lnetcdf
66!------------------------------------------------------------------------------!
67
68#if defined( __netcdf )
69    USE netcdf
70#endif
71
72    IMPLICIT NONE
73
74!
75!-- Local variables
76    CHARACTER (LEN=2)    ::  modus, model_string
77    CHARACTER (LEN=4)    ::  id_string
78    CHARACTER (LEN=10)   ::  dimname, var_name
79    CHARACTER (LEN=40)   ::  filename
80
81    CHARACTER (LEN=2000), DIMENSION(0:1) ::  var_list
82
83    INTEGER, PARAMETER ::  spk = SELECTED_REAL_KIND( 6 )
84
85    INTEGER ::  av, danz, i, id,             &
86                j, model, models, nc_stat, nxa, nxag, nxe, nxeg, nya,   &
87                nyag, nye, nyeg, nza, nzag, nze, nzeg, pos, time_step, xa, xe, &
88                ya, ye, za, ze
89
90    INTEGER, DIMENSION(0:1) ::  current_level, current_var, fanz, id_set, &
91         id_var_time, num_var
92
93    INTEGER, DIMENSION(4) ::  id_dims_loc
94
95    INTEGER, DIMENSION(0:1,4) ::  id_dims
96
97    INTEGER, DIMENSION(0:1,1000) ::  id_var, levels
98
99    LOGICAL ::  avs_output, compressed, found, iso2d_output, netcdf_output, &
100                netcdf_0, netcdf_1
101
102    REAL ::  dx, simulated_time
103    REAL, DIMENSION(:),   ALLOCATABLE   ::  eta, ho, hu
104    REAL, DIMENSION(:,:), ALLOCATABLE   ::  pf
105    REAL(spk), DIMENSION(:,:,:), ALLOCATABLE ::  pf3d
106
107    PRINT*, ''
108    PRINT*, '*** combine_plot_fields ***'
109    INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found )
110    IF ( found )  THEN
111       models = 2
112       PRINT*, '    coupled run'
113    ELSE
114       models = 1
115       PRINT*, '    uncoupled run'
116    ENDIF
117    DO model = 1, models
118       IF ( models == 2 )  THEN
119          PRINT*, ''
120          PRINT*, ''
121          PRINT*, '*** combine_plot_fields ***'
122          IF ( model == 2 )  THEN
123             model_string = "_O"
124             PRINT*, '    now combining ocean data'
125             PRINT*, '    ========================'
126          ELSE
127             model_string = ""
128             PRINT*, '    now combining atmosphere data'
129             PRINT*, '    ============================='
130          ENDIF
131       ENDIF
132!
133!--    2D-arrays for ISO2D
134!--    Main loop for the three different cross-sections, starting with
135!--    xy-section
136       modus = 'XY'
137       PRINT*, ''
138       DO  WHILE ( modus == 'XY'  .OR.  modus == 'XZ'  .OR.  modus == 'YZ' )
139!
140!--       Check, if file from PE0 exists. If it does not exist, PALM did not
141!--       create any output for this cross-section.
142          danz = 0
143          WRITE (id_string,'(I4.4)')  danz
144          INQUIRE ( &
145               FILE='PLOT2D_'//modus//TRIM( model_string )//'_'//id_string, &
146               EXIST=found )
147!
148!--       Find out the number of files (equal to the number of PEs which
149!--       have been used in PALM) and open them
150          DO  WHILE ( found )
151
152             OPEN ( danz+110, &
153                  FILE='PLOT2D_'//modus//TRIM( model_string )//'_'//id_string, &
154                  FORM='UNFORMATTED' )
155             danz = danz + 1
156             WRITE (id_string,'(I4.4)')  danz
157             INQUIRE ( &
158                  FILE='PLOT2D_'//modus//TRIM( model_string )//'_'//id_string, &
159                  EXIST=found )
160
161          ENDDO
162
163!
164!--       Inquire whether an iso2d parameter file exists
165          INQUIRE( FILE='PLOT2D_'//modus//'_GLOBAL'//TRIM( model_string ), &
166               EXIST=iso2d_output )
167
168!
169!--       Inquire whether a NetCDF file exists
170          INQUIRE( FILE='DATA_2D_'//modus//'_NETCDF'//TRIM( model_string ), &
171               EXIST=netcdf_0 )
172
173!
174!--       Inquire whether a NetCDF file for time-averaged data exists
175          INQUIRE( FILE='DATA_2D_'//modus//'_AV_NETCDF'//TRIM( model_string ),&
176               EXIST=netcdf_1 )
177
178          IF ( netcdf_0  .OR.  netcdf_1 )  THEN
179             netcdf_output = .TRUE.
180          ELSE
181             netcdf_output = .FALSE.
182          ENDIF
183
184!
185!--       Info-output
186          PRINT*, ''
187          PRINT*, '*** combine_plot_fields ***'
188#if defined( __netcdf )
189          IF ( netcdf_output )  PRINT*, '    NetCDF output enabled'
190#else
191          IF ( netcdf_output )  THEN
192             PRINT*, '--- Sorry, no NetCDF support on this host'
193             netcdf_output = .FALSE.
194          ENDIF
195#endif
196          IF ( danz /= 0 )  THEN
197             PRINT*, '    ',modus,'-section:  ', danz, ' file(s) found'
198          ELSE
199             PRINT*, '    no ', modus, '-section data available'
200          ENDIF
201
202          IF ( netcdf_output  .AND.  danz /= 0 )  THEN
203#if defined( __netcdf )
204             DO  av = 0, 1
205
206                IF ( av == 0  .AND.  .NOT.  netcdf_0 )  CYCLE
207                IF ( av == 1  .AND.  .NOT.  netcdf_1 )  CYCLE
208
209!
210!--             Open NetCDF dataset
211                IF ( av == 0 )  THEN
212                   filename = 'DATA_2D_'//modus//'_NETCDF' &
213                        //TRIM( model_string )
214                ELSE
215                   filename = 'DATA_2D_'//modus//'_AV_NETCDF' &
216                        //TRIM( model_string )
217                ENDIF
218                nc_stat = NF90_OPEN( filename, NF90_WRITE, id_set(av) )
219                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 1 )
220
221!
222!--             Get the list of variables (order of variables corresponds with
223!--             the order of data on the binary file)
224                var_list(av) = ' '    ! GET_ATT does not assign trailing blanks
225                nc_stat = NF90_GET_ATT( id_set(av), NF90_GLOBAL, 'VAR_LIST', &
226                     var_list(av) )
227                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 2 )
228
229!
230!--             Inquire id of the time coordinate variable
231                nc_stat = NF90_INQ_VARID( id_set(av), 'time', id_var_time(av) )
232                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 3 )
233
234!
235!--             Count number of variables; there is one more semicolon in the
236!--             string than variable names
237                num_var(av) = -1
238                DO  i = 1, LEN( var_list(av) )
239                   IF ( var_list(av)(i:i) == ';' )  num_var(av) = num_var(av) +1
240                ENDDO
241
242!
243!--             Extract the variable names from the list and inquire their
244!--             NetCDF IDs
245                pos = INDEX( var_list(av), ';' )
246!
247!--             Loop over all variables
248                DO  i = 1, num_var(av)
249
250!
251!--                Extract variable name from list
252                   var_list(av) = var_list(av)(pos+1:)
253                   pos = INDEX( var_list(av), ';' )
254                   var_name = var_list(av)(1:pos-1)
255
256!
257!--                Get variable ID from name
258                   nc_stat = NF90_INQ_VARID( id_set(av), TRIM( var_name ), &
259                        id_var(av,i) )
260                   IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 4 )
261
262!
263!--                Get number of x/y/z levels for that variable
264                   nc_stat = NF90_INQUIRE_VARIABLE( id_set(av), id_var(av,i), &
265                        dimids = id_dims_loc )
266                   IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 5 )
267                   id_dims(av,:) = id_dims_loc
268
269!
270!--                Inquire dimension ID
271                   DO  j = 1, 4
272                      nc_stat = NF90_INQUIRE_DIMENSION( id_set(av), &
273                           id_dims(av,j), dimname, levels(av,i) )
274                      IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 6 )
275
276                      IF ( modus == 'XY' .AND. INDEX(dimname, 'z') /= 0 )  EXIT
277                      IF ( modus == 'XZ' .AND. INDEX(dimname, 'y') /= 0 )  EXIT
278                      IF ( modus == 'YZ' .AND. INDEX(dimname, 'x') /= 0 )  EXIT
279                   ENDDO
280
281                ENDDO
282
283             ENDDO   ! av = 0, 1
284
285          ENDIF
286#endif
287
288!
289!--       Read the arrays, as long as the end of the file is reached
290          fanz          =         0
291          current_level =         1
292          current_var   = 999999999
293
294          DO  WHILE ( danz /= 0 )
295
296!
297!--          Loop over all files (reading data of the subdomains)
298             DO  id = 0, danz-1
299!
300!--             File from PE0 contains special information at the beginning,
301!--             concerning the lower and upper indices of the total-domain used
302!--             in PALM (nxag, nxeg, nyag, nyeg) and the lower and upper indices
303!--             of the array to be writte by this routine (nxa, nxe, nya,
304!--             nye). Usually in the horizontal directions nxag=-1 and nxa=0
305!--             while all other variables have the same value (i.e. nxeg=nxe).
306!--             Allocate necessary arrays, open the output file and write
307!--             the coordinate informations needed by ISO2D.
308                IF ( id == 0  .AND.  fanz(0) == 0  .AND.  fanz(1) == 0 )  THEN
309                   READ ( id+110 )  nxag, nxeg, nyag, nyeg
310                   READ ( id+110 )  nxa, nxe, nya, nye
311                   ALLOCATE ( eta(nya:nye), ho(nxa:nxe), hu(nxa:nxe), &
312                        pf(nxag:nxeg,nyag:nyeg) )
313                   READ ( id+110 )  dx, eta, hu, ho
314
315                   IF ( iso2d_output )  THEN
316                      OPEN ( 2, FILE='PLOT2D_'//modus//TRIM( model_string ), &
317                           FORM='UNFORMATTED' )
318                      WRITE ( 2 )  dx, eta, hu, ho
319                   ENDIF
320                ENDIF
321!
322!--             Read output time
323                IF ( netcdf_output  .AND.  id == 0 )  THEN
324                   IF ( netcdf_1 )  THEN
325                      READ ( id+110, END=998 )  simulated_time, time_step, av
326                   ELSE
327!
328!--                   For compatibility with earlier PALM versions
329                      READ ( id+110, END=998 )  simulated_time, time_step
330                      av = 0
331                   ENDIF
332                ENDIF
333!
334!--             Read subdomain indices
335                READ ( id+110, END=998 )  xa, xe, ya, ye
336!
337!--             IF the PE made no output (in case that no part of the
338!--             cross-section is situated on this PE), indices have the
339!--             value -1
340                IF ( .NOT. ( xa == -1  .AND.  xe == -1  .AND. &
341                             ya == -1  .AND.  ye == -1 ) )  THEN
342!
343!--                Read the subdomain grid-point values
344                   READ ( id+110 )  pf(xa:xe,ya:ye)
345                ENDIF
346                IF ( id == 0 )  fanz(av) = fanz(av) + 1
347
348             ENDDO
349!
350!--          Write the data of the total domain cross-section
351             IF ( iso2d_output )  WRITE ( 2 )  pf(nxa:nxe,nya:nye)
352       
353!
354!--          Write same data in NetCDF format
355             IF ( netcdf_output )  THEN
356#if defined( __netcdf )
357!
358!--             Check if a new time step has begun; if yes write data to time
359!--             axis
360                IF ( current_var(av) > num_var(av) )  THEN
361                   current_var(av) = 1
362                   nc_stat = NF90_PUT_VAR( id_set(av), id_var_time(av), &
363                        (/ simulated_time /),        &
364                        start = (/ time_step /),     &
365                        count = (/ 1 /) )
366                   IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 7 )
367                ENDIF
368
369!
370!--             Now write the data; this is mode dependent
371                SELECT CASE ( modus )
372
373                   CASE ( 'XY' )
374                      nc_stat = NF90_PUT_VAR( id_set(av),                      &
375                                           id_var(av,current_var(av)),         &
376                                           pf(nxa:nxe,nya:nye),                &
377                             start = (/ 1, 1, current_level(av), time_step /), &
378                                      count = (/ nxe-nxa+1, nye-nya+1, 1, 1 /) )
379                      IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 8)
380                 
381                   CASE ( 'XZ' )
382                      nc_stat = NF90_PUT_VAR( id_set(av),                      &
383                                           id_var(av,current_var(av)),         &
384                                           pf(nxa:nxe,nya:nye),                &
385                             start = (/ 1, current_level(av), 1, time_step /), &
386                                      count = (/ nxe-nxa+1, 1, nye-nya+1, 1 /) )
387                      IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 9)
388
389                   CASE ( 'YZ' )
390                      nc_stat = NF90_PUT_VAR( id_set(av),                      &
391                                           id_var(av,current_var(av)),         &
392                                           pf(nxa:nxe,nya:nye),                &
393                             start = (/ current_level(av), 1, 1, time_step /), &
394                                      count = (/ 1, nxe-nxa+1, nye-nya+1, 1 /) )
395                      IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error(10)
396
397                END SELECT
398
399!
400!--             Data is written, check if max level is reached
401                current_level(av) = current_level(av) + 1
402                IF ( current_level(av) > levels(av,current_var(av)) )  THEN
403                   current_level(av) = 1
404                   current_var(av)   = current_var(av) + 1
405                ENDIF
406
407             ENDIF
408#endif
409
410          ENDDO
411
412998       IF ( danz /= 0 )  THEN
413!
414!--          Print the number of the arrays processed
415             WRITE (*,'(16X,I4,A)')  fanz(0)+fanz(1), ' array(s) processed'
416             IF ( fanz(1) /= 0 )  THEN
417                WRITE (*,'(16X,I4,A)')  fanz(1), ' array(s) are time-averaged'
418             ENDIF
419
420!
421!--          Close all files and deallocate arrays
422             DO  id = 0, danz-1
423                CLOSE ( id+110 )
424             ENDDO
425             CLOSE ( 2 )
426             DEALLOCATE ( eta, ho, hu, pf )
427          ENDIF
428
429!
430!--       Close the NetCDF file
431          IF ( netcdf_output  .AND.  danz /= 0 )  THEN
432#if defined( __netcdf )
433             IF ( netcdf_0 )  THEN
434                nc_stat = NF90_CLOSE( id_set(0) )
435                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 11 )
436             ENDIF
437             IF ( netcdf_1 )  THEN
438                nc_stat = NF90_CLOSE( id_set(1) )
439                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 12 )
440             ENDIF
441#endif
442          ENDIF
443
444!
445!--       Choose the next cross-section
446          SELECT CASE ( modus )
447             CASE ( 'XY' )
448                modus = 'XZ'
449             CASE ( 'XZ' )
450                modus = 'YZ'
451             CASE ( 'YZ' )
452                modus = 'no'
453          END SELECT
454
455       ENDDO
456
457
458!
459!--    Combine the 3D-arrays
460
461!
462!--    Inquire whether an avs fld file exists
463       INQUIRE( FILE='PLOT3D_FLD'//TRIM( model_string ), EXIST=avs_output )
464
465!
466!--    Inquire whether a NetCDF file exists
467       INQUIRE( FILE='DATA_3D_NETCDF'//TRIM( model_string ), EXIST=netcdf_0 )
468
469!
470!--    Inquire whether a NetCDF file for time-averaged data exists
471       INQUIRE( FILE='DATA_3D_AV_NETCDF'//TRIM( model_string ), EXIST=netcdf_1 )
472
473       IF ( netcdf_0  .OR.  netcdf_1 )  THEN
474          netcdf_output = .TRUE.
475       ELSE
476          netcdf_output = .FALSE.
477       ENDIF
478
479!
480!--    Check, if file from PE0 exists
481       danz = 0
482       WRITE (id_string,'(I4.4)')  danz
483       INQUIRE ( &
484            FILE='PLOT3D_DATA'//TRIM( model_string )//'_'//TRIM( id_string ),  &
485            EXIST=found )
486
487!
488!--    Combination only works, if data are not compressed. In that case,
489!--    PALM created a flag file (PLOT3D_COMPRESSED)
490       INQUIRE ( FILE='PLOT3D_COMPRESSED'//TRIM( model_string ), &
491            EXIST=compressed )
492
493!
494!--    Find out the number of files and open them
495       DO  WHILE ( found  .AND.  .NOT. compressed )
496
497          OPEN ( danz+110, &
498               FILE='PLOT3D_DATA'//TRIM( model_string )//'_'//TRIM(id_string), &
499               FORM='UNFORMATTED')
500          danz = danz + 1
501          WRITE (id_string,'(I4.4)')  danz
502          INQUIRE ( &
503               FILE='PLOT3D_DATA'//TRIM( model_string )//'_'//TRIM(id_string), &
504               EXIST=found )
505
506       ENDDO
507
508!
509!--    Info-output
510       PRINT*, ' '
511       PRINT*, '*** combine_plot_fields ***'
512#if defined( __netcdf )
513       IF ( netcdf_output )  PRINT*, '    NetCDF output enabled'
514#else
515       IF ( netcdf_output )  THEN
516          PRINT*, '--- Sorry, no NetCDF support on this host'
517          netcdf_output = .FALSE.
518       ENDIF
519#endif
520       IF ( danz /= 0 )  THEN
521          PRINT*, '    3D-data:     ', danz, ' file(s) found'
522       ELSE
523          IF ( found .AND. compressed )  THEN
524             PRINT*, '+++ no 3D-data processing, since data are compressed'
525          ELSE
526             PRINT*, '    no 3D-data file available'
527          ENDIF
528       ENDIF
529
530       IF ( netcdf_output  .AND.  danz /= 0 )  THEN
531#if defined( __netcdf )
532          DO  av = 0, 1
533
534             IF ( av == 0  .AND.  .NOT.  netcdf_0 )  CYCLE
535             IF ( av == 1  .AND.  .NOT.  netcdf_1 )  CYCLE
536
537!
538!--          Open NetCDF dataset
539             IF ( av == 0 )  THEN
540                filename = 'DATA_3D_NETCDF'//TRIM( model_string )
541             ELSE
542                filename = 'DATA_3D_AV_NETCDF'//TRIM( model_string )
543             ENDIF
544             nc_stat = NF90_OPEN( filename, NF90_WRITE, id_set(av) )
545             IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 13 )
546
547
548!
549!--          Get the list of variables (order of variables corresponds with the
550!--          order of data on the binary file)
551             var_list(av) = ' '    ! GET_ATT does not assign trailing blanks
552             nc_stat = NF90_GET_ATT( id_set(av), NF90_GLOBAL, 'VAR_LIST', &
553                  var_list(av) )
554             IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 14 )
555
556!
557!--          Inquire id of the time coordinate variable
558             nc_stat = NF90_INQ_VARID( id_set(av), 'time', id_var_time(av) )
559             IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 15 )
560
561!
562!--          Count number of variables; there is one more semicolon in the
563!--          string than variable names
564             num_var(av) = -1
565             DO  i = 1, LEN( var_list(av) )
566                IF ( var_list(av)(i:i) == ';' )  num_var(av) = num_var(av) + 1
567             ENDDO
568
569!
570!--          Extract the variable names from the list and inquire their NetCDF
571!--          IDs
572             pos = INDEX( var_list(av), ';' )
573!
574!--          Loop over all variables
575             DO  i = 1, num_var(av)
576
577!
578!--             Extract variable name from list
579                var_list(av) = var_list(av)(pos+1:)
580                pos = INDEX( var_list(av), ';' )
581                var_name = var_list(av)(1:pos-1)
582
583!
584!--             Get variable ID from name
585!                print*, '*** find id for "',TRIM( var_name ),'" begin'
586                nc_stat = NF90_INQ_VARID( id_set(av), TRIM( var_name ), &
587                     id_var(av,i) )
588                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 16 )
589!                print*, '*** find id for "',TRIM( var_name ),'" end'
590
591             ENDDO
592
593          ENDDO    ! av=0,1
594
595       ENDIF
596#endif
597
598!
599!--    Read arrays, until the end of the file is reached
600       current_var = 999999999
601       fanz = 0
602       DO  WHILE ( danz /= 0 )
603
604!
605!--       Loop over all files
606          DO  id = 0, danz-1
607!
608!--          File from PE0 contains special information at the beginning,
609!--          concerning the lower and upper indices of the total-domain used in
610!--          PALM (nxag, nxeg, nyag, nyeg, nzag, nzeg) and the lower and upper
611!--          indices of the array to be written by this routine (nxa, nxe, nya,
612!--          nye, nza, nze). Usually nxag=-1 and nxa=0, nyag=-1 and nya=0,
613!--          nzeg=nz and nze=nz_plot3d.
614!--          Allocate necessary array and open the output file.
615             IF ( id == 0  .AND.  fanz(0) == 0  .AND.  fanz(1) == 0 )  THEN
616                READ ( id+110 )  nxag, nxeg, nyag, nyeg, nzag, nzeg
617                READ ( id+110 )  nxa, nxe, nya, nye, nza, nze
618                ALLOCATE ( pf3d(nxag:nxeg,nyag:nyeg,nzag:nzeg) )
619                IF ( avs_output )  THEN
620                   OPEN ( 2, FILE='PLOT3D_DATA'//TRIM( model_string ), &
621                        FORM='UNFORMATTED' )
622                ENDIF
623             ENDIF
624
625!
626!--          Read output time
627             IF ( netcdf_output  .AND.  id == 0 )  THEN
628                IF ( netcdf_1 )  THEN
629                   READ ( id+110, END=999 )  simulated_time, time_step, av
630                ELSE
631!
632!--                For compatibility with earlier PALM versions
633                   READ ( id+110, END=999 )  simulated_time, time_step
634                   av = 0
635                ENDIF
636             ENDIF
637
638!
639!--          Read subdomain indices and grid point values
640             READ ( id+110, END=999 )  xa, xe, ya, ye, za, ze
641             READ ( id+110 )  pf3d(xa:xe,ya:ye,za:ze)
642             IF ( id == 0 )  fanz(av) = fanz(av) + 1
643
644          ENDDO
645
646!
647!--       Write data of the total domain
648          IF ( avs_output )  WRITE ( 2 )  pf3d(nxa:nxe,nya:nye,nza:nze)
649       
650!
651!--       Write same data in NetCDF format
652          IF ( netcdf_output )  THEN
653#if defined( __netcdf )
654!
655!--          Check if a new time step has begun; if yes write data to time axis
656             IF ( current_var(av) > num_var(av) )  THEN
657                current_var(av) = 1
658                nc_stat = NF90_PUT_VAR( id_set(av), id_var_time(av), &
659                                     (/ simulated_time /),&
660                                     start = (/ time_step /), count = (/ 1 /) )
661                IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 17 )
662             ENDIF
663
664!
665!--          Now write the data
666             nc_stat = NF90_PUT_VAR( id_set(av), id_var(av,current_var(av)), &
667                                  pf3d(nxa:nxe,nya:nye,nza:nze),      &
668                                  start = (/ 1, 1, 1, time_step /),   &
669                              count = (/ nxe-nxa+1, nye-nya+1, nze-nza+1, 1 /) )
670             IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 18 )
671
672             current_var(av) = current_var(av) + 1
673
674#endif
675          ENDIF
676
677       ENDDO
678
679999    IF ( danz /= 0 )  THEN
680!
681!--       Print the number of arrays processed
682          WRITE (*,'(16X,I4,A)')  fanz(0)+fanz(1), ' array(s) processed'
683          IF ( fanz(1) /= 0 )  THEN
684             WRITE (*,'(16X,I4,A)')  fanz(1), ' array(s) are time-averaged'
685          ENDIF
686!
687!--       Close all files and deallocate array
688          DO  id = 0, danz-1
689             CLOSE ( id+110 )
690          ENDDO
691          CLOSE ( 2 )
692          DEALLOCATE ( pf3d )
693!
694!--       Close the NetCDF file
695          IF ( netcdf_output )  THEN
696#if defined( __netcdf )
697             IF ( netcdf_0 )  THEN
698                nc_stat = NF90_CLOSE( id_set(0) )
699                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 19 )
700             ENDIF
701             IF ( netcdf_1 )  THEN
702                nc_stat = NF90_CLOSE( id_set(1) )
703                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 20 )
704             ENDIF
705#endif
706          ENDIF
707       ENDIF
708
709    ENDDO  ! models
710
711
712 CONTAINS
713
714
715    SUBROUTINE handle_netcdf_error( errno )
716!
717!--    Prints out a text message corresponding to the current NetCDF status
718
719       IMPLICIT NONE
720
721       INTEGER, INTENT(IN) ::  errno
722
723       IF ( nc_stat /= NF90_NOERR )  THEN
724          PRINT*, '+++ combine_plot_fields  netcdf: ', av, errno, &
725                  TRIM( nf90_strerror( nc_stat ) )
726       ENDIF
727
728    END SUBROUTINE handle_netcdf_error
729
730
731 END PROGRAM combine_plot_fields
732
733
734
Note: See TracBrowser for help on using the repository browser.