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

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

upper bounds of cross section and 3d output changed from nx+1,ny+1 to nx,ny; no output if redundant ghost layer data to NetCDF files

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