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

Last change on this file since 2689 was 2669, checked in by raasch, 7 years ago

file attributes and activation strings in .palm.iofiles revised, file extensions for nesting, masked output, wind turbine data, etc. changed

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