source: palm/trunk/UTIL/combine_plot_fields.f90 @ 514

Last change on this file since 514 was 494, checked in by raasch, 15 years ago

last commit documented; configuration example file for netcdf4 added

  • Property svn:keywords set to Id
File size: 30.0 KB
RevLine 
[1]1 PROGRAM combine_plot_fields
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
[494]6!
[226]7!
8! Former revisions:
9! -----------------
10! $Id: combine_plot_fields.f90 494 2010-03-01 12:23:32Z heinze $
11!
[494]12! 493 2010-03-01 08:30:24Z raasch
13! Exit in case of already complete NetCDF data (due to parallel output in PALM)
14! cpu measurements included
15!
[226]16! 210 2008-11-06 08:54:02Z raasch
[210]17! Size of pf3d adjusted to the required output size (1 gridpoint less, along
18! all three dimensions), because output of a subset of the data
19! (pf3d(nxa:nxe...) in the NF90_PUT_VAR statement caused segmentation fault
20! with the INTEL compiler.
[191]21! Subdomain data are read into temporary arrays pf_tmp/pf3d_tmp in order to
22! avoid INTEL compiler warnings about (automatic) creation of temporary arrays
[179]23! Bugfix: three misplaced #endif directives
[1]24!
[139]25! 114 2007-10-10 00:03:15Z raasch
26! Bugfix: model_string needed a default value
27!
[114]28! Aug 07    Loop for processing of output by coupled runs, id_string does not
29!           contain modus any longer
30!
[1]31! 18/01/06  Output of time-averaged data
32!
33! 25/05/05  Errors removed
34!
35! 26/04/05  Output in NetCDF format, iso2d and avs output only if parameter
36!           file exists
37!
38! 31/10/01  All comments and messages translated into English
39!
40! 23/02/99  Keine Bearbeitung komprimierter 3D-Daten
41! Ursprungsversion vom 28/07/97
42!
43!
44! Description:
45! ------------
46! This routine combines data of the PALM-subdomains into one file. In PALM
47! every processor element opens its own file and writes 2D- or 3D-binary-data
48! into it (different files are opened for xy-, xz-, yz-cross-sections and
49! 3D-data). For plotting or analyzing these PE-data have to be collected and
50! to be put into single files, which is done by this routine.
51! Output format is NetCDF. Additionally, a data are output in a binary format
52! readable by ISO2D-software (cross-sections) and by AVS (3D-data).
53!------------------------------------------------------------------------------!
54
55#if defined( __netcdf )
56    USE netcdf
57#endif
58
59    IMPLICIT NONE
60
61!
62!-- Local variables
[108]63    CHARACTER (LEN=2)    ::  modus, model_string
64    CHARACTER (LEN=4)    ::  id_string
[1]65    CHARACTER (LEN=10)   ::  dimname, var_name
66    CHARACTER (LEN=40)   ::  filename
67
68    CHARACTER (LEN=2000), DIMENSION(0:1) ::  var_list
69
70    INTEGER, PARAMETER ::  spk = SELECTED_REAL_KIND( 6 )
71
72    INTEGER ::  av, danz, i, id,             &
[210]73                j, k, model, models, nc_stat, nxa, nxag, nxe, nxeg, nya,   &
[1]74                nyag, nye, nyeg, nza, nzag, nze, nzeg, pos, time_step, xa, xe, &
[210]75                xxa, xxe, ya, ye, yya, yye, za, ze, zza, zze
[1]76
[493]77#if defined( __lc ) || defined( __decalpha )
78    INTEGER(8)                  ::  count, count_rate
79#elif defined( __nec )
80    INTEGER                     ::  count, count_rate
81#elif defined( __ibm )
82    INTEGER(8)                  ::  IRTC
83#endif
84
[108]85    INTEGER, DIMENSION(0:1) ::  current_level, current_var, fanz, id_set, &
86         id_var_time, num_var
[1]87
88    INTEGER, DIMENSION(4) ::  id_dims_loc
89
90    INTEGER, DIMENSION(0:1,4) ::  id_dims
91
92    INTEGER, DIMENSION(0:1,1000) ::  id_var, levels
93
94    LOGICAL ::  avs_output, compressed, found, iso2d_output, netcdf_output, &
[493]95                netcdf_parallel, netcdf_0, netcdf_1
[1]96
[493]97    REAL ::  cpu_start_time, cpu_end_time, dx, simulated_time
[1]98    REAL, DIMENSION(:),   ALLOCATABLE   ::  eta, ho, hu
[191]99    REAL, DIMENSION(:,:), ALLOCATABLE   ::  pf, pf_tmp
100    REAL(spk), DIMENSION(:,:,:), ALLOCATABLE ::  pf3d, pf3d_tmp
[1]101
102    PRINT*, ''
[114]103    PRINT*, ''
[108]104    PRINT*, '*** combine_plot_fields ***'
[114]105
106!
107!-- Find out if a coupled run has been carried out
[108]108    INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found )
109    IF ( found )  THEN
110       models = 2
111       PRINT*, '    coupled run'
112    ELSE
113       models = 1
114       PRINT*, '    uncoupled run'
115    ENDIF
[114]116
117!
118!-- Do everything for each model
[108]119    DO model = 1, models
[114]120!
121!--    Set the model string used to identify the filenames
122       model_string = ''
[108]123       IF ( models == 2 )  THEN
124          PRINT*, ''
125          PRINT*, '*** combine_plot_fields ***'
126          IF ( model == 2 )  THEN
[114]127             model_string = '_O'
[108]128             PRINT*, '    now combining ocean data'
129             PRINT*, '    ========================'
130          ELSE
131             PRINT*, '    now combining atmosphere data'
132             PRINT*, '    ============================='
133          ENDIF
134       ENDIF
[1]135!
[108]136!--    2D-arrays for ISO2D
137!--    Main loop for the three different cross-sections, starting with
138!--    xy-section
139       modus = 'XY'
140       PRINT*, ''
141       DO  WHILE ( modus == 'XY'  .OR.  modus == 'XZ'  .OR.  modus == 'YZ' )
[1]142!
[493]143!--       Take current time
144#if defined( __lc ) || defined( __decalpha ) || defined( __nec )
145          CALL SYSTEM_CLOCK( count, count_rate )
146          cpu_start_time = REAL( count ) / REAL( count_rate )
147#elif defined( __ibm )
148          cpu_start_time = IRTC( ) * 1E-9
149#else
150          PRINT*,  '+++ INFORMATIVE: no time measurement defined on this host'
151#endif
152
153          netcdf_parallel = .FALSE.
154!
[108]155!--       Check, if file from PE0 exists. If it does not exist, PALM did not
156!--       create any output for this cross-section.
157          danz = 0
158          WRITE (id_string,'(I4.4)')  danz
159          INQUIRE ( &
160               FILE='PLOT2D_'//modus//TRIM( model_string )//'_'//id_string, &
161               EXIST=found )
162!
163!--       Find out the number of files (equal to the number of PEs which
164!--       have been used in PALM) and open them
165          DO  WHILE ( found )
[1]166
[108]167             OPEN ( danz+110, &
168                  FILE='PLOT2D_'//modus//TRIM( model_string )//'_'//id_string, &
169                  FORM='UNFORMATTED' )
170             danz = danz + 1
171             WRITE (id_string,'(I4.4)')  danz
172             INQUIRE ( &
173                  FILE='PLOT2D_'//modus//TRIM( model_string )//'_'//id_string, &
174                  EXIST=found )
[1]175
[108]176          ENDDO
[1]177
178!
[108]179!--       Inquire whether an iso2d parameter file exists
180          INQUIRE( FILE='PLOT2D_'//modus//'_GLOBAL'//TRIM( model_string ), &
181               EXIST=iso2d_output )
[1]182
183!
[108]184!--       Inquire whether a NetCDF file exists
185          INQUIRE( FILE='DATA_2D_'//modus//'_NETCDF'//TRIM( model_string ), &
186               EXIST=netcdf_0 )
[1]187
188!
[108]189!--       Inquire whether a NetCDF file for time-averaged data exists
190          INQUIRE( FILE='DATA_2D_'//modus//'_AV_NETCDF'//TRIM( model_string ),&
191               EXIST=netcdf_1 )
[1]192
[108]193          IF ( netcdf_0  .OR.  netcdf_1 )  THEN
194             netcdf_output = .TRUE.
[493]195!
196!--          Inquire whether the NetCDF file is already complete (parallel
197!--          output)
198             INQUIRE( FILE='NO_COMBINE_PLOT_FIELDS_'//modus, &
199                      EXIST=netcdf_parallel )
200             IF ( netcdf_parallel )  THEN
201                netcdf_parallel = .TRUE.
202             ELSE
203                netcdf_parallel = .FALSE.
204             ENDIF
[108]205          ELSE
206             netcdf_output = .FALSE.
207          ENDIF
[1]208
209!
[108]210!--       Info-output
211          PRINT*, ''
212          PRINT*, '*** combine_plot_fields ***'
[1]213#if defined( __netcdf )
[493]214          IF ( netcdf_output )  THEN
215             IF ( netcdf_parallel )  THEN
216             PRINT*, '    NetCDF ' // modus // '-data are in one file ', &
217                          '(NetCDF4-format) - merging not neccessary'
218             ELSE
219                PRINT*, '    NetCDF output enabled'
220             ENDIF
221          ENDIF
[1]222#else
[108]223          IF ( netcdf_output )  THEN
224             PRINT*, '--- Sorry, no NetCDF support on this host'
225             netcdf_output = .FALSE.
226          ENDIF
[1]227#endif
[493]228          IF ( .NOT. netcdf_parallel )  THEN
229             IF ( danz /= 0 )  THEN
230                PRINT*, '    ',modus,'-section:  ', danz, ' file(s) found'
231             ELSE
232                PRINT*, '    no ', modus, '-section data available'
233             ENDIF
[108]234          ENDIF
[1]235
[493]236          IF ( netcdf_output  .AND. .NOT. netcdf_parallel  .AND.  danz /= 0 ) &
237          THEN
[1]238#if defined( __netcdf )
[108]239             DO  av = 0, 1
[1]240
[108]241                IF ( av == 0  .AND.  .NOT.  netcdf_0 )  CYCLE
242                IF ( av == 1  .AND.  .NOT.  netcdf_1 )  CYCLE
[1]243
244!
[108]245!--             Open NetCDF dataset
246                IF ( av == 0 )  THEN
247                   filename = 'DATA_2D_'//modus//'_NETCDF' &
248                        //TRIM( model_string )
249                ELSE
250                   filename = 'DATA_2D_'//modus//'_AV_NETCDF' &
251                        //TRIM( model_string )
252                ENDIF
253                nc_stat = NF90_OPEN( filename, NF90_WRITE, id_set(av) )
254                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 1 )
[1]255
256!
[108]257!--             Get the list of variables (order of variables corresponds with
258!--             the order of data on the binary file)
259                var_list(av) = ' '    ! GET_ATT does not assign trailing blanks
260                nc_stat = NF90_GET_ATT( id_set(av), NF90_GLOBAL, 'VAR_LIST', &
261                     var_list(av) )
262                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 2 )
[1]263
264!
[108]265!--             Inquire id of the time coordinate variable
266                nc_stat = NF90_INQ_VARID( id_set(av), 'time', id_var_time(av) )
267                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 3 )
[1]268
269!
[108]270!--             Count number of variables; there is one more semicolon in the
271!--             string than variable names
272                num_var(av) = -1
273                DO  i = 1, LEN( var_list(av) )
274                   IF ( var_list(av)(i:i) == ';' )  num_var(av) = num_var(av) +1
275                ENDDO
[1]276
277!
[108]278!--             Extract the variable names from the list and inquire their
279!--             NetCDF IDs
280                pos = INDEX( var_list(av), ';' )
[1]281!
[108]282!--             Loop over all variables
283                DO  i = 1, num_var(av)
[1]284
285!
[108]286!--                Extract variable name from list
287                   var_list(av) = var_list(av)(pos+1:)
288                   pos = INDEX( var_list(av), ';' )
289                   var_name = var_list(av)(1:pos-1)
[1]290
291!
[108]292!--                Get variable ID from name
293                   nc_stat = NF90_INQ_VARID( id_set(av), TRIM( var_name ), &
294                        id_var(av,i) )
295                   IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 4 )
[1]296
297!
[108]298!--                Get number of x/y/z levels for that variable
299                   nc_stat = NF90_INQUIRE_VARIABLE( id_set(av), id_var(av,i), &
300                        dimids = id_dims_loc )
301                   IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 5 )
302                   id_dims(av,:) = id_dims_loc
[1]303
304!
[108]305!--                Inquire dimension ID
306                   DO  j = 1, 4
307                      nc_stat = NF90_INQUIRE_DIMENSION( id_set(av), &
308                           id_dims(av,j), dimname, levels(av,i) )
309                      IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 6 )
[1]310
[108]311                      IF ( modus == 'XY' .AND. INDEX(dimname, 'z') /= 0 )  EXIT
312                      IF ( modus == 'XZ' .AND. INDEX(dimname, 'y') /= 0 )  EXIT
313                      IF ( modus == 'YZ' .AND. INDEX(dimname, 'x') /= 0 )  EXIT
314                   ENDDO
315
[1]316                ENDDO
317
[108]318             ENDDO   ! av = 0, 1
[1]319
[179]320#endif
[108]321          ENDIF
[1]322
323!
[108]324!--       Read the arrays, as long as the end of the file is reached
[493]325          IF ( .NOT. netcdf_parallel )  THEN
[1]326
[493]327             fanz          =         0
328             current_level =         1
329             current_var   = 999999999
[1]330
[493]331             DO  WHILE ( danz /= 0 )
332
[1]333!
[493]334!--             Loop over all files (reading data of the subdomains)
335                DO  id = 0, danz-1
[1]336!
[493]337!--                File from PE0 contains special information at the beginning,
338!--                concerning the lower and upper indices of the total-domain
339!--                used in PALM (nxag, nxeg, nyag, nyeg) and the lower and
340!--                upper indices of the array to be writte by this routine
341!--                (nxa, nxe, nya, nye). Usually in the horizontal directions
342!--                nxag=-1 and nxa=0 while all other variables have the same
343!--                value (i.e. nxeg=nxe).
344!--                Allocate necessary arrays, open the output file and write
345!--                the coordinate informations needed by ISO2D.
346                   IF ( id == 0  .AND.  fanz(0) == 0  .AND.  fanz(1) == 0 ) THEN
347                      READ ( id+110 )  nxag, nxeg, nyag, nyeg
348                      READ ( id+110 )  nxa, nxe, nya, nye
349                      ALLOCATE ( eta(nya:nye), ho(nxa:nxe), hu(nxa:nxe), &
350                                 pf(nxag:nxeg,nyag:nyeg) )
351                      READ ( id+110 )  dx, eta, hu, ho
[1]352
[493]353                      IF ( iso2d_output )  THEN
354                         OPEN ( 2, FILE='PLOT2D_'//modus//TRIM( model_string ),&
355                                   FORM='UNFORMATTED' )
356                         WRITE ( 2 )  dx, eta, hu, ho
357                      ENDIF
[108]358                   ENDIF
[1]359!
[493]360!--                Read output time
361                   IF ( netcdf_output  .AND.  id == 0 )  THEN
362                      IF ( netcdf_1 )  THEN
363                         READ ( id+110, END=998 )  simulated_time, time_step, av
364                      ELSE
[1]365!
[493]366!--                      For compatibility with earlier PALM versions
367                         READ ( id+110, END=998 )  simulated_time, time_step
368                         av = 0
369                      ENDIF
[108]370                   ENDIF
[1]371!
[493]372!--                Read subdomain indices
373                   READ ( id+110, END=998 )  xa, xe, ya, ye
[1]374!
[493]375!--                IF the PE made no output (in case that no part of the
376!--                cross-section is situated on this PE), indices have the
377!--                value -1
378                   IF ( .NOT. ( xa == -1  .AND.  xe == -1  .AND. &
379                                ya == -1  .AND.  ye == -1 ) )  THEN
[1]380!
[493]381!--                   Read the subdomain grid-point values
382                      ALLOCATE( pf_tmp(xa:xe,ya:ye) )
383                      READ ( id+110 )  pf_tmp
384                      pf(xa:xe,ya:ye) = pf_tmp
385                      DEALLOCATE( pf_tmp )
386                   ENDIF
387                   IF ( id == 0 )  fanz(av) = fanz(av) + 1
[1]388
[493]389                ENDDO
[1]390!
[493]391!--             Write the data of the total domain cross-section
392                IF ( iso2d_output )  WRITE ( 2 )  pf(nxa:nxe,nya:nye)
[1]393       
394!
[493]395!--             Write same data in NetCDF format
396                IF ( netcdf_output )  THEN
[1]397#if defined( __netcdf )
398!
[493]399!--                Check if a new time step has begun; if yes write data to
400!--                time axis
401                   IF ( current_var(av) > num_var(av) )  THEN
402                      current_var(av) = 1
403                      nc_stat = NF90_PUT_VAR( id_set(av), id_var_time(av), &
404                                              (/ simulated_time /),        &
405                                              start = (/ time_step /),     &
406                                              count = (/ 1 /) )
407                      IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 7 )
408                   ENDIF
[1]409
410!
[493]411!--                Now write the data; this is mode dependent
412                   SELECT CASE ( modus )
[1]413
[493]414                      CASE ( 'XY' )
415                         nc_stat = NF90_PUT_VAR( id_set(av),                   &
[1]416                                           id_var(av,current_var(av)),         &
417                                           pf(nxa:nxe,nya:nye),                &
418                             start = (/ 1, 1, current_level(av), time_step /), &
419                                      count = (/ nxe-nxa+1, nye-nya+1, 1, 1 /) )
[493]420                         IF ( nc_stat /= NF90_NOERR )  THEN
421                            CALL handle_netcdf_error( 8 )
422                         ENDIF
[1]423                 
[493]424                      CASE ( 'XZ' )
425                         nc_stat = NF90_PUT_VAR( id_set(av),                   &
[1]426                                           id_var(av,current_var(av)),         &
427                                           pf(nxa:nxe,nya:nye),                &
428                             start = (/ 1, current_level(av), 1, time_step /), &
429                                      count = (/ nxe-nxa+1, 1, nye-nya+1, 1 /) )
[493]430                         IF ( nc_stat /= NF90_NOERR )  THEN
431                            CALL handle_netcdf_error( 9 )
432                         ENDIF
[1]433
[493]434                      CASE ( 'YZ' )
435                         nc_stat = NF90_PUT_VAR( id_set(av),                   &
[1]436                                           id_var(av,current_var(av)),         &
437                                           pf(nxa:nxe,nya:nye),                &
438                             start = (/ current_level(av), 1, 1, time_step /), &
439                                      count = (/ 1, nxe-nxa+1, nye-nya+1, 1 /) )
[493]440                         IF ( nc_stat /= NF90_NOERR )  THEN
441                            CALL handle_netcdf_error( 10 )
442                         ENDIF
[1]443
[493]444                   END SELECT
[1]445
446!
[493]447!--                Data is written, check if max level is reached
448                   current_level(av) = current_level(av) + 1
449                   IF ( current_level(av) > levels(av,current_var(av)) )  THEN
450                      current_level(av) = 1
451                      current_var(av)   = current_var(av) + 1
452                   ENDIF
[108]453
[179]454#endif
[493]455                ENDIF
[1]456
[493]457             ENDDO
[1]458
[493]459          ENDIF
460
461998       IF ( danz /= 0  .AND.  .NOT. netcdf_parallel )  THEN
[1]462!
[108]463!--          Print the number of the arrays processed
464             WRITE (*,'(16X,I4,A)')  fanz(0)+fanz(1), ' array(s) processed'
465             IF ( fanz(1) /= 0 )  THEN
466                WRITE (*,'(16X,I4,A)')  fanz(1), ' array(s) are time-averaged'
467             ENDIF
[1]468
469!
[108]470!--          Close all files and deallocate arrays
471             DO  id = 0, danz-1
472                CLOSE ( id+110 )
473             ENDDO
474             CLOSE ( 2 )
475             DEALLOCATE ( eta, ho, hu, pf )
[1]476
477!
[493]478!--          Close the NetCDF file
479             IF ( netcdf_output )  THEN
[1]480#if defined( __netcdf )
[493]481                IF ( netcdf_0 )  THEN
482                   nc_stat = NF90_CLOSE( id_set(0) )
483                   IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 11 )
484                ENDIF
485                IF ( netcdf_1 )  THEN
486                   nc_stat = NF90_CLOSE( id_set(1) )
487                   IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 12 )
488                ENDIF
489#endif
[108]490             ENDIF
[493]491          ENDIF
492
493!
494!--       Output required cpu time
495          IF ( danz /= 0  .AND.  .NOT. netcdf_parallel )  THEN
496#if defined( __lc ) || defined( __decalpha ) || defined( __nec )
497             CALL SYSTEM_CLOCK( count, count_rate )
498             cpu_end_time = REAL( count ) / REAL( count_rate )
499             WRITE (*,'(5X,A,F9.3,A)')  'Required cpu-time: ', &
500                                        cpu_end_time-cpu_start_time, ' sec'
501#elif defined( __ibm )
502             cpu_end_time = IRTC( ) * 1E-9
503             WRITE (*,'(5X,A,F9.3,A)')  'Required cpu-time: ', &
504                                        cpu_end_time-cpu_start_time, ' sec'
505#else
506             CONTINUE
[108]507#endif
[1]508          ENDIF
509
510!
[108]511!--       Choose the next cross-section
512          SELECT CASE ( modus )
513             CASE ( 'XY' )
514                modus = 'XZ'
515             CASE ( 'XZ' )
516                modus = 'YZ'
517             CASE ( 'YZ' )
518                modus = 'no'
519          END SELECT
[1]520
[108]521       ENDDO
[1]522
523
524!
[108]525!--    Combine the 3D-arrays
[493]526       netcdf_parallel = .FALSE.
[1]527
528!
[493]529!--    Info-output
530       PRINT*, ' '
531       PRINT*, '*** combine_plot_fields ***'
532
533!
534!--    Take current time
535#if defined( __lc ) || defined( __decalpha ) || defined( __nec )
536       CALL SYSTEM_CLOCK( count, count_rate )
537       cpu_start_time = REAL( count ) / REAL( count_rate )
538#elif defined( __ibm )
539       cpu_start_time = IRTC( ) * 1E-9
540#else
541       PRINT*,  '+++ INFORMATIVE: no time measurement defined on this host'
542#endif
543
544!
[108]545!--    Inquire whether an avs fld file exists
546       INQUIRE( FILE='PLOT3D_FLD'//TRIM( model_string ), EXIST=avs_output )
[1]547
548!
[108]549!--    Inquire whether a NetCDF file exists
550       INQUIRE( FILE='DATA_3D_NETCDF'//TRIM( model_string ), EXIST=netcdf_0 )
[1]551
552!
[108]553!--    Inquire whether a NetCDF file for time-averaged data exists
554       INQUIRE( FILE='DATA_3D_AV_NETCDF'//TRIM( model_string ), EXIST=netcdf_1 )
[1]555
[108]556       IF ( netcdf_0  .OR.  netcdf_1 )  THEN
557          netcdf_output = .TRUE.
[493]558!
559!--       Inquire whether the NetCDF file is already complete (parallel output)
560          INQUIRE( FILE='NO_COMBINE_PLOT_FIELDS_3D', EXIST=netcdf_parallel )
561          IF ( netcdf_parallel )  THEN
562             netcdf_parallel = .TRUE.
563          ELSE
564             netcdf_parallel = .FALSE.
565          ENDIF
[108]566       ELSE
567          netcdf_output = .FALSE.
568       ENDIF
[1]569
570!
[493]571!--    Check, if file from PE0 exists; not neccessary in case of parallel
572!--    PALM output
573       IF ( .NOT. netcdf_parallel )  THEN
574          danz = 0
575          WRITE (id_string,'(I4.4)')  danz
576          INQUIRE ( &
577               FILE='PLOT3D_DATA'//TRIM( model_string )//'_'//TRIM( id_string ),  &
578               EXIST=found )
579       ELSE
580          found = .FALSE.
581       ENDIF
[1]582
583!
[108]584!--    Combination only works, if data are not compressed. In that case,
585!--    PALM created a flag file (PLOT3D_COMPRESSED)
586       INQUIRE ( FILE='PLOT3D_COMPRESSED'//TRIM( model_string ), &
587            EXIST=compressed )
[1]588
589!
[108]590!--    Find out the number of files and open them
591       DO  WHILE ( found  .AND.  .NOT. compressed )
[1]592
[108]593          OPEN ( danz+110, &
594               FILE='PLOT3D_DATA'//TRIM( model_string )//'_'//TRIM(id_string), &
595               FORM='UNFORMATTED')
596          danz = danz + 1
597          WRITE (id_string,'(I4.4)')  danz
598          INQUIRE ( &
599               FILE='PLOT3D_DATA'//TRIM( model_string )//'_'//TRIM(id_string), &
600               EXIST=found )
[1]601
[108]602       ENDDO
[1]603
604#if defined( __netcdf )
[493]605       IF ( netcdf_output )  THEN
606          IF ( netcdf_parallel )  THEN
607             PRINT*, '    NetCDF data are in one file (NetCDF4-format)', &
608                          ' - merging not neccessary'
609          ELSE
610             PRINT*, '    NetCDF output enabled'
611          ENDIF
612       ENDIF
[1]613#else
614       IF ( netcdf_output )  THEN
615          PRINT*, '--- Sorry, no NetCDF support on this host'
616          netcdf_output = .FALSE.
617       ENDIF
618#endif
[493]619       IF ( .NOT. netcdf_parallel )  THEN
620          IF ( danz /= 0 )  THEN
621             PRINT*, '    3D-data:     ', danz, ' file(s) found'
[108]622          ELSE
[493]623             IF ( found .AND. compressed )  THEN
624                PRINT*, '+++ no 3D-data processing, since data are compressed'
625             ELSE
626                PRINT*, '    no 3D-data file available'
627             ENDIF
[108]628          ENDIF
[1]629       ENDIF
630
[493]631       IF ( netcdf_output  .AND. .NOT. netcdf_parallel  .AND.  danz /= 0 )  THEN
[1]632#if defined( __netcdf )
[108]633          DO  av = 0, 1
[1]634
[108]635             IF ( av == 0  .AND.  .NOT.  netcdf_0 )  CYCLE
636             IF ( av == 1  .AND.  .NOT.  netcdf_1 )  CYCLE
[1]637
638!
[108]639!--          Open NetCDF dataset
640             IF ( av == 0 )  THEN
641                filename = 'DATA_3D_NETCDF'//TRIM( model_string )
642             ELSE
643                filename = 'DATA_3D_AV_NETCDF'//TRIM( model_string )
644             ENDIF
645             nc_stat = NF90_OPEN( filename, NF90_WRITE, id_set(av) )
646             IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 13 )
[1]647
648
649!
[108]650!--          Get the list of variables (order of variables corresponds with the
651!--          order of data on the binary file)
652             var_list(av) = ' '    ! GET_ATT does not assign trailing blanks
653             nc_stat = NF90_GET_ATT( id_set(av), NF90_GLOBAL, 'VAR_LIST', &
654                  var_list(av) )
655             IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 14 )
[1]656
657!
[108]658!--          Inquire id of the time coordinate variable
659             nc_stat = NF90_INQ_VARID( id_set(av), 'time', id_var_time(av) )
660             IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 15 )
[1]661
662!
[108]663!--          Count number of variables; there is one more semicolon in the
664!--          string than variable names
665             num_var(av) = -1
666             DO  i = 1, LEN( var_list(av) )
667                IF ( var_list(av)(i:i) == ';' )  num_var(av) = num_var(av) + 1
668             ENDDO
[1]669
670!
[108]671!--          Extract the variable names from the list and inquire their NetCDF
672!--          IDs
673             pos = INDEX( var_list(av), ';' )
[1]674!
[108]675!--          Loop over all variables
676             DO  i = 1, num_var(av)
[1]677
678!
[108]679!--             Extract variable name from list
680                var_list(av) = var_list(av)(pos+1:)
681                pos = INDEX( var_list(av), ';' )
682                var_name = var_list(av)(1:pos-1)
[1]683
684!
[108]685!--             Get variable ID from name
686                nc_stat = NF90_INQ_VARID( id_set(av), TRIM( var_name ), &
687                     id_var(av,i) )
688                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 16 )
[1]689
[108]690             ENDDO
[1]691
[108]692          ENDDO    ! av=0,1
[1]693
[179]694#endif
[108]695       ENDIF
[1]696
697!
[108]698!--    Read arrays, until the end of the file is reached
[493]699       IF ( .NOT. netcdf_parallel )  THEN
[1]700
[493]701          current_var = 999999999
702          fanz = 0
703          DO  WHILE ( danz /= 0 )
704
[1]705!
[493]706!--          Loop over all files
707             DO  id = 0, danz-1
[1]708!
[493]709!--             File from PE0 contains special information at the beginning,
710!--             concerning the lower and upper indices of the total-domain used
711!--             in PALM (nxag, nxeg, nyag, nyeg, nzag, nzeg) and the lower and
712!--             upper indices of the array to be written by this routine (nxa,
713!--             nxe, nya, nye, nza, nze). Usually nxag=-1 and nxa=0, nyag=-1
714!--             and nya=0, nzeg=nz and nze=nz_plot3d.
715!--             Allocate necessary array and open the output file.
716                IF ( id == 0  .AND.  fanz(0) == 0  .AND.  fanz(1) == 0 )  THEN
717                   READ ( id+110 )  nxag, nxeg, nyag, nyeg, nzag, nzeg
718                   READ ( id+110 )  nxa, nxe, nya, nye, nza, nze
719                   ALLOCATE ( pf3d(nxa:nxe,nya:nye,nza:nze) )
720                   IF ( avs_output )  THEN
721                      OPEN ( 2, FILE='PLOT3D_DATA'//TRIM( model_string ), &
722                             FORM='UNFORMATTED' )
723                   ENDIF
[108]724                ENDIF
[1]725
726!
[493]727!--             Read output time
728                IF ( netcdf_output  .AND.  id == 0 )  THEN
729                   IF ( netcdf_1 )  THEN
730                      READ ( id+110, END=999 )  simulated_time, time_step, av
731                   ELSE
[1]732!
[493]733!--                   For compatibility with earlier PALM versions
734                      READ ( id+110, END=999 )  simulated_time, time_step
735                      av = 0
736                   ENDIF
[108]737                ENDIF
[1]738
739!
[493]740!--             Read subdomain indices and grid point values
741                READ ( id+110, END=999 )  xa, xe, ya, ye, za, ze
742                ALLOCATE( pf3d_tmp(xa:xe,ya:ye,za:ze) )
743                READ ( id+110 )  pf3d_tmp
[210]744
[493]745                xxa = MAX( nxa, xa )
746                xxe = MIN( nxe, xe )
747                yya = MAX( nya, ya )
748                yye = MIN( nye, ye )
749                zza = MAX( nza, za )
750                zze = MIN( nze, ze )
751                DO  k = zza, zze
752                   DO  j = yya, yye
753                      DO  i = xxa, xxe
754                         pf3d(i,j,k) = pf3d_tmp(i,j,k)
755                      ENDDO
[210]756                   ENDDO
757                ENDDO
758
[493]759                DEALLOCATE( pf3d_tmp )
760                IF ( id == 0 )  fanz(av) = fanz(av) + 1
[1]761
[493]762             ENDDO
[1]763
764!
[493]765!--          Write data of the total domain
766             IF ( avs_output )  WRITE ( 2 )  pf3d(nxa:nxe,nya:nye,nza:nze)
[1]767       
768!
[493]769!--          Write same data in NetCDF format
770             IF ( netcdf_output )  THEN
[1]771#if defined( __netcdf )
772!
[493]773!--             Check if a new time step has begun; if yes write data to time
774!--             axis
775                IF ( current_var(av) > num_var(av) )  THEN
776                   current_var(av) = 1
777                   nc_stat = NF90_PUT_VAR( id_set(av), id_var_time(av), &
778                                      (/ simulated_time /),&
779                                      start = (/ time_step /), count = (/ 1 /) )
780                   IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 17 )
781                ENDIF
[1]782
783!
[493]784!--             Now write the data
785                nc_stat = NF90_PUT_VAR( id_set(av), id_var(av,current_var(av)),&
786                                        pf3d, start = (/ 1, 1, 1, time_step /),&
[1]787                              count = (/ nxe-nxa+1, nye-nya+1, nze-nza+1, 1 /) )
[493]788                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 18 )
[1]789
[493]790                current_var(av) = current_var(av) + 1
[1]791
792#endif
[493]793             ENDIF
[1]794
[493]795          ENDDO
[1]796
[493]797       ENDIF
798
799999    IF ( danz /= 0  .AND.  .NOT. netcdf_parallel )  THEN
[1]800!
[108]801!--       Print the number of arrays processed
802          WRITE (*,'(16X,I4,A)')  fanz(0)+fanz(1), ' array(s) processed'
803          IF ( fanz(1) /= 0 )  THEN
804             WRITE (*,'(16X,I4,A)')  fanz(1), ' array(s) are time-averaged'
805          ENDIF
[1]806!
[108]807!--       Close all files and deallocate array
808          DO  id = 0, danz-1
809             CLOSE ( id+110 )
810          ENDDO
811          CLOSE ( 2 )
812          DEALLOCATE ( pf3d )
[1]813!
[108]814!--       Close the NetCDF file
815          IF ( netcdf_output )  THEN
[1]816#if defined( __netcdf )
[108]817             IF ( netcdf_0 )  THEN
818                nc_stat = NF90_CLOSE( id_set(0) )
819                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 19 )
820             ENDIF
821             IF ( netcdf_1 )  THEN
822                nc_stat = NF90_CLOSE( id_set(1) )
823                IF ( nc_stat /= NF90_NOERR )  CALL handle_netcdf_error( 20 )
824             ENDIF
825#endif
[1]826          ENDIF
[493]827
828!
829!--       Output required cpu time
830#if defined( __lc ) || defined( __decalpha ) || defined( __nec )
831          CALL SYSTEM_CLOCK( count, count_rate )
832          cpu_end_time = REAL( count ) / REAL( count_rate )
833          WRITE (*,'(5X,A,F9.3,A)')  'Required cpu-time: ', &
834                                     cpu_end_time-cpu_start_time, ' sec'
835#elif defined( __ibm )
836          cpu_end_time = IRTC( ) * 1E-9
837          WRITE (*,'(5X,A,F9.3,A)')  'Required cpu-time: ', &
838                                     cpu_end_time-cpu_start_time, ' sec'
839#endif
840
[1]841       ENDIF
842
[108]843    ENDDO  ! models
[1]844
[108]845
[1]846 CONTAINS
847
848
849    SUBROUTINE handle_netcdf_error( errno )
850!
851!--    Prints out a text message corresponding to the current NetCDF status
852
853       IMPLICIT NONE
854
855       INTEGER, INTENT(IN) ::  errno
856
[179]857#if defined( __netcdf )
[1]858       IF ( nc_stat /= NF90_NOERR )  THEN
859          PRINT*, '+++ combine_plot_fields  netcdf: ', av, errno, &
860                  TRIM( nf90_strerror( nc_stat ) )
861       ENDIF
[179]862#endif
[1]863
864    END SUBROUTINE handle_netcdf_error
865
866
867 END PROGRAM combine_plot_fields
868
869
870
Note: See TracBrowser for help on using the repository browser.