PROGRAM combine_plot_fields !------------------------------------------------------------------------------! ! Current revisions: ! ----------------- ! ! Former revisions: ! ----------------- ! $Id: combine_plot_fields.f90 692 2011-03-04 08:52:15Z letzel $ ! ! 691 2011-03-04 08:45:30Z maronga ! Bugfix: check for precursor ocean runs, removed typo ! ! 493 2010-03-01 08:30:24Z raasch ! Exit in case of already complete NetCDF data (due to parallel output in PALM) ! cpu measurements included ! ! 210 2008-11-06 08:54:02Z raasch ! Size of pf3d adjusted to the required output size (1 gridpoint less, along ! all three dimensions), because output of a subset of the data ! (pf3d(nxa:nxe...) in the NF90_PUT_VAR statement caused segmentation fault ! with the INTEL compiler. ! Subdomain data are read into temporary arrays pf_tmp/pf3d_tmp in order to ! avoid INTEL compiler warnings about (automatic) creation of temporary arrays ! Bugfix: three misplaced #endif directives ! ! 114 2007-10-10 00:03:15Z raasch ! Bugfix: model_string needed a default value ! ! Aug 07 Loop for processing of output by coupled runs, id_string does not ! contain modus any longer ! ! 18/01/06 Output of time-averaged data ! ! 25/05/05 Errors removed ! ! 26/04/05 Output in NetCDF format, iso2d and avs output only if parameter ! file exists ! ! 31/10/01 All comments and messages translated into English ! ! 23/02/99 Keine Bearbeitung komprimierter 3D-Daten ! Ursprungsversion vom 28/07/97 ! ! ! Description: ! ------------ ! This routine combines data of the PALM-subdomains into one file. In PALM ! every processor element opens its own file and writes 2D- or 3D-binary-data ! into it (different files are opened for xy-, xz-, yz-cross-sections and ! 3D-data). For plotting or analyzing these PE-data have to be collected and ! to be put into single files, which is done by this routine. ! Output format is NetCDF. Additionally, a data are output in a binary format ! readable by ISO2D-software (cross-sections) and by AVS (3D-data). !------------------------------------------------------------------------------! #if defined( __netcdf ) USE netcdf #endif IMPLICIT NONE ! !-- Local variables CHARACTER (LEN=2) :: modus, model_string CHARACTER (LEN=4) :: id_string CHARACTER (LEN=10) :: dimname, var_name CHARACTER (LEN=40) :: filename CHARACTER (LEN=2000), DIMENSION(0:1) :: var_list INTEGER, PARAMETER :: spk = SELECTED_REAL_KIND( 6 ) INTEGER :: av, danz, i, id, & j, k, model, models, nc_stat, nxa, nxag, nxe, nxeg, nya, & nyag, nye, nyeg, nza, nzag, nze, nzeg, pos, time_step, xa, xe, & xxa, xxe, ya, ye, yya, yye, za, ze, zza, zze #if defined( __lc ) || defined( __decalpha ) INTEGER(8) :: count, count_rate #elif defined( __nec ) INTEGER :: count, count_rate #elif defined( __ibm ) INTEGER(8) :: IRTC #endif INTEGER, DIMENSION(0:1) :: current_level, current_var, fanz, id_set, & id_var_time, num_var INTEGER, DIMENSION(4) :: id_dims_loc INTEGER, DIMENSION(0:1,4) :: id_dims INTEGER, DIMENSION(0:1,1000) :: id_var, levels LOGICAL :: avs_output, compressed, found, iso2d_output, netcdf_output, & netcdf_parallel, netcdf_0, netcdf_1 REAL :: cpu_start_time, cpu_end_time, dx, simulated_time REAL, DIMENSION(:), ALLOCATABLE :: eta, ho, hu REAL, DIMENSION(:,:), ALLOCATABLE :: pf, pf_tmp REAL(spk), DIMENSION(:,:,:), ALLOCATABLE :: pf3d, pf3d_tmp PRINT*, '' PRINT*, '' PRINT*, '*** combine_plot_fields ***' ! !-- Find out if a coupled run has been carried out INQUIRE( FILE='COUPLING_PORT_OPENED', EXIST=found ) IF ( found ) THEN models = 2 PRINT*, ' coupled run' ELSE models = 1 PRINT*, ' uncoupled run' ENDIF ! !-- Find out if a precursor ocean run has been carried out INQUIRE( FILE='PRECURSOR_OCEAN', EXIST=found ) IF ( found ) THEN model_string = '_O' PRINT*, ' precursor ocean run' ELSE model_string = '' ENDIF ! !-- Do everything for each model DO model = 1, models ! !-- Set the model string used to identify the filenames IF ( models == 2 ) THEN PRINT*, '' PRINT*, '*** combine_plot_fields ***' IF ( model == 2 ) THEN model_string = '_O' PRINT*, ' now combining ocean data' PRINT*, ' ========================' ELSE PRINT*, ' now combining atmosphere data' PRINT*, ' =============================' ENDIF ENDIF ! !-- 2D-arrays for ISO2D !-- Main loop for the three different cross-sections, starting with !-- xy-section modus = 'XY' PRINT*, '' DO WHILE ( modus == 'XY' .OR. modus == 'XZ' .OR. modus == 'YZ' ) ! !-- Take current time #if defined( __lc ) || defined( __decalpha ) || defined( __nec ) CALL SYSTEM_CLOCK( count, count_rate ) cpu_start_time = REAL( count ) / REAL( count_rate ) #elif defined( __ibm ) cpu_start_time = IRTC( ) * 1E-9 #else PRINT*, '+++ INFORMATIVE: no time measurement defined on this host' #endif netcdf_parallel = .FALSE. ! !-- Check, if file from PE0 exists. If it does not exist, PALM did not !-- create any output for this cross-section. danz = 0 WRITE (id_string,'(I4.4)') danz INQUIRE ( & FILE='PLOT2D_'//modus//TRIM( model_string )//'_'//id_string, & EXIST=found ) ! !-- Find out the number of files (equal to the number of PEs which !-- have been used in PALM) and open them DO WHILE ( found ) OPEN ( danz+110, & FILE='PLOT2D_'//modus//TRIM( model_string )//'_'//id_string, & FORM='UNFORMATTED' ) danz = danz + 1 WRITE (id_string,'(I4.4)') danz INQUIRE ( & FILE='PLOT2D_'//modus//TRIM( model_string )//'_'//id_string, & EXIST=found ) ENDDO ! !-- Inquire whether an iso2d parameter file exists INQUIRE( FILE='PLOT2D_'//modus//'_GLOBAL'//TRIM( model_string ), & EXIST=iso2d_output ) ! !-- Inquire whether a NetCDF file exists INQUIRE( FILE='DATA_2D_'//modus//'_NETCDF'//TRIM( model_string ), & EXIST=netcdf_0 ) ! !-- Inquire whether a NetCDF file for time-averaged data exists INQUIRE( FILE='DATA_2D_'//modus//'_AV_NETCDF'//TRIM( model_string ),& EXIST=netcdf_1 ) IF ( netcdf_0 .OR. netcdf_1 ) THEN netcdf_output = .TRUE. ! !-- Inquire whether the NetCDF file is already complete (parallel !-- output) INQUIRE( FILE='NO_COMBINE_PLOT_FIELDS_'//modus, & EXIST=netcdf_parallel ) IF ( netcdf_parallel ) THEN netcdf_parallel = .TRUE. ELSE netcdf_parallel = .FALSE. ENDIF ELSE netcdf_output = .FALSE. ENDIF ! !-- Info-output PRINT*, '' PRINT*, '*** combine_plot_fields ***' #if defined( __netcdf ) IF ( netcdf_output ) THEN IF ( netcdf_parallel ) THEN PRINT*, ' NetCDF ' // modus // '-data are in one file ', & '(NetCDF4-format) - merging not neccessary' ELSE PRINT*, ' NetCDF output enabled' ENDIF ENDIF #else IF ( netcdf_output ) THEN PRINT*, '--- Sorry, no NetCDF support on this host' netcdf_output = .FALSE. ENDIF #endif IF ( .NOT. netcdf_parallel ) THEN IF ( danz /= 0 ) THEN PRINT*, ' ',modus,'-section: ', danz, ' file(s) found' ELSE PRINT*, ' no ', modus, '-section data available' ENDIF ENDIF IF ( netcdf_output .AND. .NOT. netcdf_parallel .AND. danz /= 0 ) & THEN #if defined( __netcdf ) DO av = 0, 1 IF ( av == 0 .AND. .NOT. netcdf_0 ) CYCLE IF ( av == 1 .AND. .NOT. netcdf_1 ) CYCLE ! !-- Open NetCDF dataset IF ( av == 0 ) THEN filename = 'DATA_2D_'//modus//'_NETCDF' & //TRIM( model_string ) ELSE filename = 'DATA_2D_'//modus//'_AV_NETCDF' & //TRIM( model_string ) ENDIF nc_stat = NF90_OPEN( filename, NF90_WRITE, id_set(av) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 1 ) ! !-- Get the list of variables (order of variables corresponds with !-- the order of data on the binary file) var_list(av) = ' ' ! GET_ATT does not assign trailing blanks nc_stat = NF90_GET_ATT( id_set(av), NF90_GLOBAL, 'VAR_LIST', & var_list(av) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 2 ) ! !-- Inquire id of the time coordinate variable nc_stat = NF90_INQ_VARID( id_set(av), 'time', id_var_time(av) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 3 ) ! !-- Count number of variables; there is one more semicolon in the !-- string than variable names num_var(av) = -1 DO i = 1, LEN( var_list(av) ) IF ( var_list(av)(i:i) == ';' ) num_var(av) = num_var(av) +1 ENDDO ! !-- Extract the variable names from the list and inquire their !-- NetCDF IDs pos = INDEX( var_list(av), ';' ) ! !-- Loop over all variables DO i = 1, num_var(av) ! !-- Extract variable name from list var_list(av) = var_list(av)(pos+1:) pos = INDEX( var_list(av), ';' ) var_name = var_list(av)(1:pos-1) ! !-- Get variable ID from name nc_stat = NF90_INQ_VARID( id_set(av), TRIM( var_name ), & id_var(av,i) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 4 ) ! !-- Get number of x/y/z levels for that variable nc_stat = NF90_INQUIRE_VARIABLE( id_set(av), id_var(av,i), & dimids = id_dims_loc ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 5 ) id_dims(av,:) = id_dims_loc ! !-- Inquire dimension ID DO j = 1, 4 nc_stat = NF90_INQUIRE_DIMENSION( id_set(av), & id_dims(av,j), dimname, levels(av,i) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 6 ) IF ( modus == 'XY' .AND. INDEX(dimname, 'z') /= 0 ) EXIT IF ( modus == 'XZ' .AND. INDEX(dimname, 'y') /= 0 ) EXIT IF ( modus == 'YZ' .AND. INDEX(dimname, 'x') /= 0 ) EXIT ENDDO ENDDO ENDDO ! av = 0, 1 #endif ENDIF ! !-- Read the arrays, as long as the end of the file is reached IF ( .NOT. netcdf_parallel ) THEN fanz = 0 current_level = 1 current_var = 999999999 DO WHILE ( danz /= 0 ) ! !-- Loop over all files (reading data of the subdomains) DO id = 0, danz-1 ! !-- File from PE0 contains special information at the beginning, !-- concerning the lower and upper indices of the total-domain !-- used in PALM (nxag, nxeg, nyag, nyeg) and the lower and !-- upper indices of the array to be writte by this routine !-- (nxa, nxe, nya, nye). Usually in the horizontal directions !-- nxag=-1 and nxa=0 while all other variables have the same !-- value (i.e. nxeg=nxe). !-- Allocate necessary arrays, open the output file and write !-- the coordinate informations needed by ISO2D. IF ( id == 0 .AND. fanz(0) == 0 .AND. fanz(1) == 0 ) THEN READ ( id+110 ) nxag, nxeg, nyag, nyeg READ ( id+110 ) nxa, nxe, nya, nye ALLOCATE ( eta(nya:nye), ho(nxa:nxe), hu(nxa:nxe), & pf(nxag:nxeg,nyag:nyeg) ) READ ( id+110 ) dx, eta, hu, ho IF ( iso2d_output ) THEN OPEN ( 2, FILE='PLOT2D_'//modus//TRIM( model_string ),& FORM='UNFORMATTED' ) WRITE ( 2 ) dx, eta, hu, ho ENDIF ENDIF ! !-- Read output time IF ( netcdf_output .AND. id == 0 ) THEN IF ( netcdf_1 ) THEN READ ( id+110, END=998 ) simulated_time, time_step, av ELSE ! !-- For compatibility with earlier PALM versions READ ( id+110, END=998 ) simulated_time, time_step av = 0 ENDIF ENDIF ! !-- Read subdomain indices READ ( id+110, END=998 ) xa, xe, ya, ye ! !-- IF the PE made no output (in case that no part of the !-- cross-section is situated on this PE), indices have the !-- value -1 IF ( .NOT. ( xa == -1 .AND. xe == -1 .AND. & ya == -1 .AND. ye == -1 ) ) THEN ! !-- Read the subdomain grid-point values ALLOCATE( pf_tmp(xa:xe,ya:ye) ) READ ( id+110 ) pf_tmp pf(xa:xe,ya:ye) = pf_tmp DEALLOCATE( pf_tmp ) ENDIF IF ( id == 0 ) fanz(av) = fanz(av) + 1 ENDDO ! !-- Write the data of the total domain cross-section IF ( iso2d_output ) WRITE ( 2 ) pf(nxa:nxe,nya:nye) ! !-- Write same data in NetCDF format IF ( netcdf_output ) THEN #if defined( __netcdf ) ! !-- Check if a new time step has begun; if yes write data to !-- time axis IF ( current_var(av) > num_var(av) ) THEN current_var(av) = 1 nc_stat = NF90_PUT_VAR( id_set(av), id_var_time(av), & (/ simulated_time /), & start = (/ time_step /), & count = (/ 1 /) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 7 ) ENDIF ! !-- Now write the data; this is mode dependent SELECT CASE ( modus ) CASE ( 'XY' ) nc_stat = NF90_PUT_VAR( id_set(av), & id_var(av,current_var(av)), & pf(nxa:nxe,nya:nye), & start = (/ 1, 1, current_level(av), time_step /), & count = (/ nxe-nxa+1, nye-nya+1, 1, 1 /) ) IF ( nc_stat /= NF90_NOERR ) THEN CALL handle_netcdf_error( 8 ) ENDIF CASE ( 'XZ' ) nc_stat = NF90_PUT_VAR( id_set(av), & id_var(av,current_var(av)), & pf(nxa:nxe,nya:nye), & start = (/ 1, current_level(av), 1, time_step /), & count = (/ nxe-nxa+1, 1, nye-nya+1, 1 /) ) IF ( nc_stat /= NF90_NOERR ) THEN CALL handle_netcdf_error( 9 ) ENDIF CASE ( 'YZ' ) nc_stat = NF90_PUT_VAR( id_set(av), & id_var(av,current_var(av)), & pf(nxa:nxe,nya:nye), & start = (/ current_level(av), 1, 1, time_step /), & count = (/ 1, nxe-nxa+1, nye-nya+1, 1 /) ) IF ( nc_stat /= NF90_NOERR ) THEN CALL handle_netcdf_error( 10 ) ENDIF END SELECT ! !-- Data is written, check if max level is reached current_level(av) = current_level(av) + 1 IF ( current_level(av) > levels(av,current_var(av)) ) THEN current_level(av) = 1 current_var(av) = current_var(av) + 1 ENDIF #endif ENDIF ENDDO ENDIF 998 IF ( danz /= 0 .AND. .NOT. netcdf_parallel ) THEN ! !-- Print the number of the arrays processed WRITE (*,'(16X,I4,A)') fanz(0)+fanz(1), ' array(s) processed' IF ( fanz(1) /= 0 ) THEN WRITE (*,'(16X,I4,A)') fanz(1), ' array(s) are time-averaged' ENDIF ! !-- Close all files and deallocate arrays DO id = 0, danz-1 CLOSE ( id+110 ) ENDDO CLOSE ( 2 ) DEALLOCATE ( eta, ho, hu, pf ) ! !-- Close the NetCDF file IF ( netcdf_output ) THEN #if defined( __netcdf ) IF ( netcdf_0 ) THEN nc_stat = NF90_CLOSE( id_set(0) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 11 ) ENDIF IF ( netcdf_1 ) THEN nc_stat = NF90_CLOSE( id_set(1) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 12 ) ENDIF #endif ENDIF ENDIF ! !-- Output required cpu time IF ( danz /= 0 .AND. .NOT. netcdf_parallel ) THEN #if defined( __lc ) || defined( __decalpha ) || defined( __nec ) CALL SYSTEM_CLOCK( count, count_rate ) cpu_end_time = REAL( count ) / REAL( count_rate ) WRITE (*,'(5X,A,F9.3,A)') 'Required cpu-time: ', & cpu_end_time-cpu_start_time, ' sec' #elif defined( __ibm ) cpu_end_time = IRTC( ) * 1E-9 WRITE (*,'(5X,A,F9.3,A)') 'Required cpu-time: ', & cpu_end_time-cpu_start_time, ' sec' #else CONTINUE #endif ENDIF ! !-- Choose the next cross-section SELECT CASE ( modus ) CASE ( 'XY' ) modus = 'XZ' CASE ( 'XZ' ) modus = 'YZ' CASE ( 'YZ' ) modus = 'no' END SELECT ENDDO ! !-- Combine the 3D-arrays netcdf_parallel = .FALSE. ! !-- Info-output PRINT*, ' ' PRINT*, '*** combine_plot_fields ***' ! !-- Take current time #if defined( __lc ) || defined( __decalpha ) || defined( __nec ) CALL SYSTEM_CLOCK( count, count_rate ) cpu_start_time = REAL( count ) / REAL( count_rate ) #elif defined( __ibm ) cpu_start_time = IRTC( ) * 1E-9 #else PRINT*, '+++ INFORMATIVE: no time measurement defined on this host' #endif ! !-- Inquire whether an avs fld file exists INQUIRE( FILE='PLOT3D_FLD'//TRIM( model_string ), EXIST=avs_output ) ! !-- Inquire whether a NetCDF file exists INQUIRE( FILE='DATA_3D_NETCDF'//TRIM( model_string ), EXIST=netcdf_0 ) ! !-- Inquire whether a NetCDF file for time-averaged data exists INQUIRE( FILE='DATA_3D_AV_NETCDF'//TRIM( model_string ), EXIST=netcdf_1 ) IF ( netcdf_0 .OR. netcdf_1 ) THEN netcdf_output = .TRUE. ! !-- Inquire whether the NetCDF file is already complete (parallel output) INQUIRE( FILE='NO_COMBINE_PLOT_FIELDS_3D', EXIST=netcdf_parallel ) IF ( netcdf_parallel ) THEN netcdf_parallel = .TRUE. ELSE netcdf_parallel = .FALSE. ENDIF ELSE netcdf_output = .FALSE. ENDIF ! !-- Check, if file from PE0 exists; not neccessary in case of parallel !-- PALM output IF ( .NOT. netcdf_parallel ) THEN danz = 0 WRITE (id_string,'(I4.4)') danz INQUIRE ( & FILE='PLOT3D_DATA'//TRIM( model_string )//'_'//TRIM( id_string ), & EXIST=found ) ELSE found = .FALSE. ENDIF ! !-- Combination only works, if data are not compressed. In that case, !-- PALM created a flag file (PLOT3D_COMPRESSED) INQUIRE ( FILE='PLOT3D_COMPRESSED'//TRIM( model_string ), & EXIST=compressed ) ! !-- Find out the number of files and open them DO WHILE ( found .AND. .NOT. compressed ) OPEN ( danz+110, & FILE='PLOT3D_DATA'//TRIM( model_string )//'_'//TRIM(id_string), & FORM='UNFORMATTED') danz = danz + 1 WRITE (id_string,'(I4.4)') danz INQUIRE ( & FILE='PLOT3D_DATA'//TRIM( model_string )//'_'//TRIM(id_string), & EXIST=found ) ENDDO #if defined( __netcdf ) IF ( netcdf_output ) THEN IF ( netcdf_parallel ) THEN PRINT*, ' NetCDF data are in one file (NetCDF4-format)', & ' - merging not neccessary' ELSE PRINT*, ' NetCDF output enabled' ENDIF ENDIF #else IF ( netcdf_output ) THEN PRINT*, '--- Sorry, no NetCDF support on this host' netcdf_output = .FALSE. ENDIF #endif IF ( .NOT. netcdf_parallel ) THEN IF ( danz /= 0 ) THEN PRINT*, ' 3D-data: ', danz, ' file(s) found' ELSE IF ( found .AND. compressed ) THEN PRINT*, '+++ no 3D-data processing, since data are compressed' ELSE PRINT*, ' no 3D-data file available' ENDIF ENDIF ENDIF IF ( netcdf_output .AND. .NOT. netcdf_parallel .AND. danz /= 0 ) THEN #if defined( __netcdf ) DO av = 0, 1 IF ( av == 0 .AND. .NOT. netcdf_0 ) CYCLE IF ( av == 1 .AND. .NOT. netcdf_1 ) CYCLE ! !-- Open NetCDF dataset IF ( av == 0 ) THEN filename = 'DATA_3D_NETCDF'//TRIM( model_string ) ELSE filename = 'DATA_3D_AV_NETCDF'//TRIM( model_string ) ENDIF nc_stat = NF90_OPEN( filename, NF90_WRITE, id_set(av) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 13 ) ! !-- Get the list of variables (order of variables corresponds with the !-- order of data on the binary file) var_list(av) = ' ' ! GET_ATT does not assign trailing blanks nc_stat = NF90_GET_ATT( id_set(av), NF90_GLOBAL, 'VAR_LIST', & var_list(av) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 14 ) ! !-- Inquire id of the time coordinate variable nc_stat = NF90_INQ_VARID( id_set(av), 'time', id_var_time(av) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 15 ) ! !-- Count number of variables; there is one more semicolon in the !-- string than variable names num_var(av) = -1 DO i = 1, LEN( var_list(av) ) IF ( var_list(av)(i:i) == ';' ) num_var(av) = num_var(av) + 1 ENDDO ! !-- Extract the variable names from the list and inquire their NetCDF !-- IDs pos = INDEX( var_list(av), ';' ) ! !-- Loop over all variables DO i = 1, num_var(av) ! !-- Extract variable name from list var_list(av) = var_list(av)(pos+1:) pos = INDEX( var_list(av), ';' ) var_name = var_list(av)(1:pos-1) ! !-- Get variable ID from name nc_stat = NF90_INQ_VARID( id_set(av), TRIM( var_name ), & id_var(av,i) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 16 ) ENDDO ENDDO ! av=0,1 #endif ENDIF ! !-- Read arrays, until the end of the file is reached IF ( .NOT. netcdf_parallel ) THEN current_var = 999999999 fanz = 0 DO WHILE ( danz /= 0 ) ! !-- Loop over all files DO id = 0, danz-1 ! !-- File from PE0 contains special information at the beginning, !-- concerning the lower and upper indices of the total-domain used !-- in PALM (nxag, nxeg, nyag, nyeg, nzag, nzeg) and the lower and !-- upper indices of the array to be written by this routine (nxa, !-- nxe, nya, nye, nza, nze). Usually nxag=-1 and nxa=0, nyag=-1 !-- and nya=0, nzeg=nz and nze=nz_plot3d. !-- Allocate necessary array and open the output file. IF ( id == 0 .AND. fanz(0) == 0 .AND. fanz(1) == 0 ) THEN READ ( id+110 ) nxag, nxeg, nyag, nyeg, nzag, nzeg READ ( id+110 ) nxa, nxe, nya, nye, nza, nze ALLOCATE ( pf3d(nxa:nxe,nya:nye,nza:nze) ) IF ( avs_output ) THEN OPEN ( 2, FILE='PLOT3D_DATA'//TRIM( model_string ), & FORM='UNFORMATTED' ) ENDIF ENDIF ! !-- Read output time IF ( netcdf_output .AND. id == 0 ) THEN IF ( netcdf_1 ) THEN READ ( id+110, END=999 ) simulated_time, time_step, av ELSE ! !-- For compatibility with earlier PALM versions READ ( id+110, END=999 ) simulated_time, time_step av = 0 ENDIF ENDIF ! !-- Read subdomain indices and grid point values READ ( id+110, END=999 ) xa, xe, ya, ye, za, ze ALLOCATE( pf3d_tmp(xa:xe,ya:ye,za:ze) ) READ ( id+110 ) pf3d_tmp xxa = MAX( nxa, xa ) xxe = MIN( nxe, xe ) yya = MAX( nya, ya ) yye = MIN( nye, ye ) zza = MAX( nza, za ) zze = MIN( nze, ze ) DO k = zza, zze DO j = yya, yye DO i = xxa, xxe pf3d(i,j,k) = pf3d_tmp(i,j,k) ENDDO ENDDO ENDDO DEALLOCATE( pf3d_tmp ) IF ( id == 0 ) fanz(av) = fanz(av) + 1 ENDDO ! !-- Write data of the total domain IF ( avs_output ) WRITE ( 2 ) pf3d(nxa:nxe,nya:nye,nza:nze) ! !-- Write same data in NetCDF format IF ( netcdf_output ) THEN #if defined( __netcdf ) ! !-- Check if a new time step has begun; if yes write data to time !-- axis IF ( current_var(av) > num_var(av) ) THEN current_var(av) = 1 nc_stat = NF90_PUT_VAR( id_set(av), id_var_time(av), & (/ simulated_time /),& start = (/ time_step /), count = (/ 1 /) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 17 ) ENDIF ! !-- Now write the data nc_stat = NF90_PUT_VAR( id_set(av), id_var(av,current_var(av)),& pf3d, start = (/ 1, 1, 1, time_step /),& count = (/ nxe-nxa+1, nye-nya+1, nze-nza+1, 1 /) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 18 ) current_var(av) = current_var(av) + 1 #endif ENDIF ENDDO ENDIF 999 IF ( danz /= 0 .AND. .NOT. netcdf_parallel ) THEN ! !-- Print the number of arrays processed WRITE (*,'(16X,I4,A)') fanz(0)+fanz(1), ' array(s) processed' IF ( fanz(1) /= 0 ) THEN WRITE (*,'(16X,I4,A)') fanz(1), ' array(s) are time-averaged' ENDIF ! !-- Close all files and deallocate array DO id = 0, danz-1 CLOSE ( id+110 ) ENDDO CLOSE ( 2 ) DEALLOCATE ( pf3d ) ! !-- Close the NetCDF file IF ( netcdf_output ) THEN #if defined( __netcdf ) IF ( netcdf_0 ) THEN nc_stat = NF90_CLOSE( id_set(0) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 19 ) ENDIF IF ( netcdf_1 ) THEN nc_stat = NF90_CLOSE( id_set(1) ) IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 20 ) ENDIF #endif ENDIF ! !-- Output required cpu time #if defined( __lc ) || defined( __decalpha ) || defined( __nec ) CALL SYSTEM_CLOCK( count, count_rate ) cpu_end_time = REAL( count ) / REAL( count_rate ) WRITE (*,'(5X,A,F9.3,A)') 'Required cpu-time: ', & cpu_end_time-cpu_start_time, ' sec' #elif defined( __ibm ) cpu_end_time = IRTC( ) * 1E-9 WRITE (*,'(5X,A,F9.3,A)') 'Required cpu-time: ', & cpu_end_time-cpu_start_time, ' sec' #endif ENDIF ENDDO ! models CONTAINS SUBROUTINE handle_netcdf_error( errno ) ! !-- Prints out a text message corresponding to the current NetCDF status IMPLICIT NONE INTEGER, INTENT(IN) :: errno #if defined( __netcdf ) IF ( nc_stat /= NF90_NOERR ) THEN PRINT*, '+++ combine_plot_fields netcdf: ', av, errno, & TRIM( nf90_strerror( nc_stat ) ) ENDIF #endif END SUBROUTINE handle_netcdf_error END PROGRAM combine_plot_fields