SUBROUTINE data_output_3d( av ) !--------------------------------------------------------------------------------! ! This file is part of PALM. ! ! PALM is free software: you can redistribute it and/or modify it under the terms ! of the GNU General Public License as published by the Free Software Foundation, ! either version 3 of the License, or (at your option) any later version. ! ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. ! ! You should have received a copy of the GNU General Public License along with ! PALM. If not, see . ! ! Copyright 1997-2012 Leibniz University Hannover !--------------------------------------------------------------------------------! ! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: data_output_3d.f90 1037 2012-10-22 14:10:22Z banzhafs $ ! ! 1036 2012-10-22 13:43:42Z raasch ! code put under GPL (PALM 3.9) ! ! 1031 2012-10-19 14:35:30Z raasch ! netCDF4 without parallel file support implemented ! ! 1007 2012-09-19 14:30:36Z franke ! Bugfix: missing calculation of ql_vp added ! ! 790 2011-11-29 03:11:20Z raasch ! bugfix: calculation of 'pr' must depend on the particle weighting factor, ! nzt+1 replaced by nz_do3d for 'pr' ! ! 771 2011-10-27 10:56:21Z heinze ! +lpt ! ! 759 2011-09-15 13:58:31Z raasch ! Splitting of parallel I/O ! ! 727 2011-04-20 20:05:25Z suehring ! Exchange ghost layers also for p_av. ! ! 725 2011-04-11 09:37:01Z suehring ! Exchange ghost layers for p regardless of used pressure solver (except SOR). ! ! 691 2011-03-04 08:45:30Z maronga ! Replaced simulated_time by time_since_reference_point ! ! 673 2011-01-18 16:19:48Z suehring ! When using Multigrid or SOR solver an additional CALL exchange_horiz is ! is needed for pressure output. ! ! 667 2010-12-23 12:06:00Z suehring/gryschka ! nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng in loops and ! allocation of arrays. Calls of exchange_horiz are modified. ! Skip-value skip_do_avs changed to a dynamic adaption of ghost points. ! ! 646 2010-12-15 13:03:52Z raasch ! bugfix: missing define statements for netcdf added ! ! 493 2010-03-01 08:30:24Z raasch ! netCDF4 support (parallel output) ! ! 355 2009-07-17 01:03:01Z letzel ! simulated_time in netCDF output replaced by time_since_reference_point. ! Output of netCDF messages with aid of message handling routine. ! Output of messages replaced by message handling routine. ! Bugfix: to_be_resorted => s_av for time-averaged scalars ! ! 96 2007-06-04 08:07:41Z raasch ! Output of density and salinity ! ! 75 2007-03-22 09:54:05Z raasch ! 2nd+3rd argument removed from exchange horiz ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.3 2006/06/02 15:18:59 raasch ! +argument "found", -argument grid in call of routine user_data_output_3d ! ! Revision 1.2 2006/02/23 10:23:07 raasch ! Former subroutine plot_3d renamed data_output_3d, pl.. renamed do.., ! .._anz renamed .._n, ! output extended to (almost) all quantities, output of user-defined quantities ! ! Revision 1.1 1997/09/03 06:29:36 raasch ! Initial revision ! ! ! Description: ! ------------ ! Output of the 3D-arrays in netCDF and/or AVS format. !------------------------------------------------------------------------------! USE array_kind USE arrays_3d USE averaging USE cloud_parameters USE control_parameters USE cpulog USE indices USE interfaces USE netcdf_control USE particle_attributes USE pegrid IMPLICIT NONE CHARACTER (LEN=9) :: simulated_time_mod INTEGER :: av, i, if, j, k, n, pos, prec, psi LOGICAL :: found, resorted REAL :: mean_r, s_r3, s_r4 REAL(spk), DIMENSION(:,:,:), ALLOCATABLE :: local_pf REAL, DIMENSION(:,:,:), POINTER :: to_be_resorted ! !-- Return, if nothing to output IF ( do3d_no(av) == 0 ) RETURN CALL cpu_log (log_point(14),'data_output_3d','start') ! !-- Open output file. !-- Also creates coordinate and fld-file for AVS. !-- For classic or 64bit netCDF output or output of other (old) data formats, !-- for a run on more than one PE, each PE opens its own file and !-- writes the data of its subdomain in binary format (regardless of the format !-- the user has requested). After the run, these files are combined to one !-- file by combine_plot_fields in the format requested by the user (netcdf !-- and/or avs). !-- For netCDF4/HDF5 output, data is written in parallel into one file. IF ( netcdf_output ) THEN IF ( netcdf_data_format < 5 ) THEN CALL check_open( 30 ) IF ( myid == 0 ) CALL check_open( 106+av*10 ) ELSE CALL check_open( 106+av*10 ) ENDIF ELSE IF ( avs_output .OR. ( numprocs > 1 ) ) CALL check_open( 30 ) ENDIF ! !-- Allocate a temporary array with the desired output dimensions. ALLOCATE( local_pf(nxlg:nxrg,nysg:nyng,nzb:nz_do3d) ) ! !-- Update the netCDF time axis #if defined( __netcdf ) IF ( myid == 0 .OR. netcdf_data_format > 4 ) THEN do3d_time_count(av) = do3d_time_count(av) + 1 IF ( netcdf_output ) THEN nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_time_3d(av), & (/ time_since_reference_point /), & start = (/ do3d_time_count(av) /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'data_output_3d', 376 ) ENDIF ENDIF #endif ! !-- Loop over all variables to be written. if = 1 DO WHILE ( do3d(av,if)(1:1) /= ' ' ) ! !-- Set the precision for data compression. IF ( do3d_compress ) THEN DO i = 1, 100 IF ( plot_3d_precision(i)%variable == do3d(av,if) ) THEN prec = plot_3d_precision(i)%precision EXIT ENDIF ENDDO ENDIF ! !-- Store the array chosen on the temporary array. resorted = .FALSE. SELECT CASE ( TRIM( do3d(av,if) ) ) CASE ( 'e' ) IF ( av == 0 ) THEN to_be_resorted => e ELSE to_be_resorted => e_av ENDIF CASE ( 'lpt' ) IF ( av == 0 ) THEN to_be_resorted => pt ELSE to_be_resorted => lpt_av ENDIF CASE ( 'p' ) IF ( av == 0 ) THEN IF ( psolver /= 'sor' ) CALL exchange_horiz( p, nbgp ) to_be_resorted => p ELSE IF ( psolver /= 'sor' ) CALL exchange_horiz( p_av, nbgp ) to_be_resorted => p_av ENDIF CASE ( 'pc' ) ! particle concentration (requires ghostpoint exchange) IF ( av == 0 ) THEN tend = prt_count CALL exchange_horiz( tend, nbgp ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nz_do3d local_pf(i,j,k) = tend(k,j,i) ENDDO ENDDO ENDDO resorted = .TRUE. ELSE CALL exchange_horiz( pc_av, nbgp ) to_be_resorted => pc_av ENDIF CASE ( 'pr' ) ! mean particle radius IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb, nz_do3d psi = prt_start_index(k,j,i) s_r3 = 0.0 s_r4 = 0.0 DO n = psi, psi+prt_count(k,j,i)-1 s_r3 = s_r3 + particles(n)%radius**3 * & particles(n)%weight_factor s_r4 = s_r4 + particles(n)%radius**4 * & particles(n)%weight_factor ENDDO IF ( s_r3 /= 0.0 ) THEN mean_r = s_r4 / s_r3 ELSE mean_r = 0.0 ENDIF tend(k,j,i) = mean_r ENDDO ENDDO ENDDO CALL exchange_horiz( tend, nbgp ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nz_do3d local_pf(i,j,k) = tend(k,j,i) ENDDO ENDDO ENDDO resorted = .TRUE. ELSE CALL exchange_horiz( pr_av, nbgp ) to_be_resorted => pr_av ENDIF CASE ( 'pt' ) IF ( av == 0 ) THEN IF ( .NOT. cloud_physics ) THEN to_be_resorted => pt ELSE DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nz_do3d local_pf(i,j,k) = pt(k,j,i) + l_d_cp * & pt_d_t(k) * & ql(k,j,i) ENDDO ENDDO ENDDO resorted = .TRUE. ENDIF ELSE to_be_resorted => pt_av ENDIF CASE ( 'q' ) IF ( av == 0 ) THEN to_be_resorted => q ELSE to_be_resorted => q_av ENDIF CASE ( 'ql' ) IF ( av == 0 ) THEN to_be_resorted => ql ELSE to_be_resorted => ql_av ENDIF CASE ( 'ql_c' ) IF ( av == 0 ) THEN to_be_resorted => ql_c ELSE to_be_resorted => ql_c_av ENDIF CASE ( 'ql_v' ) IF ( av == 0 ) THEN to_be_resorted => ql_v ELSE to_be_resorted => ql_v_av ENDIF CASE ( 'ql_vp' ) IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb, nz_do3d psi = prt_start_index(k,j,i) DO n = psi, psi+prt_count(k,j,i)-1 tend(k,j,i) = tend(k,j,i) + & particles(n)%weight_factor / & prt_count(k,j,i) ENDDO ENDDO ENDDO ENDDO CALL exchange_horiz( tend, nbgp ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nz_do3d local_pf(i,j,k) = tend(k,j,i) ENDDO ENDDO ENDDO resorted = .TRUE. ELSE CALL exchange_horiz( ql_vp_av, nbgp ) to_be_resorted => ql_vp_av ENDIF CASE ( 'qv' ) IF ( av == 0 ) THEN DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nz_do3d local_pf(i,j,k) = q(k,j,i) - ql(k,j,i) ENDDO ENDDO ENDDO resorted = .TRUE. ELSE to_be_resorted => qv_av ENDIF CASE ( 'rho' ) IF ( av == 0 ) THEN to_be_resorted => rho ELSE to_be_resorted => rho_av ENDIF CASE ( 's' ) IF ( av == 0 ) THEN to_be_resorted => q ELSE to_be_resorted => s_av ENDIF CASE ( 'sa' ) IF ( av == 0 ) THEN to_be_resorted => sa ELSE to_be_resorted => sa_av ENDIF CASE ( 'u' ) IF ( av == 0 ) THEN to_be_resorted => u ELSE to_be_resorted => u_av ENDIF CASE ( 'v' ) IF ( av == 0 ) THEN to_be_resorted => v ELSE to_be_resorted => v_av ENDIF CASE ( 'vpt' ) IF ( av == 0 ) THEN to_be_resorted => vpt ELSE to_be_resorted => vpt_av ENDIF CASE ( 'w' ) IF ( av == 0 ) THEN to_be_resorted => w ELSE to_be_resorted => w_av ENDIF CASE DEFAULT ! !-- User defined quantity CALL user_data_output_3d( av, do3d(av,if), found, local_pf, & nz_do3d ) resorted = .TRUE. IF ( .NOT. found ) THEN message_string = 'no output available for: ' // & TRIM( do3d(av,if) ) CALL message( 'data_output_3d', 'PA0182', 0, 0, 0, 6, 0 ) ENDIF END SELECT ! !-- Resort the array to be output, if not done above IF ( .NOT. resorted ) THEN DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nz_do3d local_pf(i,j,k) = to_be_resorted(k,j,i) ENDDO ENDDO ENDDO ENDIF ! !-- Output of the volume data information for the AVS-FLD-file. do3d_avs_n = do3d_avs_n + 1 IF ( myid == 0 .AND. avs_output ) THEN ! !-- AVS-labels must not contain any colons. Hence they must be removed !-- from the time character string. simulated_time_mod = simulated_time_chr DO WHILE ( SCAN( simulated_time_mod, ':' ) /= 0 ) pos = SCAN( simulated_time_mod, ':' ) simulated_time_mod(pos:pos) = '/' ENDDO IF ( av == 0 ) THEN WRITE ( 33, 3300 ) do3d_avs_n, TRIM( avs_data_file ), & skip_do_avs, TRIM( do3d(av,if) ), & TRIM( simulated_time_mod ) ELSE WRITE ( 33, 3300 ) do3d_avs_n, TRIM( avs_data_file ), & skip_do_avs, TRIM( do3d(av,if) ) // & ' averaged', TRIM( simulated_time_mod ) ENDIF ! !-- Determine the Skip-value for the next array. Record end and start !-- require 4 byte each. skip_do_avs = skip_do_avs + ( ( ( nx+2*nbgp ) * ( ny+2*nbgp ) * & ( nz_do3d+1 ) ) * 4 + 8 ) ENDIF ! !-- Output of the 3D-array. (compressed/uncompressed) IF ( do3d_compress ) THEN ! !-- Compression, output of compression information on FLD-file and output !-- of compressed data. CALL write_compressed( local_pf, 30, 33, myid, nxl, nxr, nyn, nys, & nzb, nz_do3d, prec, nbgp ) ELSE ! !-- Uncompressed output. #if defined( __parallel ) IF ( netcdf_output ) THEN IF ( netcdf_data_format < 5 ) THEN ! !-- Non-parallel netCDF output. Data is output in parallel in !-- FORTRAN binary format here, and later collected into one file by !-- combine_plot_fields IF ( myid == 0 ) THEN WRITE ( 30 ) time_since_reference_point, & do3d_time_count(av), av ENDIF DO i = 0, io_blocks-1 IF ( i == io_group ) THEN WRITE ( 30 ) nxlg, nxrg, nysg, nyng, nzb, nz_do3d WRITE ( 30 ) local_pf ENDIF #if defined( __parallel ) CALL MPI_BARRIER( comm2d, ierr ) #endif ENDDO ELSE #if defined( __netcdf ) ! !-- Parallel output in netCDF4/HDF5 format. !-- Do not output redundant ghost point data except for the !-- boundaries of the total domain. IF ( nxr == nx .AND. nyn /= ny ) THEN nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & local_pf(nxl:nxrg,nys:nyn,nzb:nz_do3d), & start = (/ nxl+1, nys+1, nzb+1, do3d_time_count(av) /), & count = (/ nxr-nxl+1+nbgp, nyn-nys+1, nz_do3d-nzb+1, 1 /) ) ELSEIF ( nxr /= nx .AND. nyn == ny ) THEN nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & local_pf(nxl:nxr,nys:nyng,nzb:nz_do3d), & start = (/ nxl+1, nys+1, nzb+1, do3d_time_count(av) /), & count = (/ nxr-nxl+1, nyn-nys+1+nbgp, nz_do3d-nzb+1, 1 /) ) ELSEIF ( nxr == nx .AND. nyn == ny ) THEN nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & local_pf(nxl:nxrg,nys:nyng,nzb:nz_do3d), & start = (/ nxl+1, nys+1, nzb+1, do3d_time_count(av) /), & count = (/ nxr-nxl+1+nbgp, nyn-nys+1+nbgp, nz_do3d-nzb+1, 1 /) ) ELSE nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & local_pf(nxl:nxr,nys:nyn,nzb:nz_do3d), & start = (/ nxl+1, nys+1, nzb+1, do3d_time_count(av) /), & count = (/ nxr-nxl+1, nyn-nys+1, nz_do3d-nzb+1, 1 /) ) ENDIF CALL handle_netcdf_error( 'data_output_3d', 386 ) #endif ENDIF ENDIF #else IF ( avs_output ) THEN WRITE ( 30 ) local_pf(nxl:nxr+1,nys:nyn+1,:) ENDIF #if defined( __netcdf ) IF ( netcdf_output ) THEN nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & local_pf(nxl:nxr+1,nys:nyn+1,nzb:nz_do3d), & start = (/ 1, 1, 1, do3d_time_count(av) /), & count = (/ nx+2, ny+2, nz_do3d-nzb+1, 1 /) ) CALL handle_netcdf_error( 'data_output_3d', 446 ) ENDIF #endif #endif ENDIF if = if + 1 ENDDO ! !-- Deallocate temporary array. DEALLOCATE( local_pf ) CALL cpu_log (log_point(14),'data_output_3d','stop','nobarrier') ! !-- Formats. 3300 FORMAT ('variable ',I4,' file=',A,' filetype=unformatted skip=',I12/ & 'label = ',A,A) END SUBROUTINE data_output_3d