#if defined( __ibmy_special ) @PROCESS NOOPTimize #endif SUBROUTINE define_netcdf_header( callmode, extend, av ) !------------------------------------------------------------------------------! ! Current revisions: ! ------------------ ! ! ! Former revisions: ! ----------------- ! $Id: netcdf.f90 601 2010-11-24 16:42:38Z helmke $ ! ! 600 2010-11-24 16:10:51Z raasch ! bugfix concerning check of cross-section levels on netcdf-files to be ! extended (xz,yz) ! ! 564 2010-09-30 13:18:59Z helmke ! nc_precision changed from 40 masks to 1 mask, start number of mask output ! files changed to 201, netcdf message identifiers of masked output changed ! ! 519 2010-03-19 05:30:02Z raasch ! particle number defined as unlimited dimension in case of NetCDF4 output, ! special characters like * and " are now allowed for NetCDF variable names, ! replacement of these characters removed, routine clean_netcdf_varname ! removed ! ! 493 2010-03-01 08:30:24Z raasch ! Extensions for NetCDF4 output ! ! 410 2009-12-04 17:05:40Z letzel ! masked data output ! ! 359 2009-08-19 16:56:44Z letzel ! for extended NetCDF files, the updated title attribute includes an update of ! time_average_text where appropriate. ! Bugfix for extended NetCDF files: In order to avoid 'data mode' errors if ! updated attributes are larger than their original size, NF90_PUT_ATT is called ! in 'define mode' enclosed by NF90_REDEF and NF90_ENDDEF calls. This implies a ! possible performance loss; an alternative strategy would be to ensure equal ! attribute size in a job chain. ! NetCDF unit attribute in timeseries output in case of statistic ! regions added. ! Output of NetCDF messages with aid of message handling routine. ! Output of messages replaced by message handling routine. ! Typographical errors fixed. ! ! 216 2008-11-25 07:12:43Z raasch ! Origin of the xy-coordinate system shifted from the center of the first ! grid cell (indices i=0, j=0) to the south-left corner of this cell. ! ! 189 2008-08-13 17:09:26Z letzel ! consistently allow 100 spectra levels instead of 10 ! bug fix in the determination of the number of output heights for spectra, ! +user-defined spectra ! ! 97 2007-06-21 08:23:15Z raasch ! Grids defined for rho and sa ! ! 48 2007-03-06 12:28:36Z raasch ! Output topography height information (zu_s_inner, zw_s_inner) to 2d-xy and 3d ! datasets ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.12 2006/09/26 19:35:16 raasch ! Bugfix yv coordinates for yz cross sections ! ! Revision 1.1 2005/05/18 15:37:16 raasch ! Initial revision ! ! ! Description: ! ------------ ! In case of extend = .FALSE.: ! Define all necessary dimensions, axes and variables for the different ! NetCDF datasets. This subroutine is called from check_open after a new ! dataset is created. It leaves the open NetCDF files ready to write. ! ! In case of extend = .TRUE.: ! Find out if dimensions and variables of an existing file match the values ! of the actual run. If so, get all necessary informations (ids, etc.) from ! this file. ! ! Parameter av can assume values 0 (non-averaged data) and 1 (time averaged ! data) !------------------------------------------------------------------------------! #if defined( __netcdf ) USE arrays_3d USE constants USE control_parameters USE grid_variables USE indices USE netcdf_control USE pegrid USE particle_attributes USE profil_parameter USE spectrum USE statistics IMPLICIT NONE CHARACTER (LEN=2) :: suffix CHARACTER (LEN=2), INTENT (IN) :: callmode CHARACTER (LEN=3) :: suffix1 CHARACTER (LEN=4) :: grid_x, grid_y, grid_z CHARACTER (LEN=6) :: mode CHARACTER (LEN=10) :: netcdf_var_name, precision, var CHARACTER (LEN=80) :: time_average_text CHARACTER (LEN=2000) :: var_list, var_list_old INTEGER :: av, file_id, i, id_x, id_y, id_z, j, ns, ns_old, nz_old INTEGER, DIMENSION(1) :: id_dim_time_old, id_dim_x_yz_old, & id_dim_y_xz_old, id_dim_zu_sp_old, & id_dim_zu_xy_old, id_dim_zu_3d_old, & id_dim_zu_mask_old LOGICAL :: found LOGICAL, INTENT (INOUT) :: extend LOGICAL, SAVE :: init_netcdf = .FALSE. REAL, DIMENSION(1) :: last_time_coordinate REAL, DIMENSION(:), ALLOCATABLE :: netcdf_data REAL, DIMENSION(:,:), ALLOCATABLE :: netcdf_data_2d ! !-- Initializing actions (return to calling routine check_parameters afterwards) IF ( .NOT. init_netcdf ) THEN ! !-- Check and set accuracy for NetCDF output. First set default value nc_precision = NF90_REAL4 i = 1 DO WHILE ( netcdf_precision(i) /= ' ' ) j = INDEX( netcdf_precision(i), '_' ) IF ( j == 0 ) THEN WRITE ( message_string, * ) 'netcdf_precision must contain a ', & '"_"netcdf_precision(', i, ')="', & TRIM( netcdf_precision(i) ),'"' CALL message( 'define_netcdf_header', 'PA0241', 1, 2, 0, 6, 0 ) ENDIF var = netcdf_precision(i)(1:j-1) precision = netcdf_precision(i)(j+1:) IF ( precision == 'NF90_REAL4' ) THEN j = NF90_REAL4 ELSEIF ( precision == 'NF90_REAL8' ) THEN j = NF90_REAL8 ELSE WRITE ( message_string, * ) 'illegal netcdf precision: ', & 'netcdf_precision(', i, ')="', & TRIM( netcdf_precision(i) ),'"' CALL message( 'define_netcdf_header', 'PA0242', 1, 2, 0, 6, 0 ) ENDIF SELECT CASE ( var ) CASE ( 'xy' ) nc_precision(1) = j CASE ( 'xz' ) nc_precision(2) = j CASE ( 'yz' ) nc_precision(3) = j CASE ( '2d' ) nc_precision(1:3) = j CASE ( '3d' ) nc_precision(4) = j CASE ( 'pr' ) nc_precision(5) = j CASE ( 'ts' ) nc_precision(6) = j CASE ( 'sp' ) nc_precision(7) = j CASE ( 'prt' ) nc_precision(8) = j CASE ( 'masks' ) nc_precision(11) = j CASE ( 'all' ) nc_precision = j CASE DEFAULT WRITE ( message_string, * ) 'unknown variable in inipar ', & 'assignment: netcdf_precision(', i, ')="', & TRIM( netcdf_precision(i) ),'"' CALL message( 'define_netcdf_header', 'PA0243', 1, 2, 0, 6, 0 ) END SELECT i = i + 1 IF ( i > 50 ) EXIT ENDDO init_netcdf = .TRUE. RETURN ENDIF ! !-- Determine the mode to be processed IF ( extend ) THEN mode = callmode // '_ext' ELSE mode = callmode // '_new' ENDIF ! !-- Select the mode to be processed. Possibilities are 3d, mask, xy, xz, yz, !-- pr and ts. SELECT CASE ( mode ) CASE ( 'ma_new' ) ! !-- decompose actual parameter file_id (=formal parameter av) into !-- mid and av file_id = av IF ( file_id <= 200+max_masks ) THEN mid = file_id - 200 av = 0 ELSE mid = file_id - (200+max_masks) av = 1 ENDIF ! !-- Define some global attributes of the dataset nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), NF90_GLOBAL, & 'Conventions', 'COARDS' ) CALL handle_netcdf_error( 'netcdf', 464 ) IF ( av == 0 ) THEN time_average_text = ' ' ELSE WRITE (time_average_text, '('', '',F7.1,'' s average'')') & averaging_interval ENDIF nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 465 ) IF ( av == 1 ) THEN WRITE ( time_average_text,'(F7.1,'' s avg'')' ) averaging_interval nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), NF90_GLOBAL, & 'time_avg', TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 466 ) ENDIF ! !-- Define time coordinate for volume data (unlimited dimension) nc_stat = NF90_DEF_DIM( id_set_mask(mid,av), 'time', NF90_UNLIMITED, & id_dim_time_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 467 ) nc_stat = NF90_DEF_VAR( id_set_mask(mid,av), 'time', NF90_DOUBLE, & id_dim_time_mask(mid,av), & id_var_time_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 468 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), & id_var_time_mask(mid,av), 'units', & 'seconds') CALL handle_netcdf_error( 'netcdf', 469 ) ! !-- Define spatial dimensions and coordinates: !-- Define vertical coordinate grid (zu grid) nc_stat = NF90_DEF_DIM( id_set_mask(mid,av), 'zu_3d', & mask_size(mid,3), id_dim_zu_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 470 ) nc_stat = NF90_DEF_VAR( id_set_mask(mid,av), 'zu_3d', NF90_DOUBLE, & id_dim_zu_mask(mid,av), & id_var_zu_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 471 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), id_var_zu_mask(mid,av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 472 ) ! !-- Define vertical coordinate grid (zw grid) nc_stat = NF90_DEF_DIM( id_set_mask(mid,av), 'zw_3d', & mask_size(mid,3), id_dim_zw_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 473 ) nc_stat = NF90_DEF_VAR( id_set_mask(mid,av), 'zw_3d', NF90_DOUBLE, & id_dim_zw_mask(mid,av), & id_var_zw_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 474 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), id_var_zw_mask(mid,av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 475 ) ! !-- Define x-axis (for scalar position) nc_stat = NF90_DEF_DIM( id_set_mask(mid,av), 'x', & mask_size(mid,1), id_dim_x_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 476 ) nc_stat = NF90_DEF_VAR( id_set_mask(mid,av), 'x', NF90_DOUBLE, & id_dim_x_mask(mid,av), id_var_x_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 477 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), id_var_x_mask(mid,av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 478 ) ! !-- Define x-axis (for u position) nc_stat = NF90_DEF_DIM( id_set_mask(mid,av), 'xu', & mask_size(mid,1), id_dim_xu_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 479 ) nc_stat = NF90_DEF_VAR( id_set_mask(mid,av), 'xu', NF90_DOUBLE, & id_dim_xu_mask(mid,av), & id_var_xu_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 480 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), id_var_xu_mask(mid,av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 481 ) ! !-- Define y-axis (for scalar position) nc_stat = NF90_DEF_DIM( id_set_mask(mid,av), 'y', & mask_size(mid,2), id_dim_y_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 482 ) nc_stat = NF90_DEF_VAR( id_set_mask(mid,av), 'y', NF90_DOUBLE, & id_dim_y_mask(mid,av), id_var_y_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 483 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), id_var_y_mask(mid,av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 484 ) ! !-- Define y-axis (for v position) nc_stat = NF90_DEF_DIM( id_set_mask(mid,av), 'yv', & mask_size(mid,2), id_dim_yv_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 485 ) nc_stat = NF90_DEF_VAR( id_set_mask(mid,av), 'yv', NF90_DOUBLE, & id_dim_yv_mask(mid,av), & id_var_yv_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 486 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), id_var_yv_mask(mid,av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 487 ) ! !-- In case of non-flat topography define 2d-arrays containing the height !-- informations IF ( TRIM( topography ) /= 'flat' ) THEN ! !-- Define zusi = zu(nzb_s_inner) nc_stat = NF90_DEF_VAR( id_set_mask(mid,av), 'zusi', NF90_DOUBLE, & (/ id_dim_x_mask(mid,av), & id_dim_y_mask(mid,av) /), & id_var_zusi_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 488 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), & id_var_zusi_mask(mid,av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 489 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), & id_var_zusi_mask(mid,av), & 'long_name', 'zu(nzb_s_inner)' ) CALL handle_netcdf_error( 'netcdf', 490 ) ! !-- Define zwwi = zw(nzb_w_inner) nc_stat = NF90_DEF_VAR( id_set_mask(mid,av), 'zwwi', NF90_DOUBLE, & (/ id_dim_x_mask(mid,av), & id_dim_y_mask(mid,av) /), & id_var_zwwi_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 491 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), & id_var_zwwi_mask(mid,av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 492 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), & id_var_zwwi_mask(mid,av), & 'long_name', 'zw(nzb_w_inner)' ) CALL handle_netcdf_error( 'netcdf', 493 ) ENDIF ! !-- Define the variables var_list = ';' i = 1 DO WHILE ( domask(mid,av,i)(1:1) /= ' ' ) ! !-- Check for the grid found = .TRUE. SELECT CASE ( domask(mid,av,i) ) ! !-- Most variables are defined on the scalar grid CASE ( 'e', 'p', 'pc', 'pr', 'pt', 'q', 'ql', 'ql_c', 'ql_v', & 'ql_vp', 'qv', 'rho', 's', 'sa', 'vpt' ) grid_x = 'x' grid_y = 'y' grid_z = 'zu' ! !-- u grid CASE ( 'u' ) grid_x = 'xu' grid_y = 'y' grid_z = 'zu' ! !-- v grid CASE ( 'v' ) grid_x = 'x' grid_y = 'yv' grid_z = 'zu' ! !-- w grid CASE ( 'w' ) grid_x = 'x' grid_y = 'y' grid_z = 'zw' CASE DEFAULT ! !-- Check for user-defined quantities CALL user_define_netcdf_grid( domask(mid,av,i), found, & grid_x, grid_y, grid_z ) END SELECT ! !-- Select the respective dimension ids IF ( grid_x == 'x' ) THEN id_x = id_dim_x_mask(mid,av) ELSEIF ( grid_x == 'xu' ) THEN id_x = id_dim_xu_mask(mid,av) ENDIF IF ( grid_y == 'y' ) THEN id_y = id_dim_y_mask(mid,av) ELSEIF ( grid_y == 'yv' ) THEN id_y = id_dim_yv_mask(mid,av) ENDIF IF ( grid_z == 'zu' ) THEN id_z = id_dim_zu_mask(mid,av) ELSEIF ( grid_z == 'zw' ) THEN id_z = id_dim_zw_mask(mid,av) ENDIF ! !-- Define the grid nc_stat = NF90_DEF_VAR( id_set_mask(mid,av), domask(mid,av,i), & nc_precision(11), & (/ id_x, id_y, id_z, & id_dim_time_mask(mid,av) /), & id_var_domask(mid,av,i) ) IF ( .NOT. found ) THEN WRITE ( message_string, * ) 'no grid defined for', & ' variable ', TRIM( domask(mid,av,i) ) CALL message( 'define_netcdf_header', 'PA0244', 0, 1, 0, 6, 0 ) ENDIF var_list = TRIM( var_list ) // TRIM( domask(mid,av,i) ) // ';' CALL handle_netcdf_error( 'netcdf', 494 ) ! !-- Store the 'real' name of the variable (with *, for example) !-- in the long_name attribute. This is evaluated by Ferret, !-- for example. nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), & id_var_domask(mid,av,i), & 'long_name', domask(mid,av,i) ) CALL handle_netcdf_error( 'netcdf', 495 ) ! !-- Define the variable's unit nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), & id_var_domask(mid,av,i), & 'units', TRIM( domask_unit(mid,av,i) ) ) CALL handle_netcdf_error( 'netcdf', 496 ) i = i + 1 ENDDO ! !-- No arrays to output IF ( i == 1 ) RETURN ! !-- Write the list of variables as global attribute (this is used by !-- restart runs and by combine_plot_fields) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), NF90_GLOBAL, & 'VAR_LIST', var_list ) CALL handle_netcdf_error( 'netcdf', 497 ) ! !-- Leave NetCDF define mode nc_stat = NF90_ENDDEF( id_set_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 498 ) ! !-- Write data for x (shifted by +dx/2) and xu axis ALLOCATE( netcdf_data(mask_size(mid,1)) ) netcdf_data = ( mask_i_global(mid,:mask_size(mid,1)) + 0.5 ) * dx nc_stat = NF90_PUT_VAR( id_set_mask(mid,av), id_var_x_mask(mid,av), & netcdf_data, start = (/ 1 /), & count = (/ mask_size(mid,1) /) ) CALL handle_netcdf_error( 'netcdf', 499 ) netcdf_data = mask_i_global(mid,:mask_size(mid,1)) * dx nc_stat = NF90_PUT_VAR( id_set_mask(mid,av), id_var_xu_mask(mid,av),& netcdf_data, start = (/ 1 /), & count = (/ mask_size(mid,1) /) ) CALL handle_netcdf_error( 'netcdf', 500 ) DEALLOCATE( netcdf_data ) ! !-- Write data for y (shifted by +dy/2) and yv axis ALLOCATE( netcdf_data(mask_size(mid,2)) ) netcdf_data = ( mask_j_global(mid,:mask_size(mid,2)) + 0.5 ) * dy nc_stat = NF90_PUT_VAR( id_set_mask(mid,av), id_var_y_mask(mid,av), & netcdf_data, start = (/ 1 /), & count = (/ mask_size(mid,2) /)) CALL handle_netcdf_error( 'netcdf', 501 ) netcdf_data = mask_j_global(mid,:mask_size(mid,2)) * dy nc_stat = NF90_PUT_VAR( id_set_mask(mid,av), id_var_yv_mask(mid,av), & netcdf_data, start = (/ 1 /), & count = (/ mask_size(mid,2) /)) CALL handle_netcdf_error( 'netcdf', 502 ) DEALLOCATE( netcdf_data ) ! !-- Write zu and zw data (vertical axes) ALLOCATE( netcdf_data(mask_size(mid,3)) ) netcdf_data = zu( mask_k_global(mid,:mask_size(mid,3)) ) nc_stat = NF90_PUT_VAR( id_set_mask(mid,av), id_var_zu_mask(mid,av), & netcdf_data, start = (/ 1 /), & count = (/ mask_size(mid,3) /) ) CALL handle_netcdf_error( 'netcdf', 503 ) netcdf_data = zw( mask_k_global(mid,:mask_size(mid,3)) ) nc_stat = NF90_PUT_VAR( id_set_mask(mid,av), id_var_zw_mask(mid,av), & netcdf_data, start = (/ 1 /), & count = (/ mask_size(mid,3) /) ) CALL handle_netcdf_error( 'netcdf', 504 ) DEALLOCATE( netcdf_data ) ! !-- In case of non-flat topography write height information IF ( TRIM( topography ) /= 'flat' ) THEN ALLOCATE( netcdf_data_2d(mask_size(mid,1),mask_size(mid,2)) ) netcdf_data_2d = zu_s_inner( mask_i_global(mid,:mask_size(mid,1)),& mask_j_global(mid,:mask_size(mid,2)) ) nc_stat = NF90_PUT_VAR( id_set_mask(mid,av), & id_var_zusi_mask(mid,av), & netcdf_data_2d, & start = (/ 1, 1 /), & count = (/ mask_size(mid,1), & mask_size(mid,2) /) ) CALL handle_netcdf_error( 'netcdf', 505 ) netcdf_data_2d = zw_w_inner( mask_i_global(mid,:mask_size(mid,1)),& mask_j_global(mid,:mask_size(mid,2)) ) nc_stat = NF90_PUT_VAR( id_set_mask(mid,av), & id_var_zwwi_mask(mid,av), & netcdf_data_2d, & start = (/ 1, 1 /), & count = (/ mask_size(mid,1), & mask_size(mid,2) /) ) CALL handle_netcdf_error( 'netcdf', 506 ) DEALLOCATE( netcdf_data_2d ) ENDIF ! !-- restore original parameter file_id (=formal parameter av) into av av = file_id CASE ( 'ma_ext' ) ! !-- decompose actual parameter file_id (=formal parameter av) into !-- mid and av file_id = av IF ( file_id <= 200+max_masks ) THEN mid = file_id - 200 av = 0 ELSE mid = file_id - (200+max_masks) av = 1 ENDIF ! !-- Get the list of variables and compare with the actual run. !-- First var_list_old has to be reset, since GET_ATT does not assign !-- trailing blanks. var_list_old = ' ' nc_stat = NF90_GET_ATT( id_set_mask(mid,av), NF90_GLOBAL, 'VAR_LIST',& var_list_old ) CALL handle_netcdf_error( 'netcdf', 507 ) var_list = ';' i = 1 DO WHILE ( domask(mid,av,i)(1:1) /= ' ' ) var_list = TRIM(var_list) // TRIM( domask(mid,av,i) ) // ';' i = i + 1 ENDDO IF ( av == 0 ) THEN var = '(mask)' ELSE var = '(mask_av)' ENDIF IF ( TRIM( var_list ) /= TRIM( var_list_old ) ) THEN WRITE ( message_string, * ) 'NetCDF file for ', TRIM( var ), & ' data for mask', mid, ' from previous run found,', & '&but this file cannot be extended due to variable ', & 'mismatch.&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0335', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Get and compare the number of vertical gridpoints nc_stat = NF90_INQ_VARID( id_set_mask(mid,av), 'zu', & id_var_zu_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 508 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_mask(mid,av), & id_var_zu_mask(mid,av), & dimids = id_dim_zu_mask_old ) CALL handle_netcdf_error( 'netcdf', 509 ) id_dim_zu_mask(mid,av) = id_dim_zu_mask_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_mask(mid,av), & id_dim_zu_mask(mid,av), & len = nz_old ) CALL handle_netcdf_error( 'netcdf', 510 ) IF ( mask_size(mid,3) /= nz_old ) THEN WRITE ( message_string, * ) 'NetCDF file for ', TRIM( var ), & ' data for mask', mid, ' from previous run found,', & '&but this file cannot be extended due to mismatch in ', & ' number of&vertical grid points.', & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0336', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Get the id of the time coordinate (unlimited coordinate) and its !-- last index on the file. The next time level is plmask..count+1. !-- The current time must be larger than the last output time !-- on the file. nc_stat = NF90_INQ_VARID( id_set_mask(mid,av), 'time', & id_var_time_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 511 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_mask(mid,av), & id_var_time_mask(mid,av), & dimids = id_dim_time_old ) CALL handle_netcdf_error( 'netcdf', 512 ) id_dim_time_mask(mid,av) = id_dim_time_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_mask(mid,av), & id_dim_time_mask(mid,av), & len = domask_time_count(mid,av) ) CALL handle_netcdf_error( 'netcdf', 513 ) nc_stat = NF90_GET_VAR( id_set_mask(mid,av), & id_var_time_mask(mid,av), & last_time_coordinate, & start = (/ domask_time_count(mid,av) /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'netcdf', 514 ) IF ( last_time_coordinate(1) >= simulated_time ) THEN WRITE ( message_string, * ) 'NetCDF file for ', TRIM( var ), & ' data for mask', mid, ' from previous run found,', & '&but this file cannot be extended because the current ', & 'output time&is less or equal than the last output time ', & 'on this file.&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0337', 0, 1, 0, 6, 0 ) domask_time_count(mid,av) = 0 extend = .FALSE. RETURN ENDIF ! !-- Dataset seems to be extendable. !-- Now get the variable ids. i = 1 DO WHILE ( domask(mid,av,i)(1:1) /= ' ' ) nc_stat = NF90_INQ_VARID( id_set_mask(mid,av), & TRIM( domask(mid,av,i) ), & id_var_domask(mid,av,i) ) CALL handle_netcdf_error( 'netcdf', 515 ) i = i + 1 ENDDO ! !-- Update the title attribute on file !-- In order to avoid 'data mode' errors if updated attributes are larger !-- than their original size, NF90_PUT_ATT is called in 'define mode' !-- enclosed by NF90_REDEF and NF90_ENDDEF calls. This implies a possible !-- performance loss due to data copying; an alternative strategy would be !-- to ensure equal attribute size in a job chain. Maybe revise later. IF ( av == 0 ) THEN time_average_text = ' ' ELSE WRITE (time_average_text, '('', '',F7.1,'' s average'')') & averaging_interval ENDIF nc_stat = NF90_REDEF( id_set_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 516 ) nc_stat = NF90_PUT_ATT( id_set_mask(mid,av), NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 517 ) nc_stat = NF90_ENDDEF( id_set_mask(mid,av) ) CALL handle_netcdf_error( 'netcdf', 518 ) WRITE ( message_string, * ) 'NetCDF file for ', TRIM( var ), & ' data for mask', mid, ' from previous run found.', & '&This file will be extended.' CALL message( 'define_netcdf_header', 'PA0338', 0, 0, 0, 6, 0 ) ! !-- restore original parameter file_id (=formal parameter av) into av av = file_id CASE ( '3d_new' ) ! !-- Define some global attributes of the dataset nc_stat = NF90_PUT_ATT( id_set_3d(av), NF90_GLOBAL, 'Conventions', & 'COARDS' ) CALL handle_netcdf_error( 'netcdf', 62 ) IF ( av == 0 ) THEN time_average_text = ' ' ELSE WRITE (time_average_text, '('', '',F7.1,'' s average'')') & averaging_interval ENDIF nc_stat = NF90_PUT_ATT( id_set_3d(av), NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 63 ) IF ( av == 1 ) THEN WRITE ( time_average_text,'(F7.1,'' s avg'')' ) averaging_interval nc_stat = NF90_PUT_ATT( id_set_3d(av), NF90_GLOBAL, 'time_avg', & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 63 ) ENDIF ! !-- Define time coordinate for volume data (unlimited dimension) nc_stat = NF90_DEF_DIM( id_set_3d(av), 'time', NF90_UNLIMITED, & id_dim_time_3d(av) ) CALL handle_netcdf_error( 'netcdf', 64 ) nc_stat = NF90_DEF_VAR( id_set_3d(av), 'time', NF90_DOUBLE, & id_dim_time_3d(av), id_var_time_3d(av) ) CALL handle_netcdf_error( 'netcdf', 65 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_time_3d(av), 'units', & 'seconds') CALL handle_netcdf_error( 'netcdf', 66 ) ! !-- Define spatial dimensions and coordinates: !-- Define vertical coordinate grid (zu grid) nc_stat = NF90_DEF_DIM( id_set_3d(av), 'zu_3d', nz_do3d-nzb+1, & id_dim_zu_3d(av) ) CALL handle_netcdf_error( 'netcdf', 67 ) nc_stat = NF90_DEF_VAR( id_set_3d(av), 'zu_3d', NF90_DOUBLE, & id_dim_zu_3d(av), id_var_zu_3d(av) ) CALL handle_netcdf_error( 'netcdf', 68 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_zu_3d(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 69 ) ! !-- Define vertical coordinate grid (zw grid) nc_stat = NF90_DEF_DIM( id_set_3d(av), 'zw_3d', nz_do3d-nzb+1, & id_dim_zw_3d(av) ) CALL handle_netcdf_error( 'netcdf', 70 ) nc_stat = NF90_DEF_VAR( id_set_3d(av), 'zw_3d', NF90_DOUBLE, & id_dim_zw_3d(av), id_var_zw_3d(av) ) CALL handle_netcdf_error( 'netcdf', 71 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_zw_3d(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 72 ) ! !-- Define x-axis (for scalar position) nc_stat = NF90_DEF_DIM( id_set_3d(av), 'x', nx+2, id_dim_x_3d(av) ) CALL handle_netcdf_error( 'netcdf', 73 ) nc_stat = NF90_DEF_VAR( id_set_3d(av), 'x', NF90_DOUBLE, & id_dim_x_3d(av), id_var_x_3d(av) ) CALL handle_netcdf_error( 'netcdf', 74 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_x_3d(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 75 ) ! !-- Define x-axis (for u position) nc_stat = NF90_DEF_DIM( id_set_3d(av), 'xu', nx+2, id_dim_xu_3d(av) ) CALL handle_netcdf_error( 'netcdf', 358 ) nc_stat = NF90_DEF_VAR( id_set_3d(av), 'xu', NF90_DOUBLE, & id_dim_xu_3d(av), id_var_xu_3d(av) ) CALL handle_netcdf_error( 'netcdf', 359 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_xu_3d(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 360 ) ! !-- Define y-axis (for scalar position) nc_stat = NF90_DEF_DIM( id_set_3d(av), 'y', ny+2, id_dim_y_3d(av) ) CALL handle_netcdf_error( 'netcdf', 76 ) nc_stat = NF90_DEF_VAR( id_set_3d(av), 'y', NF90_DOUBLE, & id_dim_y_3d(av), id_var_y_3d(av) ) CALL handle_netcdf_error( 'netcdf', 77 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_y_3d(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 78 ) ! !-- Define y-axis (for v position) nc_stat = NF90_DEF_DIM( id_set_3d(av), 'yv', ny+2, id_dim_yv_3d(av) ) CALL handle_netcdf_error( 'netcdf', 361 ) nc_stat = NF90_DEF_VAR( id_set_3d(av), 'yv', NF90_DOUBLE, & id_dim_yv_3d(av), id_var_yv_3d(av) ) CALL handle_netcdf_error( 'netcdf', 362 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_yv_3d(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 363 ) ! !-- In case of non-flat topography define 2d-arrays containing the height !-- informations IF ( TRIM( topography ) /= 'flat' ) THEN ! !-- Define zusi = zu(nzb_s_inner) nc_stat = NF90_DEF_VAR( id_set_3d(av), 'zusi', NF90_DOUBLE, & (/ id_dim_x_3d(av), id_dim_y_3d(av) /), & id_var_zusi_3d(av) ) CALL handle_netcdf_error( 'netcdf', 413 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_zusi_3d(av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 414 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_zusi_3d(av), & 'long_name', 'zu(nzb_s_inner)' ) CALL handle_netcdf_error( 'netcdf', 415 ) ! !-- Define zwwi = zw(nzb_w_inner) nc_stat = NF90_DEF_VAR( id_set_3d(av), 'zwwi', NF90_DOUBLE, & (/ id_dim_x_3d(av), id_dim_y_3d(av) /), & id_var_zwwi_3d(av) ) CALL handle_netcdf_error( 'netcdf', 416 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_zwwi_3d(av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 417 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_zwwi_3d(av), & 'long_name', 'zw(nzb_w_inner)' ) CALL handle_netcdf_error( 'netcdf', 418 ) ENDIF ! !-- Define the variables var_list = ';' i = 1 DO WHILE ( do3d(av,i)(1:1) /= ' ' ) ! !-- Check for the grid found = .TRUE. SELECT CASE ( do3d(av,i) ) ! !-- Most variables are defined on the scalar grid CASE ( 'e', 'p', 'pc', 'pr', 'pt', 'q', 'ql', 'ql_c', 'ql_v', & 'ql_vp', 'qv', 'rho', 's', 'sa', 'vpt' ) grid_x = 'x' grid_y = 'y' grid_z = 'zu' ! !-- u grid CASE ( 'u' ) grid_x = 'xu' grid_y = 'y' grid_z = 'zu' ! !-- v grid CASE ( 'v' ) grid_x = 'x' grid_y = 'yv' grid_z = 'zu' ! !-- w grid CASE ( 'w' ) grid_x = 'x' grid_y = 'y' grid_z = 'zw' CASE DEFAULT ! !-- Check for user-defined quantities CALL user_define_netcdf_grid( do3d(av,i), found, grid_x, & grid_y, grid_z ) END SELECT ! !-- Select the respective dimension ids IF ( grid_x == 'x' ) THEN id_x = id_dim_x_3d(av) ELSEIF ( grid_x == 'xu' ) THEN id_x = id_dim_xu_3d(av) ENDIF IF ( grid_y == 'y' ) THEN id_y = id_dim_y_3d(av) ELSEIF ( grid_y == 'yv' ) THEN id_y = id_dim_yv_3d(av) ENDIF IF ( grid_z == 'zu' ) THEN id_z = id_dim_zu_3d(av) ELSEIF ( grid_z == 'zw' ) THEN id_z = id_dim_zw_3d(av) ENDIF ! !-- Define the grid nc_stat = NF90_DEF_VAR( id_set_3d(av), do3d(av,i), & nc_precision(4), & (/ id_x, id_y, id_z, id_dim_time_3d(av) /), & id_var_do3d(av,i) ) IF ( .NOT. found ) THEN WRITE ( message_string, * ) 'no grid defined for', & ' variable ', TRIM( do3d(av,i) ) CALL message( 'define_netcdf_header', 'PA0244', 0, 1, 0, 6, 0 ) ENDIF var_list = TRIM( var_list ) // TRIM( do3d(av,i) ) // ';' CALL handle_netcdf_error( 'netcdf', 79 ) ! !-- Store the 'real' name of the variable (with *, for example) !-- in the long_name attribute. This is evaluated by Ferret, !-- for example. nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_do3d(av,i), & 'long_name', do3d(av,i) ) CALL handle_netcdf_error( 'netcdf', 80 ) ! !-- Define the variable's unit nc_stat = NF90_PUT_ATT( id_set_3d(av), id_var_do3d(av,i), & 'units', TRIM( do3d_unit(av,i) ) ) CALL handle_netcdf_error( 'netcdf', 357 ) #if defined( __netcdf4 ) ! !-- Set collective io operations for parallel io IF ( netcdf_data_format > 2 ) THEN nc_stat = NF90_VAR_PAR_ACCESS( id_set_3d(av), & id_var_do3d(av,i), & NF90_COLLECTIVE ) CALL handle_netcdf_error( 'netcdf', 445 ) ENDIF #endif i = i + 1 ENDDO ! !-- No arrays to output IF ( i == 1 ) RETURN ! !-- Write the list of variables as global attribute (this is used by !-- restart runs and by combine_plot_fields) nc_stat = NF90_PUT_ATT( id_set_3d(av), NF90_GLOBAL, 'VAR_LIST', & var_list ) CALL handle_netcdf_error( 'netcdf', 81 ) ! !-- Leave NetCDF define mode nc_stat = NF90_ENDDEF( id_set_3d(av) ) CALL handle_netcdf_error( 'netcdf', 82 ) ! !-- Write data for x (shifted by +dx/2) and xu axis ALLOCATE( netcdf_data(0:nx+1) ) DO i = 0, nx+1 netcdf_data(i) = ( i + 0.5 ) * dx ENDDO nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_x_3d(av), netcdf_data, & start = (/ 1 /), count = (/ nx+2 /) ) CALL handle_netcdf_error( 'netcdf', 83 ) DO i = 0, nx+1 netcdf_data(i) = i * dx ENDDO nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_xu_3d(av), & netcdf_data, start = (/ 1 /), & count = (/ nx+2 /) ) CALL handle_netcdf_error( 'netcdf', 385 ) DEALLOCATE( netcdf_data ) ! !-- Write data for y (shifted by +dy/2) and yv axis ALLOCATE( netcdf_data(0:ny+1) ) DO i = 0, ny+1 netcdf_data(i) = ( i + 0.5 ) * dy ENDDO nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_y_3d(av), netcdf_data, & start = (/ 1 /), count = (/ ny+2 /)) CALL handle_netcdf_error( 'netcdf', 84 ) DO i = 0, ny+1 netcdf_data(i) = i * dy ENDDO nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_yv_3d(av), & netcdf_data, start = (/ 1 /), & count = (/ ny+2 /)) CALL handle_netcdf_error( 'netcdf', 387 ) DEALLOCATE( netcdf_data ) ! !-- Write zu and zw data (vertical axes) nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_zu_3d(av), & zu(nzb:nz_do3d), start = (/ 1 /), & count = (/ nz_do3d-nzb+1 /) ) CALL handle_netcdf_error( 'netcdf', 85 ) nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_zw_3d(av), & zw(nzb:nz_do3d), start = (/ 1 /), & count = (/ nz_do3d-nzb+1 /) ) CALL handle_netcdf_error( 'netcdf', 86 ) ! !-- In case of non-flat topography write height information IF ( TRIM( topography ) /= 'flat' ) THEN nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_zusi_3d(av), & zu_s_inner(0:nx+1,0:ny+1), & start = (/ 1, 1 /), & count = (/ nx+2, ny+2 /) ) CALL handle_netcdf_error( 'netcdf', 419 ) nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_zwwi_3d(av), & zw_w_inner(0:nx+1,0:ny+1), & start = (/ 1, 1 /), & count = (/ nx+2, ny+2 /) ) CALL handle_netcdf_error( 'netcdf', 420 ) ENDIF CASE ( '3d_ext' ) ! !-- Get the list of variables and compare with the actual run. !-- First var_list_old has to be reset, since GET_ATT does not assign !-- trailing blanks. var_list_old = ' ' nc_stat = NF90_GET_ATT( id_set_3d(av), NF90_GLOBAL, 'VAR_LIST', & var_list_old ) CALL handle_netcdf_error( 'netcdf', 87 ) var_list = ';' i = 1 DO WHILE ( do3d(av,i)(1:1) /= ' ' ) var_list = TRIM(var_list) // TRIM( do3d(av,i) ) // ';' i = i + 1 ENDDO IF ( av == 0 ) THEN var = '(3d)' ELSE var = '(3d_av)' ENDIF IF ( TRIM( var_list ) /= TRIM( var_list_old ) ) THEN message_string = 'NetCDF file for volume data ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended due to' // & ' variable mismatch.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0245', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Get and compare the number of vertical gridpoints nc_stat = NF90_INQ_VARID( id_set_3d(av), 'zu_3d', id_var_zu_3d(av) ) CALL handle_netcdf_error( 'netcdf', 88 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_3d(av), id_var_zu_3d(av), & dimids = id_dim_zu_3d_old ) CALL handle_netcdf_error( 'netcdf', 89 ) id_dim_zu_3d(av) = id_dim_zu_3d_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_3d(av), id_dim_zu_3d(av), & len = nz_old ) CALL handle_netcdf_error( 'netcdf', 90 ) IF ( nz_do3d-nzb+1 /= nz_old ) THEN message_string = 'NetCDF file for volume data ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended due to' // & ' mismatch in number of' // & '&vertical grid points (nz_do3d).' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0246', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Get the id of the time coordinate (unlimited coordinate) and its !-- last index on the file. The next time level is pl3d..count+1. !-- The current time must be larger than the last output time !-- on the file. nc_stat = NF90_INQ_VARID( id_set_3d(av), 'time', id_var_time_3d(av) ) CALL handle_netcdf_error( 'netcdf', 91 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_3d(av), id_var_time_3d(av), & dimids = id_dim_time_old ) CALL handle_netcdf_error( 'netcdf', 92 ) id_dim_time_3d(av) = id_dim_time_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_3d(av), id_dim_time_3d(av), & len = do3d_time_count(av) ) CALL handle_netcdf_error( 'netcdf', 93 ) nc_stat = NF90_GET_VAR( id_set_3d(av), id_var_time_3d(av), & last_time_coordinate, & start = (/ do3d_time_count(av) /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'netcdf', 94 ) IF ( last_time_coordinate(1) >= simulated_time ) THEN message_string = 'NetCDF file for volume data ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended becaus' // & 'e the current output time' // & '&is less or equal than the last output t' // & 'ime on this file.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0247', 0, 1, 0, 6, 0 ) do3d_time_count(av) = 0 extend = .FALSE. RETURN ENDIF ! !-- Dataset seems to be extendable. !-- Now get the variable ids. i = 1 DO WHILE ( do3d(av,i)(1:1) /= ' ' ) nc_stat = NF90_INQ_VARID( id_set_3d(av), TRIM( do3d(av,i) ), & id_var_do3d(av,i) ) CALL handle_netcdf_error( 'netcdf', 95 ) #if defined( __netcdf4 ) ! !-- Set collective io operations for parallel io IF ( netcdf_data_format > 2 ) THEN nc_stat = NF90_VAR_PAR_ACCESS( id_set_3d(av), & id_var_do3d(av,i), & NF90_COLLECTIVE ) CALL handle_netcdf_error( 'netcdf', 453 ) ENDIF #endif i = i + 1 ENDDO ! !-- Update the title attribute on file !-- In order to avoid 'data mode' errors if updated attributes are larger !-- than their original size, NF90_PUT_ATT is called in 'define mode' !-- enclosed by NF90_REDEF and NF90_ENDDEF calls. This implies a possible !-- performance loss due to data copying; an alternative strategy would be !-- to ensure equal attribute size. Maybe revise later. IF ( av == 0 ) THEN time_average_text = ' ' ELSE WRITE (time_average_text, '('', '',F7.1,'' s average'')') & averaging_interval ENDIF nc_stat = NF90_REDEF( id_set_3d(av) ) CALL handle_netcdf_error( 'netcdf', 429 ) nc_stat = NF90_PUT_ATT( id_set_3d(av), NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 96 ) nc_stat = NF90_ENDDEF( id_set_3d(av) ) CALL handle_netcdf_error( 'netcdf', 430 ) message_string = 'NetCDF file for volume data ' // & TRIM( var ) // ' from previous run found.' // & '&This file will be extended.' CALL message( 'define_netcdf_header', 'PA0248', 0, 0, 0, 6, 0 ) CASE ( 'xy_new' ) ! !-- Define some global attributes of the dataset nc_stat = NF90_PUT_ATT( id_set_xy(av), NF90_GLOBAL, 'Conventions', & 'COARDS' ) CALL handle_netcdf_error( 'netcdf', 97 ) IF ( av == 0 ) THEN time_average_text = ' ' ELSE WRITE (time_average_text, '('', '',F7.1,'' s average'')') & averaging_interval ENDIF nc_stat = NF90_PUT_ATT( id_set_xy(av), NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 98 ) IF ( av == 1 ) THEN WRITE ( time_average_text,'(F7.1,'' s avg'')' ) averaging_interval nc_stat = NF90_PUT_ATT( id_set_xy(av), NF90_GLOBAL, 'time_avg', & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 98 ) ENDIF ! !-- Define time coordinate for xy sections (unlimited dimension) nc_stat = NF90_DEF_DIM( id_set_xy(av), 'time', NF90_UNLIMITED, & id_dim_time_xy(av) ) CALL handle_netcdf_error( 'netcdf', 99 ) nc_stat = NF90_DEF_VAR( id_set_xy(av), 'time', NF90_DOUBLE, & id_dim_time_xy(av), id_var_time_xy(av) ) CALL handle_netcdf_error( 'netcdf', 100 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_time_xy(av), 'units', & 'seconds') CALL handle_netcdf_error( 'netcdf', 101 ) ! !-- Define the spatial dimensions and coordinates for xy-sections. !-- First, determine the number of horizontal sections to be written. IF ( section(1,1) == -9999 ) THEN RETURN ELSE ns = 1 DO WHILE ( section(ns,1) /= -9999 .AND. ns <= 100 ) ns = ns + 1 ENDDO ns = ns - 1 ENDIF ! !-- Define vertical coordinate grid (zu grid) nc_stat = NF90_DEF_DIM( id_set_xy(av), 'zu_xy', ns, id_dim_zu_xy(av) ) CALL handle_netcdf_error( 'netcdf', 102 ) nc_stat = NF90_DEF_VAR( id_set_xy(av), 'zu_xy', NF90_DOUBLE, & id_dim_zu_xy(av), id_var_zu_xy(av) ) CALL handle_netcdf_error( 'netcdf', 103 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_zu_xy(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 104 ) ! !-- Define vertical coordinate grid (zw grid) nc_stat = NF90_DEF_DIM( id_set_xy(av), 'zw_xy', ns, id_dim_zw_xy(av) ) CALL handle_netcdf_error( 'netcdf', 105 ) nc_stat = NF90_DEF_VAR( id_set_xy(av), 'zw_xy', NF90_DOUBLE, & id_dim_zw_xy(av), id_var_zw_xy(av) ) CALL handle_netcdf_error( 'netcdf', 106 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_zw_xy(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 107 ) ! !-- Define a pseudo vertical coordinate grid for the surface variables !-- u* and t* to store their height level nc_stat = NF90_DEF_DIM( id_set_xy(av), 'zu1_xy', 1, & id_dim_zu1_xy(av) ) CALL handle_netcdf_error( 'netcdf', 108 ) nc_stat = NF90_DEF_VAR( id_set_xy(av), 'zu1_xy', NF90_DOUBLE, & id_dim_zu1_xy(av), id_var_zu1_xy(av) ) CALL handle_netcdf_error( 'netcdf', 109 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_zu1_xy(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 110 ) ! !-- Define a variable to store the layer indices of the horizontal cross !-- sections, too nc_stat = NF90_DEF_VAR( id_set_xy(av), 'ind_z_xy', NF90_DOUBLE, & id_dim_zu_xy(av), id_var_ind_z_xy(av) ) CALL handle_netcdf_error( 'netcdf', 111 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_ind_z_xy(av), 'units', & 'gridpoints') CALL handle_netcdf_error( 'netcdf', 112 ) ! !-- Define x-axis (for scalar position) nc_stat = NF90_DEF_DIM( id_set_xy(av), 'x', nx+2, id_dim_x_xy(av) ) CALL handle_netcdf_error( 'netcdf', 113 ) nc_stat = NF90_DEF_VAR( id_set_xy(av), 'x', NF90_DOUBLE, & id_dim_x_xy(av), id_var_x_xy(av) ) CALL handle_netcdf_error( 'netcdf', 114 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_x_xy(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 115 ) ! !-- Define x-axis (for u position) nc_stat = NF90_DEF_DIM( id_set_xy(av), 'xu', nx+2, id_dim_xu_xy(av) ) CALL handle_netcdf_error( 'netcdf', 388 ) nc_stat = NF90_DEF_VAR( id_set_xy(av), 'xu', NF90_DOUBLE, & id_dim_xu_xy(av), id_var_xu_xy(av) ) CALL handle_netcdf_error( 'netcdf', 389 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_xu_xy(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 390 ) ! !-- Define y-axis (for scalar position) nc_stat = NF90_DEF_DIM( id_set_xy(av), 'y', ny+2, id_dim_y_xy(av) ) CALL handle_netcdf_error( 'netcdf', 116 ) nc_stat = NF90_DEF_VAR( id_set_xy(av), 'y', NF90_DOUBLE, & id_dim_y_xy(av), id_var_y_xy(av) ) CALL handle_netcdf_error( 'netcdf', 117 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_y_xy(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 118 ) ! !-- Define y-axis (for scalar position) nc_stat = NF90_DEF_DIM( id_set_xy(av), 'yv', ny+2, id_dim_yv_xy(av) ) CALL handle_netcdf_error( 'netcdf', 364 ) nc_stat = NF90_DEF_VAR( id_set_xy(av), 'yv', NF90_DOUBLE, & id_dim_yv_xy(av), id_var_yv_xy(av) ) CALL handle_netcdf_error( 'netcdf', 365 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_yv_xy(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 366 ) ! !-- In case of non-flat topography define 2d-arrays containing the height !-- informations IF ( TRIM( topography ) /= 'flat' ) THEN ! !-- Define zusi = zu(nzb_s_inner) nc_stat = NF90_DEF_VAR( id_set_xy(av), 'zusi', NF90_DOUBLE, & (/ id_dim_x_xy(av), id_dim_y_xy(av) /), & id_var_zusi_xy(av) ) CALL handle_netcdf_error( 'netcdf', 421 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_zusi_xy(av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 422 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_zusi_xy(av), & 'long_name', 'zu(nzb_s_inner)' ) CALL handle_netcdf_error( 'netcdf', 423 ) ! !-- Define zwwi = zw(nzb_w_inner) nc_stat = NF90_DEF_VAR( id_set_xy(av), 'zwwi', NF90_DOUBLE, & (/ id_dim_x_xy(av), id_dim_y_xy(av) /), & id_var_zwwi_xy(av) ) CALL handle_netcdf_error( 'netcdf', 424 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_zwwi_xy(av), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 425 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_zwwi_xy(av), & 'long_name', 'zw(nzb_w_inner)' ) CALL handle_netcdf_error( 'netcdf', 426 ) ENDIF ! !-- Define the variables var_list = ';' i = 1 DO WHILE ( do2d(av,i)(1:1) /= ' ' ) IF ( INDEX( do2d(av,i), 'xy' ) /= 0 ) THEN ! !-- If there is a star in the variable name (u* or t*), it is a !-- surface variable. Define it with id_dim_zu1_xy. IF ( INDEX( do2d(av,i), '*' ) /= 0 ) THEN nc_stat = NF90_DEF_VAR( id_set_xy(av), do2d(av,i), & nc_precision(1), & (/ id_dim_x_xy(av), id_dim_y_xy(av),& id_dim_zu1_xy(av), & id_dim_time_xy(av) /), & id_var_do2d(av,i) ) var_list = TRIM(var_list) // TRIM( do2d(av,i) ) // ';' ELSE ! !-- Check for the grid found = .TRUE. SELECT CASE ( do2d(av,i) ) ! !-- Most variables are defined on the zu grid CASE ( 'e_xy', 'p_xy', 'pc_xy', 'pr_xy', 'pt_xy', 'q_xy',& 'ql_xy', 'ql_c_xy', 'ql_v_xy', 'ql_vp_xy', & 'qv_xy', 'rho_xy', 's_xy', 'sa_xy', 'vpt_xy' ) grid_x = 'x' grid_y = 'y' grid_z = 'zu' ! !-- u grid CASE ( 'u_xy' ) grid_x = 'xu' grid_y = 'y' grid_z = 'zu' ! !-- v grid CASE ( 'v_xy' ) grid_x = 'x' grid_y = 'yv' grid_z = 'zu' ! !-- w grid CASE ( 'w_xy' ) grid_x = 'x' grid_y = 'y' grid_z = 'zw' CASE DEFAULT ! !-- Check for user-defined quantities CALL user_define_netcdf_grid( do2d(av,i), found, & grid_x, grid_y, grid_z ) END SELECT ! !-- Select the respective dimension ids IF ( grid_x == 'x' ) THEN id_x = id_dim_x_xy(av) ELSEIF ( grid_x == 'xu' ) THEN id_x = id_dim_xu_xy(av) ENDIF IF ( grid_y == 'y' ) THEN id_y = id_dim_y_xy(av) ELSEIF ( grid_y == 'yv' ) THEN id_y = id_dim_yv_xy(av) ENDIF IF ( grid_z == 'zu' ) THEN id_z = id_dim_zu_xy(av) ELSEIF ( grid_z == 'zw' ) THEN id_z = id_dim_zw_xy(av) ENDIF ! !-- Define the grid nc_stat = NF90_DEF_VAR( id_set_xy(av), do2d(av,i), & nc_precision(1), & (/ id_x, id_y, id_z, id_dim_time_xy(av) /), & id_var_do2d(av,i) ) IF ( .NOT. found ) THEN WRITE ( message_string, * ) 'no grid defined for', & ' variable ', TRIM( do2d(av,i) ) CALL message( 'define_netcdf_header', 'PA0244',& 0, 1, 0, 6, 0 ) ENDIF var_list = TRIM( var_list ) // TRIM( do2d(av,i) ) // ';' ENDIF CALL handle_netcdf_error( 'netcdf', 119 ) ! !-- Store the 'real' name of the variable (with *, for example) !-- in the long_name attribute. This is evaluated by Ferret, !-- for example. nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_do2d(av,i), & 'long_name', do2d(av,i) ) CALL handle_netcdf_error( 'netcdf', 120 ) ! !-- Define the variable's unit nc_stat = NF90_PUT_ATT( id_set_xy(av), id_var_do2d(av,i), & 'units', TRIM( do2d_unit(av,i) ) ) CALL handle_netcdf_error( 'netcdf', 354 ) #if defined( __netcdf4 ) ! !-- Set collective io operations for parallel io IF ( netcdf_data_format > 2 ) THEN nc_stat = NF90_VAR_PAR_ACCESS( id_set_xy(av), & id_var_do2d(av,i), & NF90_COLLECTIVE ) CALL handle_netcdf_error( 'netcdf', 448 ) ENDIF #endif ENDIF i = i + 1 ENDDO ! !-- No arrays to output. Close the netcdf file and return. IF ( i == 1 ) RETURN ! !-- Write the list of variables as global attribute (this is used by !-- restart runs and by combine_plot_fields) nc_stat = NF90_PUT_ATT( id_set_xy(av), NF90_GLOBAL, 'VAR_LIST', & var_list ) CALL handle_netcdf_error( 'netcdf', 121 ) ! !-- Leave NetCDF define mode nc_stat = NF90_ENDDEF( id_set_xy(av) ) CALL handle_netcdf_error( 'netcdf', 122 ) ! !-- Write axis data: z_xy, x, y ALLOCATE( netcdf_data(1:ns) ) ! !-- Write zu data DO i = 1, ns IF( section(i,1) == -1 ) THEN netcdf_data(i) = -1.0 ! section averaged along z ELSE netcdf_data(i) = zu( section(i,1) ) ENDIF ENDDO nc_stat = NF90_PUT_VAR( id_set_xy(av), id_var_zu_xy(av), & netcdf_data, start = (/ 1 /), & count = (/ ns /) ) CALL handle_netcdf_error( 'netcdf', 123 ) ! !-- Write zw data DO i = 1, ns IF( section(i,1) == -1 ) THEN netcdf_data(i) = -1.0 ! section averaged along z ELSE netcdf_data(i) = zw( section(i,1) ) ENDIF ENDDO nc_stat = NF90_PUT_VAR( id_set_xy(av), id_var_zw_xy(av), & netcdf_data, start = (/ 1 /), & count = (/ ns /) ) CALL handle_netcdf_error( 'netcdf', 124 ) ! !-- Write gridpoint number data netcdf_data(1:ns) = section(1:ns,1) nc_stat = NF90_PUT_VAR( id_set_xy(av), id_var_ind_z_xy(av), & netcdf_data, start = (/ 1 /), & count = (/ ns /) ) CALL handle_netcdf_error( 'netcdf', 125 ) DEALLOCATE( netcdf_data ) ! !-- Write the cross section height u*, t* nc_stat = NF90_PUT_VAR( id_set_xy(av), id_var_zu1_xy(av), & (/ zu(nzb+1) /), start = (/ 1 /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'netcdf', 126 ) ! !-- Write data for x (shifted by +dx/2) and xu axis ALLOCATE( netcdf_data(0:nx+1) ) DO i = 0, nx+1 netcdf_data(i) = ( i + 0.5 ) * dx ENDDO nc_stat = NF90_PUT_VAR( id_set_xy(av), id_var_x_xy(av), netcdf_data, & start = (/ 1 /), count = (/ nx+2 /) ) CALL handle_netcdf_error( 'netcdf', 127 ) DO i = 0, nx+1 netcdf_data(i) = i * dx ENDDO nc_stat = NF90_PUT_VAR( id_set_xy(av), id_var_xu_xy(av), & netcdf_data, start = (/ 1 /), & count = (/ nx+2 /) ) CALL handle_netcdf_error( 'netcdf', 367 ) DEALLOCATE( netcdf_data ) ! !-- Write data for y (shifted by +dy/2) and yv axis ALLOCATE( netcdf_data(0:ny+1) ) DO i = 0, ny+1 netcdf_data(i) = ( i + 0.5 ) * dy ENDDO nc_stat = NF90_PUT_VAR( id_set_xy(av), id_var_y_xy(av), netcdf_data, & start = (/ 1 /), count = (/ ny+2 /)) CALL handle_netcdf_error( 'netcdf', 128 ) DO i = 0, ny+1 netcdf_data(i) = i * dy ENDDO nc_stat = NF90_PUT_VAR( id_set_xy(av), id_var_yv_xy(av), & netcdf_data, start = (/ 1 /), & count = (/ ny+2 /)) CALL handle_netcdf_error( 'netcdf', 368 ) DEALLOCATE( netcdf_data ) ! !-- In case of non-flat topography write height information IF ( TRIM( topography ) /= 'flat' ) THEN nc_stat = NF90_PUT_VAR( id_set_xy(av), id_var_zusi_xy(av), & zu_s_inner(0:nx+1,0:ny+1), & start = (/ 1, 1 /), & count = (/ nx+2, ny+2 /) ) CALL handle_netcdf_error( 'netcdf', 427 ) nc_stat = NF90_PUT_VAR( id_set_xy(av), id_var_zwwi_xy(av), & zw_w_inner(0:nx+1,0:ny+1), & start = (/ 1, 1 /), & count = (/ nx+2, ny+2 /) ) CALL handle_netcdf_error( 'netcdf', 428 ) ENDIF CASE ( 'xy_ext' ) ! !-- Get the list of variables and compare with the actual run. !-- First var_list_old has to be reset, since GET_ATT does not assign !-- trailing blanks. var_list_old = ' ' nc_stat = NF90_GET_ATT( id_set_xy(av), NF90_GLOBAL, 'VAR_LIST', & var_list_old ) CALL handle_netcdf_error( 'netcdf', 129 ) var_list = ';' i = 1 DO WHILE ( do2d(av,i)(1:1) /= ' ' ) IF ( INDEX( do2d(av,i), 'xy' ) /= 0 ) THEN var_list = TRIM( var_list ) // TRIM( do2d(av,i) ) // ';' ENDIF i = i + 1 ENDDO IF ( av == 0 ) THEN var = '(xy)' ELSE var = '(xy_av)' ENDIF IF ( TRIM( var_list ) /= TRIM( var_list_old ) ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '& but this file cannot be extended due to' // & ' variable mismatch.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0249', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Calculate the number of current sections ns = 1 DO WHILE ( section(ns,1) /= -9999 .AND. ns <= 100 ) ns = ns + 1 ENDDO ns = ns - 1 ! !-- Get and compare the number of horizontal cross sections nc_stat = NF90_INQ_VARID( id_set_xy(av), 'zu_xy', id_var_zu_xy(av) ) CALL handle_netcdf_error( 'netcdf', 130 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_xy(av), id_var_zu_xy(av), & dimids = id_dim_zu_xy_old ) CALL handle_netcdf_error( 'netcdf', 131 ) id_dim_zu_xy(av) = id_dim_zu_xy_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_xy(av), id_dim_zu_xy(av), & len = ns_old ) CALL handle_netcdf_error( 'netcdf', 132 ) IF ( ns /= ns_old ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended due to' // & ' mismatch in number of' // & '&cross sections.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0250', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Get and compare the heights of the cross sections ALLOCATE( netcdf_data(1:ns_old) ) nc_stat = NF90_GET_VAR( id_set_xy(av), id_var_zu_xy(av), netcdf_data ) CALL handle_netcdf_error( 'netcdf', 133 ) DO i = 1, ns IF ( section(i,1) /= -1 ) THEN IF ( zu(section(i,1)) /= netcdf_data(i) ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended' // & ' due to mismatch in cross' // & '§ion levels.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0251', & 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ELSE IF ( -1.0 /= netcdf_data(i) ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended' // & ' due to mismatch in cross' // & '§ion levels.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0251',& 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ENDIF ENDDO DEALLOCATE( netcdf_data ) ! !-- Get the id of the time coordinate (unlimited coordinate) and its !-- last index on the file. The next time level is do2d..count+1. !-- The current time must be larger than the last output time !-- on the file. nc_stat = NF90_INQ_VARID( id_set_xy(av), 'time', id_var_time_xy(av) ) CALL handle_netcdf_error( 'netcdf', 134 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_xy(av), id_var_time_xy(av), & dimids = id_dim_time_old ) CALL handle_netcdf_error( 'netcdf', 135 ) id_dim_time_xy(av) = id_dim_time_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_xy(av), id_dim_time_xy(av), & len = do2d_xy_time_count(av) ) CALL handle_netcdf_error( 'netcdf', 136 ) nc_stat = NF90_GET_VAR( id_set_xy(av), id_var_time_xy(av), & last_time_coordinate, & start = (/ do2d_xy_time_count(av) /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'netcdf', 137 ) IF ( last_time_coordinate(1) >= simulated_time ) THEN message_string = 'NetCDF file for cross sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended becaus' // & 'e the current output time' // & '&is less or equal than the last output t' // & 'ime on this file.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0252', 0, 1, 0, 6, 0 ) do2d_xy_time_count(av) = 0 extend = .FALSE. RETURN ENDIF ! !-- Dataset seems to be extendable. !-- Now get the variable ids. i = 1 DO WHILE ( do2d(av,i)(1:1) /= ' ' ) IF ( INDEX( do2d(av,i), 'xy' ) /= 0 ) THEN nc_stat = NF90_INQ_VARID( id_set_xy(av), do2d(av,i), & id_var_do2d(av,i) ) CALL handle_netcdf_error( 'netcdf', 138 ) #if defined( __netcdf4 ) ! !-- Set collective io operations for parallel io IF ( netcdf_data_format > 2 ) THEN nc_stat = NF90_VAR_PAR_ACCESS( id_set_xy(av), & id_var_do2d(av,i), & NF90_COLLECTIVE ) CALL handle_netcdf_error( 'netcdf', 454 ) ENDIF #endif ENDIF i = i + 1 ENDDO ! !-- Update the title attribute on file !-- In order to avoid 'data mode' errors if updated attributes are larger !-- than their original size, NF90_PUT_ATT is called in 'define mode' !-- enclosed by NF90_REDEF and NF90_ENDDEF calls. This implies a possible !-- performance loss due to data copying; an alternative strategy would be !-- to ensure equal attribute size in a job chain. Maybe revise later. IF ( av == 0 ) THEN time_average_text = ' ' ELSE WRITE (time_average_text, '('', '',F7.1,'' s average'')') & averaging_interval ENDIF nc_stat = NF90_REDEF( id_set_xy(av) ) CALL handle_netcdf_error( 'netcdf', 431 ) nc_stat = NF90_PUT_ATT( id_set_xy(av), NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 139 ) nc_stat = NF90_ENDDEF( id_set_xy(av) ) CALL handle_netcdf_error( 'netcdf', 432 ) message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found.' // & '&This file will be extended.' CALL message( 'define_netcdf_header', 'PA0253', 0, 0, 0, 6, 0 ) CASE ( 'xz_new' ) ! !-- Define some global attributes of the dataset nc_stat = NF90_PUT_ATT( id_set_xz(av), NF90_GLOBAL, 'Conventions', & 'COARDS' ) CALL handle_netcdf_error( 'netcdf', 140 ) IF ( av == 0 ) THEN time_average_text = ' ' ELSE WRITE (time_average_text, '('', '',F7.1,'' s average'')') & averaging_interval ENDIF nc_stat = NF90_PUT_ATT( id_set_xz(av), NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 141 ) IF ( av == 1 ) THEN WRITE ( time_average_text,'(F7.1,'' s avg'')' ) averaging_interval nc_stat = NF90_PUT_ATT( id_set_xz(av), NF90_GLOBAL, 'time_avg', & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 141 ) ENDIF ! !-- Define time coordinate for xz sections (unlimited dimension) nc_stat = NF90_DEF_DIM( id_set_xz(av), 'time', NF90_UNLIMITED, & id_dim_time_xz(av) ) CALL handle_netcdf_error( 'netcdf', 142 ) nc_stat = NF90_DEF_VAR( id_set_xz(av), 'time', NF90_DOUBLE, & id_dim_time_xz(av), id_var_time_xz(av) ) CALL handle_netcdf_error( 'netcdf', 143 ) nc_stat = NF90_PUT_ATT( id_set_xz(av), id_var_time_xz(av), 'units', & 'seconds') CALL handle_netcdf_error( 'netcdf', 144 ) ! !-- Define the spatial dimensions and coordinates for xz-sections. !-- First, determine the number of vertical sections to be written. IF ( section(1,2) == -9999 ) THEN RETURN ELSE ns = 1 DO WHILE ( section(ns,2) /= -9999 .AND. ns <= 100 ) ns = ns + 1 ENDDO ns = ns - 1 ENDIF ! !-- Define y-axis (for scalar position) nc_stat = NF90_DEF_DIM( id_set_xz(av), 'y_xz', ns, id_dim_y_xz(av) ) CALL handle_netcdf_error( 'netcdf', 145 ) nc_stat = NF90_DEF_VAR( id_set_xz(av), 'y_xz', NF90_DOUBLE, & id_dim_y_xz(av), id_var_y_xz(av) ) CALL handle_netcdf_error( 'netcdf', 146 ) nc_stat = NF90_PUT_ATT( id_set_xz(av), id_var_y_xz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 147 ) ! !-- Define y-axis (for v position) nc_stat = NF90_DEF_DIM( id_set_xz(av), 'yv_xz', ns, id_dim_yv_xz(av) ) CALL handle_netcdf_error( 'netcdf', 369 ) nc_stat = NF90_DEF_VAR( id_set_xz(av), 'yv_xz', NF90_DOUBLE, & id_dim_yv_xz(av), id_var_yv_xz(av) ) CALL handle_netcdf_error( 'netcdf', 370 ) nc_stat = NF90_PUT_ATT( id_set_xz(av), id_var_yv_xz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 371 ) ! !-- Define a variable to store the layer indices of the vertical cross !-- sections nc_stat = NF90_DEF_VAR( id_set_xz(av), 'ind_y_xz', NF90_DOUBLE, & id_dim_y_xz(av), id_var_ind_y_xz(av) ) CALL handle_netcdf_error( 'netcdf', 148 ) nc_stat = NF90_PUT_ATT( id_set_xz(av), id_var_ind_y_xz(av), 'units', & 'gridpoints') CALL handle_netcdf_error( 'netcdf', 149 ) ! !-- Define x-axis (for scalar position) nc_stat = NF90_DEF_DIM( id_set_xz(av), 'x', nx+2, id_dim_x_xz(av) ) CALL handle_netcdf_error( 'netcdf', 150 ) nc_stat = NF90_DEF_VAR( id_set_xz(av), 'x', NF90_DOUBLE, & id_dim_x_xz(av), id_var_x_xz(av) ) CALL handle_netcdf_error( 'netcdf', 151 ) nc_stat = NF90_PUT_ATT( id_set_xz(av), id_var_x_xz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 152 ) ! !-- Define x-axis (for u position) nc_stat = NF90_DEF_DIM( id_set_xz(av), 'xu', nx+2, id_dim_xu_xz(av) ) CALL handle_netcdf_error( 'netcdf', 372 ) nc_stat = NF90_DEF_VAR( id_set_xz(av), 'xu', NF90_DOUBLE, & id_dim_xu_xz(av), id_var_xu_xz(av) ) CALL handle_netcdf_error( 'netcdf', 373 ) nc_stat = NF90_PUT_ATT( id_set_xz(av), id_var_xu_xz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 374 ) ! !-- Define the two z-axes (zu and zw) nc_stat = NF90_DEF_DIM( id_set_xz(av), 'zu', nz+2, id_dim_zu_xz(av) ) CALL handle_netcdf_error( 'netcdf', 153 ) nc_stat = NF90_DEF_VAR( id_set_xz(av), 'zu', NF90_DOUBLE, & id_dim_zu_xz(av), id_var_zu_xz(av) ) CALL handle_netcdf_error( 'netcdf', 154 ) nc_stat = NF90_PUT_ATT( id_set_xz(av), id_var_zu_xz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 155 ) nc_stat = NF90_DEF_DIM( id_set_xz(av), 'zw', nz+2, id_dim_zw_xz(av) ) CALL handle_netcdf_error( 'netcdf', 156 ) nc_stat = NF90_DEF_VAR( id_set_xz(av), 'zw', NF90_DOUBLE, & id_dim_zw_xz(av), id_var_zw_xz(av) ) CALL handle_netcdf_error( 'netcdf', 157 ) nc_stat = NF90_PUT_ATT( id_set_xz(av), id_var_zw_xz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 158 ) ! !-- Define the variables var_list = ';' i = 1 DO WHILE ( do2d(av,i)(1:1) /= ' ' ) IF ( INDEX( do2d(av,i), 'xz' ) /= 0 ) THEN ! !-- Check for the grid found = .TRUE. SELECT CASE ( do2d(av,i) ) ! !-- Most variables are defined on the zu grid CASE ( 'e_xz', 'p_xz', 'pc_xz', 'pr_xz', 'pt_xz', 'q_xz', & 'ql_xz', 'ql_c_xz', 'ql_v_xz', 'ql_vp_xz', 'qv_xz', & 'rho_xz', 's_xz', 'sa_xz', 'vpt_xz' ) grid_x = 'x' grid_y = 'y' grid_z = 'zu' ! !-- u grid CASE ( 'u_xz' ) grid_x = 'xu' grid_y = 'y' grid_z = 'zu' ! !-- v grid CASE ( 'v_xz' ) grid_x = 'x' grid_y = 'yv' grid_z = 'zu' ! !-- w grid CASE ( 'w_xz' ) grid_x = 'x' grid_y = 'y' grid_z = 'zw' CASE DEFAULT ! !-- Check for user-defined quantities CALL user_define_netcdf_grid( do2d(av,i), found, & grid_x, grid_y, grid_z ) END SELECT ! !-- Select the respective dimension ids IF ( grid_x == 'x' ) THEN id_x = id_dim_x_xz(av) ELSEIF ( grid_x == 'xu' ) THEN id_x = id_dim_xu_xz(av) ENDIF IF ( grid_y == 'y' ) THEN id_y = id_dim_y_xz(av) ELSEIF ( grid_y == 'yv' ) THEN id_y = id_dim_yv_xz(av) ENDIF IF ( grid_z == 'zu' ) THEN id_z = id_dim_zu_xz(av) ELSEIF ( grid_z == 'zw' ) THEN id_z = id_dim_zw_xz(av) ENDIF ! !-- Define the grid nc_stat = NF90_DEF_VAR( id_set_xz(av), do2d(av,i), & nc_precision(2), & (/ id_x, id_y, id_z, id_dim_time_xz(av) /), & id_var_do2d(av,i) ) IF ( .NOT. found ) THEN WRITE ( message_string, * ) 'no grid defined for', & ' variable ', TRIM( do2d(av,i) ) CALL message( 'define_netcdf_header', 'PA0244',& 0, 1, 0, 6, 0 ) ENDIF var_list = TRIM( var_list ) // TRIM( do2d(av,i) ) // ';' CALL handle_netcdf_error( 'netcdf', 159 ) ! !-- Store the 'real' name of the variable (with *, for example) !-- in the long_name attribute. This is evaluated by Ferret, !-- for example. nc_stat = NF90_PUT_ATT( id_set_xz(av), id_var_do2d(av,i), & 'long_name', do2d(av,i) ) CALL handle_netcdf_error( 'netcdf', 160 ) ! !-- Define the variable's unit nc_stat = NF90_PUT_ATT( id_set_xz(av), id_var_do2d(av,i), & 'units', TRIM( do2d_unit(av,i) ) ) CALL handle_netcdf_error( 'netcdf', 355 ) #if defined( __netcdf4 ) ! !-- Set independent io operations for parallel io. Collective io !-- is only allowed in case of a 1d-decomposition along x, because !-- otherwise, not all PEs have output data. IF ( netcdf_data_format > 2 ) THEN IF ( npey == 1 ) THEN nc_stat = NF90_VAR_PAR_ACCESS( id_set_xz(av), & id_var_do2d(av,i), & NF90_COLLECTIVE ) ELSE ! !-- ATTENTION: Due to a probable bug in the NetCDF4 !-- installation, independet output does not work !-- A workaround is used in data_output_2d on those !-- PEs having no data nc_stat = NF90_VAR_PAR_ACCESS( id_set_xz(av), & id_var_do2d(av,i), & NF90_COLLECTIVE ) ! nc_stat = NF90_VAR_PAR_ACCESS( id_set_xz(av), & ! id_var_do2d(av,i), & ! NF90_INDEPENDENT ) ENDIF CALL handle_netcdf_error( 'netcdf', 449 ) ENDIF #endif ENDIF i = i + 1 ENDDO ! !-- No arrays to output. Close the netcdf file and return. IF ( i == 1 ) RETURN ! !-- Write the list of variables as global attribute (this is used by !-- restart runs and by combine_plot_fields) nc_stat = NF90_PUT_ATT( id_set_xz(av), NF90_GLOBAL, 'VAR_LIST', & var_list ) CALL handle_netcdf_error( 'netcdf', 161 ) ! !-- Leave NetCDF define mode nc_stat = NF90_ENDDEF( id_set_xz(av) ) CALL handle_netcdf_error( 'netcdf', 162 ) ! !-- Write axis data: y_xz, x, zu, zw ALLOCATE( netcdf_data(1:ns) ) ! !-- Write y_xz data (shifted by +dy/2) DO i = 1, ns IF( section(i,2) == -1 ) THEN netcdf_data(i) = -1.0 ! section averaged along y ELSE netcdf_data(i) = ( section(i,2) + 0.5 ) * dy ENDIF ENDDO nc_stat = NF90_PUT_VAR( id_set_xz(av), id_var_y_xz(av), netcdf_data, & start = (/ 1 /), count = (/ ns /) ) CALL handle_netcdf_error( 'netcdf', 163 ) ! !-- Write yv_xz data DO i = 1, ns IF( section(i,2) == -1 ) THEN netcdf_data(i) = -1.0 ! section averaged along y ELSE netcdf_data(i) = section(i,2) * dy ENDIF ENDDO nc_stat = NF90_PUT_VAR( id_set_xz(av), id_var_yv_xz(av), & netcdf_data, start = (/ 1 /), & count = (/ ns /) ) CALL handle_netcdf_error( 'netcdf', 375 ) ! !-- Write gridpoint number data netcdf_data(1:ns) = section(1:ns,2) nc_stat = NF90_PUT_VAR( id_set_xz(av), id_var_ind_y_xz(av), & netcdf_data, start = (/ 1 /), & count = (/ ns /) ) CALL handle_netcdf_error( 'netcdf', 164 ) DEALLOCATE( netcdf_data ) ! !-- Write data for x (shifted by +dx/2) and xu axis ALLOCATE( netcdf_data(0:nx+1) ) DO i = 0, nx+1 netcdf_data(i) = ( i + 0.5 ) * dx ENDDO nc_stat = NF90_PUT_VAR( id_set_xz(av), id_var_x_xz(av), netcdf_data, & start = (/ 1 /), count = (/ nx+2 /) ) CALL handle_netcdf_error( 'netcdf', 165 ) DO i = 0, nx+1 netcdf_data(i) = i * dx ENDDO nc_stat = NF90_PUT_VAR( id_set_xz(av), id_var_xu_xz(av), & netcdf_data, start = (/ 1 /), & count = (/ nx+2 /) ) CALL handle_netcdf_error( 'netcdf', 377 ) DEALLOCATE( netcdf_data ) ! !-- Write zu and zw data (vertical axes) ALLOCATE( netcdf_data(0:nz+1) ) netcdf_data(0:nz+1) = zu(nzb:nzt+1) nc_stat = NF90_PUT_VAR( id_set_xz(av), id_var_zu_xz(av), & netcdf_data, start = (/ 1 /), & count = (/ nz+2 /) ) CALL handle_netcdf_error( 'netcdf', 166 ) netcdf_data(0:nz+1) = zw(nzb:nzt+1) nc_stat = NF90_PUT_VAR( id_set_xz(av), id_var_zw_xz(av), & netcdf_data, start = (/ 1 /), & count = (/ nz+2 /) ) CALL handle_netcdf_error( 'netcdf', 167 ) DEALLOCATE( netcdf_data ) CASE ( 'xz_ext' ) ! !-- Get the list of variables and compare with the actual run. !-- First var_list_old has to be reset, since GET_ATT does not assign !-- trailing blanks. var_list_old = ' ' nc_stat = NF90_GET_ATT( id_set_xz(av), NF90_GLOBAL, 'VAR_LIST', & var_list_old ) CALL handle_netcdf_error( 'netcdf', 168 ) var_list = ';' i = 1 DO WHILE ( do2d(av,i)(1:1) /= ' ' ) IF ( INDEX( do2d(av,i), 'xz' ) /= 0 ) THEN var_list = TRIM( var_list ) // TRIM( do2d(av,i) ) // ';' ENDIF i = i + 1 ENDDO IF ( av == 0 ) THEN var = '(xz)' ELSE var = '(xz_av)' ENDIF IF ( TRIM( var_list ) /= TRIM( var_list_old ) ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '& but this file cannot be extended due to' // & ' variable mismatch.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0249', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Calculate the number of current sections ns = 1 DO WHILE ( section(ns,2) /= -9999 .AND. ns <= 100 ) ns = ns + 1 ENDDO ns = ns - 1 ! !-- Get and compare the number of vertical cross sections nc_stat = NF90_INQ_VARID( id_set_xz(av), 'y_xz', id_var_y_xz(av) ) CALL handle_netcdf_error( 'netcdf', 169 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_xz(av), id_var_y_xz(av), & dimids = id_dim_y_xz_old ) CALL handle_netcdf_error( 'netcdf', 170 ) id_dim_y_xz(av) = id_dim_y_xz_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_xz(av), id_dim_y_xz(av), & len = ns_old ) CALL handle_netcdf_error( 'netcdf', 171 ) IF ( ns /= ns_old ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended due to' // & ' mismatch in number of' // & '&cross sections.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0250', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Get and compare the heights of the cross sections ALLOCATE( netcdf_data(1:ns_old) ) nc_stat = NF90_GET_VAR( id_set_xz(av), id_var_y_xz(av), netcdf_data ) CALL handle_netcdf_error( 'netcdf', 172 ) DO i = 1, ns IF ( section(i,2) /= -1 ) THEN IF ( ( ( section(i,2) + 0.5 ) * dy ) /= netcdf_data(i) ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended' // & ' due to mismatch in cross' // & '§ion levels.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0251', & 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ELSE IF ( -1.0 /= netcdf_data(i) ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended' // & ' due to mismatch in cross' // & '§ion levels.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0251', & 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ENDIF ENDDO DEALLOCATE( netcdf_data ) ! !-- Get the id of the time coordinate (unlimited coordinate) and its !-- last index on the file. The next time level is do2d..count+1. !-- The current time must be larger than the last output time !-- on the file. nc_stat = NF90_INQ_VARID( id_set_xz(av), 'time', id_var_time_xz(av) ) CALL handle_netcdf_error( 'netcdf', 173 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_xz(av), id_var_time_xz(av), & dimids = id_dim_time_old ) CALL handle_netcdf_error( 'netcdf', 174 ) id_dim_time_xz(av) = id_dim_time_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_xz(av), id_dim_time_xz(av), & len = do2d_xz_time_count(av) ) CALL handle_netcdf_error( 'netcdf', 175 ) nc_stat = NF90_GET_VAR( id_set_xz(av), id_var_time_xz(av), & last_time_coordinate, & start = (/ do2d_xz_time_count(av) /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'netcdf', 176 ) IF ( last_time_coordinate(1) >= simulated_time ) THEN message_string = 'NetCDF file for cross sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended becaus' // & 'e the current output time' // & '&is less or equal than the last output t' // & 'ime on this file.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0252', 0, 1, 0, 6, 0 ) do2d_xz_time_count(av) = 0 extend = .FALSE. RETURN ENDIF ! !-- Dataset seems to be extendable. !-- Now get the variable ids. i = 1 DO WHILE ( do2d(av,i)(1:1) /= ' ' ) IF ( INDEX( do2d(av,i), 'xz' ) /= 0 ) THEN nc_stat = NF90_INQ_VARID( id_set_xz(av), do2d(av,i), & id_var_do2d(av,i) ) CALL handle_netcdf_error( 'netcdf', 177 ) #if defined( __netcdf4 ) ! !-- Set independent io operations for parallel io. Collective io !-- is only allowed in case of a 1d-decomposition along x, because !-- otherwise, not all PEs have output data. IF ( netcdf_data_format > 2 ) THEN IF ( npey == 1 ) THEN nc_stat = NF90_VAR_PAR_ACCESS( id_set_xz(av), & id_var_do2d(av,i), & NF90_COLLECTIVE ) ELSE ! !-- ATTENTION: Due to a probable bug in the NetCDF4 !-- installation, independet output does not work !-- A workaround is used in data_output_2d on those !-- PEs having no data nc_stat = NF90_VAR_PAR_ACCESS( id_set_xz(av), & id_var_do2d(av,i), & NF90_COLLECTIVE ) ! nc_stat = NF90_VAR_PAR_ACCESS( id_set_xz(av), & ! id_var_do2d(av,i), & ! NF90_INDEPENDENT ) ENDIF CALL handle_netcdf_error( 'netcdf', 455 ) ENDIF #endif ENDIF i = i + 1 ENDDO ! !-- Update the title attribute on file !-- In order to avoid 'data mode' errors if updated attributes are larger !-- than their original size, NF90_PUT_ATT is called in 'define mode' !-- enclosed by NF90_REDEF and NF90_ENDDEF calls. This implies a possible !-- performance loss due to data copying; an alternative strategy would be !-- to ensure equal attribute size in a job chain. Maybe revise later. IF ( av == 0 ) THEN time_average_text = ' ' ELSE WRITE (time_average_text, '('', '',F7.1,'' s average'')') & averaging_interval ENDIF nc_stat = NF90_REDEF( id_set_xz(av) ) CALL handle_netcdf_error( 'netcdf', 433 ) nc_stat = NF90_PUT_ATT( id_set_xz(av), NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 178 ) nc_stat = NF90_ENDDEF( id_set_xz(av) ) CALL handle_netcdf_error( 'netcdf', 434 ) message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found.' // & '&This file will be extended.' CALL message( 'define_netcdf_header', 'PA0253', 0, 0, 0, 6, 0 ) CASE ( 'yz_new' ) ! !-- Define some global attributes of the dataset nc_stat = NF90_PUT_ATT( id_set_yz(av), NF90_GLOBAL, 'Conventions', & 'COARDS' ) CALL handle_netcdf_error( 'netcdf', 179 ) IF ( av == 0 ) THEN time_average_text = ' ' ELSE WRITE (time_average_text, '('', '',F7.1,'' s average'')') & averaging_interval ENDIF nc_stat = NF90_PUT_ATT( id_set_yz(av), NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 180 ) IF ( av == 1 ) THEN WRITE ( time_average_text,'(F7.1,'' s avg'')' ) averaging_interval nc_stat = NF90_PUT_ATT( id_set_yz(av), NF90_GLOBAL, 'time_avg', & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 180 ) ENDIF ! !-- Define time coordinate for yz sections (unlimited dimension) nc_stat = NF90_DEF_DIM( id_set_yz(av), 'time', NF90_UNLIMITED, & id_dim_time_yz(av) ) CALL handle_netcdf_error( 'netcdf', 181 ) nc_stat = NF90_DEF_VAR( id_set_yz(av), 'time', NF90_DOUBLE, & id_dim_time_yz(av), id_var_time_yz(av) ) CALL handle_netcdf_error( 'netcdf', 182 ) nc_stat = NF90_PUT_ATT( id_set_yz(av), id_var_time_yz(av), 'units', & 'seconds') CALL handle_netcdf_error( 'netcdf', 183 ) ! !-- Define the spatial dimensions and coordinates for yz-sections. !-- First, determine the number of vertical sections to be written. IF ( section(1,3) == -9999 ) THEN RETURN ELSE ns = 1 DO WHILE ( section(ns,3) /= -9999 .AND. ns <= 100 ) ns = ns + 1 ENDDO ns = ns - 1 ENDIF ! !-- Define x axis (for scalar position) nc_stat = NF90_DEF_DIM( id_set_yz(av), 'x_yz', ns, id_dim_x_yz(av) ) CALL handle_netcdf_error( 'netcdf', 184 ) nc_stat = NF90_DEF_VAR( id_set_yz(av), 'x_yz', NF90_DOUBLE, & id_dim_x_yz(av), id_var_x_yz(av) ) CALL handle_netcdf_error( 'netcdf', 185 ) nc_stat = NF90_PUT_ATT( id_set_yz(av), id_var_x_yz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 186 ) ! !-- Define x axis (for u position) nc_stat = NF90_DEF_DIM( id_set_yz(av), 'xu_yz', ns, id_dim_xu_yz(av) ) CALL handle_netcdf_error( 'netcdf', 377 ) nc_stat = NF90_DEF_VAR( id_set_yz(av), 'xu_yz', NF90_DOUBLE, & id_dim_xu_yz(av), id_var_xu_yz(av) ) CALL handle_netcdf_error( 'netcdf', 378 ) nc_stat = NF90_PUT_ATT( id_set_yz(av), id_var_xu_yz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 379 ) ! !-- Define a variable to store the layer indices of the vertical cross !-- sections nc_stat = NF90_DEF_VAR( id_set_yz(av), 'ind_x_yz', NF90_DOUBLE, & id_dim_x_yz(av), id_var_ind_x_yz(av) ) CALL handle_netcdf_error( 'netcdf', 187 ) nc_stat = NF90_PUT_ATT( id_set_yz(av), id_var_ind_x_yz(av), 'units', & 'gridpoints') CALL handle_netcdf_error( 'netcdf', 188 ) ! !-- Define y-axis (for scalar position) nc_stat = NF90_DEF_DIM( id_set_yz(av), 'y', ny+2, id_dim_y_yz(av) ) CALL handle_netcdf_error( 'netcdf', 189 ) nc_stat = NF90_DEF_VAR( id_set_yz(av), 'y', NF90_DOUBLE, & id_dim_y_yz(av), id_var_y_yz(av) ) CALL handle_netcdf_error( 'netcdf', 190 ) nc_stat = NF90_PUT_ATT( id_set_yz(av), id_var_y_yz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 191 ) ! !-- Define y-axis (for v position) nc_stat = NF90_DEF_DIM( id_set_yz(av), 'yv', ny+2, id_dim_yv_yz(av) ) CALL handle_netcdf_error( 'netcdf', 380 ) nc_stat = NF90_DEF_VAR( id_set_yz(av), 'yv', NF90_DOUBLE, & id_dim_yv_yz(av), id_var_yv_yz(av) ) CALL handle_netcdf_error( 'netcdf', 381 ) nc_stat = NF90_PUT_ATT( id_set_yz(av), id_var_yv_yz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 382 ) ! !-- Define the two z-axes (zu and zw) nc_stat = NF90_DEF_DIM( id_set_yz(av), 'zu', nz+2, id_dim_zu_yz(av) ) CALL handle_netcdf_error( 'netcdf', 192 ) nc_stat = NF90_DEF_VAR( id_set_yz(av), 'zu', NF90_DOUBLE, & id_dim_zu_yz(av), id_var_zu_yz(av) ) CALL handle_netcdf_error( 'netcdf', 193 ) nc_stat = NF90_PUT_ATT( id_set_yz(av), id_var_zu_yz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 194 ) nc_stat = NF90_DEF_DIM( id_set_yz(av), 'zw', nz+2, id_dim_zw_yz(av) ) CALL handle_netcdf_error( 'netcdf', 195 ) nc_stat = NF90_DEF_VAR( id_set_yz(av), 'zw', NF90_DOUBLE, & id_dim_zw_yz(av), id_var_zw_yz(av) ) CALL handle_netcdf_error( 'netcdf', 196 ) nc_stat = NF90_PUT_ATT( id_set_yz(av), id_var_zw_yz(av), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 197 ) ! !-- Define the variables var_list = ';' i = 1 DO WHILE ( do2d(av,i)(1:1) /= ' ' ) IF ( INDEX( do2d(av,i), 'yz' ) /= 0 ) THEN ! !-- Check for the grid found = .TRUE. SELECT CASE ( do2d(av,i) ) ! !-- Most variables are defined on the zu grid CASE ( 'e_yz', 'p_yz', 'pc_yz', 'pr_yz', 'pt_yz', 'q_yz', & 'ql_yz', 'ql_c_yz', 'ql_v_yz', 'ql_vp_yz', 'qv_yz', & 'rho_yz', 's_yz', 'sa_yz', 'vpt_yz' ) grid_x = 'x' grid_y = 'y' grid_z = 'zu' ! !-- u grid CASE ( 'u_yz' ) grid_x = 'xu' grid_y = 'y' grid_z = 'zu' ! !-- v grid CASE ( 'v_yz' ) grid_x = 'x' grid_y = 'yv' grid_z = 'zu' ! !-- w grid CASE ( 'w_yz' ) grid_x = 'x' grid_y = 'y' grid_z = 'zw' CASE DEFAULT ! !-- Check for user-defined quantities CALL user_define_netcdf_grid( do2d(av,i), found, & grid_x, grid_y, grid_z ) END SELECT ! !-- Select the respective dimension ids IF ( grid_x == 'x' ) THEN id_x = id_dim_x_yz(av) ELSEIF ( grid_x == 'xu' ) THEN id_x = id_dim_xu_yz(av) ENDIF IF ( grid_y == 'y' ) THEN id_y = id_dim_y_yz(av) ELSEIF ( grid_y == 'yv' ) THEN id_y = id_dim_yv_yz(av) ENDIF IF ( grid_z == 'zu' ) THEN id_z = id_dim_zu_yz(av) ELSEIF ( grid_z == 'zw' ) THEN id_z = id_dim_zw_yz(av) ENDIF ! !-- Define the grid nc_stat = NF90_DEF_VAR( id_set_yz(av), do2d(av,i), & nc_precision(3), & (/ id_x, id_y, id_z, id_dim_time_yz(av) /), & id_var_do2d(av,i) ) IF ( .NOT. found ) THEN WRITE ( message_string, * ) 'no grid defined for', & ' variable ', TRIM( do2d(av,i) ) CALL message( 'define_netcdf_header', 'PA0244',& 0, 1, 0, 6, 0 ) ENDIF var_list = TRIM( var_list ) // TRIM( do2d(av,i) ) // ';' CALL handle_netcdf_error( 'netcdf', 198 ) ! !-- Store the 'real' name of the variable (with *, for example) !-- in the long_name attribute. This is evaluated by Ferret, !-- for example. nc_stat = NF90_PUT_ATT( id_set_yz(av), id_var_do2d(av,i), & 'long_name', do2d(av,i) ) CALL handle_netcdf_error( 'netcdf', 199 ) ! !-- Define the variable's unit nc_stat = NF90_PUT_ATT( id_set_yz(av), id_var_do2d(av,i), & 'units', TRIM( do2d_unit(av,i) ) ) CALL handle_netcdf_error( 'netcdf', 356 ) #if defined( __netcdf4 ) ! !-- Set independent io operations for parallel io. Collective io !-- is only allowed in case of a 1d-decomposition along y, because !-- otherwise, not all PEs have output data. IF ( netcdf_data_format > 2 ) THEN IF ( npex == 1 ) THEN nc_stat = NF90_VAR_PAR_ACCESS( id_set_yz(av), & id_var_do2d(av,i), & NF90_COLLECTIVE ) ELSE ! !-- ATTENTION: Due to a probable bug in the NetCDF4 !-- installation, independet output does not work !-- A workaround is used in data_output_2d on those !-- PEs having no data nc_stat = NF90_VAR_PAR_ACCESS( id_set_yz(av), & id_var_do2d(av,i), & NF90_COLLECTIVE ) ! nc_stat = NF90_VAR_PAR_ACCESS( id_set_yz(av), & ! id_var_do2d(av,i), & ! NF90_INDEPENDENT ) ENDIF CALL handle_netcdf_error( 'netcdf', 450 ) ENDIF #endif ENDIF i = i + 1 ENDDO ! !-- No arrays to output. Close the netcdf file and return. IF ( i == 1 ) RETURN ! !-- Write the list of variables as global attribute (this is used by !-- restart runs and by combine_plot_fields) nc_stat = NF90_PUT_ATT( id_set_yz(av), NF90_GLOBAL, 'VAR_LIST', & var_list ) CALL handle_netcdf_error( 'netcdf', 200 ) ! !-- Leave NetCDF define mode nc_stat = NF90_ENDDEF( id_set_yz(av) ) CALL handle_netcdf_error( 'netcdf', 201 ) ! !-- Write axis data: x_yz, y, zu, zw ALLOCATE( netcdf_data(1:ns) ) ! !-- Write x_yz data (shifted by +dx/2) DO i = 1, ns IF( section(i,3) == -1 ) THEN netcdf_data(i) = -1.0 ! section averaged along x ELSE netcdf_data(i) = ( section(i,3) + 0.5 ) * dx ENDIF ENDDO nc_stat = NF90_PUT_VAR( id_set_yz(av), id_var_x_yz(av), netcdf_data, & start = (/ 1 /), count = (/ ns /) ) CALL handle_netcdf_error( 'netcdf', 202 ) ! !-- Write x_yz data (xu grid) DO i = 1, ns IF( section(i,3) == -1 ) THEN netcdf_data(i) = -1.0 ! section averaged along x ELSE netcdf_data(i) = section(i,3) * dx ENDIF ENDDO nc_stat = NF90_PUT_VAR( id_set_yz(av), id_var_xu_yz(av), netcdf_data, & start = (/ 1 /), count = (/ ns /) ) CALL handle_netcdf_error( 'netcdf', 383 ) ! !-- Write gridpoint number data netcdf_data(1:ns) = section(1:ns,3) nc_stat = NF90_PUT_VAR( id_set_yz(av), id_var_ind_x_yz(av), & netcdf_data, start = (/ 1 /), & count = (/ ns /) ) CALL handle_netcdf_error( 'netcdf', 203 ) DEALLOCATE( netcdf_data ) ! !-- Write data for y (shifted by +dy/2) and yv axis ALLOCATE( netcdf_data(0:ny+1) ) DO j = 0, ny+1 netcdf_data(j) = ( j + 0.5 ) * dy ENDDO nc_stat = NF90_PUT_VAR( id_set_yz(av), id_var_y_yz(av), netcdf_data, & start = (/ 1 /), count = (/ ny+2 /) ) CALL handle_netcdf_error( 'netcdf', 204 ) DO j = 0, ny+1 netcdf_data(j) = j * dy ENDDO nc_stat = NF90_PUT_VAR( id_set_yz(av), id_var_yv_yz(av), & netcdf_data, start = (/ 1 /), & count = (/ ny+2 /) ) CALL handle_netcdf_error( 'netcdf', 384 ) DEALLOCATE( netcdf_data ) ! !-- Write zu and zw data (vertical axes) ALLOCATE( netcdf_data(0:nz+1) ) netcdf_data(0:nz+1) = zu(nzb:nzt+1) nc_stat = NF90_PUT_VAR( id_set_yz(av), id_var_zu_yz(av), & netcdf_data, start = (/ 1 /), & count = (/ nz+2 /) ) CALL handle_netcdf_error( 'netcdf', 205 ) netcdf_data(0:nz+1) = zw(nzb:nzt+1) nc_stat = NF90_PUT_VAR( id_set_yz(av), id_var_zw_yz(av), & netcdf_data, start = (/ 1 /), & count = (/ nz+2 /) ) CALL handle_netcdf_error( 'netcdf', 206 ) DEALLOCATE( netcdf_data ) CASE ( 'yz_ext' ) ! !-- Get the list of variables and compare with the actual run. !-- First var_list_old has to be reset, since GET_ATT does not assign !-- trailing blanks. var_list_old = ' ' nc_stat = NF90_GET_ATT( id_set_yz(av), NF90_GLOBAL, 'VAR_LIST', & var_list_old ) CALL handle_netcdf_error( 'netcdf', 207 ) var_list = ';' i = 1 DO WHILE ( do2d(av,i)(1:1) /= ' ' ) IF ( INDEX( do2d(av,i), 'yz' ) /= 0 ) THEN var_list = TRIM( var_list ) // TRIM( do2d(av,i) ) // ';' ENDIF i = i + 1 ENDDO IF ( av == 0 ) THEN var = '(yz)' ELSE var = '(yz_av)' ENDIF IF ( TRIM( var_list ) /= TRIM( var_list_old ) ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '& but this file cannot be extended due to' // & ' variable mismatch.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0249', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Calculate the number of current sections ns = 1 DO WHILE ( section(ns,3) /= -9999 .AND. ns <= 100 ) ns = ns + 1 ENDDO ns = ns - 1 ! !-- Get and compare the number of vertical cross sections nc_stat = NF90_INQ_VARID( id_set_yz(av), 'x_yz', id_var_x_yz(av) ) CALL handle_netcdf_error( 'netcdf', 208 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_yz(av), id_var_x_yz(av), & dimids = id_dim_x_yz_old ) CALL handle_netcdf_error( 'netcdf', 209 ) id_dim_x_yz(av) = id_dim_x_yz_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_yz(av), id_dim_x_yz(av), & len = ns_old ) CALL handle_netcdf_error( 'netcdf', 210 ) IF ( ns /= ns_old ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended due to' // & ' mismatch in number of' // & '&cross sections.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0250', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Get and compare the heights of the cross sections ALLOCATE( netcdf_data(1:ns_old) ) nc_stat = NF90_GET_VAR( id_set_yz(av), id_var_x_yz(av), netcdf_data ) CALL handle_netcdf_error( 'netcdf', 211 ) DO i = 1, ns IF ( section(i,3) /= -1 ) THEN IF ( ( ( section(i,3) + 0.5 ) * dx ) /= netcdf_data(i) ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended' // & ' due to mismatch in cross' // & '§ion levels.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0251', & 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ELSE IF ( -1.0 /= netcdf_data(i) ) THEN message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended' // & ' due to mismatch in cross' // & '§ion levels.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0251', & 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ENDIF ENDDO DEALLOCATE( netcdf_data ) ! !-- Get the id of the time coordinate (unlimited coordinate) and its !-- last index on the file. The next time level is pl2d..count+1. !-- The current time must be larger than the last output time !-- on the file. nc_stat = NF90_INQ_VARID( id_set_yz(av), 'time', id_var_time_yz(av) ) CALL handle_netcdf_error( 'netcdf', 212 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_yz(av), id_var_time_yz(av), & dimids = id_dim_time_old ) CALL handle_netcdf_error( 'netcdf', 213 ) id_dim_time_yz(av) = id_dim_time_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_yz(av), id_dim_time_yz(av), & len = do2d_yz_time_count(av) ) CALL handle_netcdf_error( 'netcdf', 214 ) nc_stat = NF90_GET_VAR( id_set_yz(av), id_var_time_yz(av), & last_time_coordinate, & start = (/ do2d_yz_time_count(av) /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'netcdf', 215 ) IF ( last_time_coordinate(1) >= simulated_time ) THEN message_string = 'NetCDF file for cross sections ' // & TRIM( var ) // ' from previous run found,' // & '&but this file cannot be extended becaus' // & 'e the current output time' // & '&is less or equal than the last output t' // & 'ime on this file.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0252', 0, 1, 0, 6, 0 ) do2d_yz_time_count(av) = 0 extend = .FALSE. RETURN ENDIF ! !-- Dataset seems to be extendable. !-- Now get the variable ids. i = 1 DO WHILE ( do2d(av,i)(1:1) /= ' ' ) IF ( INDEX( do2d(av,i), 'yz' ) /= 0 ) THEN nc_stat = NF90_INQ_VARID( id_set_yz(av), do2d(av,i), & id_var_do2d(av,i) ) CALL handle_netcdf_error( 'netcdf', 216 ) #if defined( __netcdf4 ) ! !-- Set independent io operations for parallel io. Collective io !-- is only allowed in case of a 1d-decomposition along y, because !-- otherwise, not all PEs have output data. IF ( netcdf_data_format > 2 ) THEN IF ( npex == 1 ) THEN nc_stat = NF90_VAR_PAR_ACCESS( id_set_yz(av), & id_var_do2d(av,i), & NF90_COLLECTIVE ) ELSE ! !-- ATTENTION: Due to a probable bug in the NetCDF4 !-- installation, independet output does not work !-- A workaround is used in data_output_2d on those !-- PEs having no data nc_stat = NF90_VAR_PAR_ACCESS( id_set_yz(av), & id_var_do2d(av,i), & NF90_COLLECTIVE ) ! nc_stat = NF90_VAR_PAR_ACCESS( id_set_yz(av), & ! id_var_do2d(av,i), & ! NF90_INDEPENDENT ) ENDIF CALL handle_netcdf_error( 'netcdf', 450 ) ENDIF #endif ENDIF i = i + 1 ENDDO ! !-- Update the title attribute on file !-- In order to avoid 'data mode' errors if updated attributes are larger !-- than their original size, NF90_PUT_ATT is called in 'define mode' !-- enclosed by NF90_REDEF and NF90_ENDDEF calls. This implies a possible !-- performance loss due to data copying; an alternative strategy would be !-- to ensure equal attribute size in a job chain. Maybe revise later. IF ( av == 0 ) THEN time_average_text = ' ' ELSE WRITE (time_average_text, '('', '',F7.1,'' s average'')') & averaging_interval ENDIF nc_stat = NF90_REDEF( id_set_yz(av) ) CALL handle_netcdf_error( 'netcdf', 435 ) nc_stat = NF90_PUT_ATT( id_set_yz(av), NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 217 ) nc_stat = NF90_ENDDEF( id_set_yz(av) ) CALL handle_netcdf_error( 'netcdf', 436 ) message_string = 'NetCDF file for cross-sections ' // & TRIM( var ) // ' from previous run found.' // & '&This file will be extended.' CALL message( 'define_netcdf_header', 'PA0253', 0, 0, 0, 6, 0 ) CASE ( 'pr_new' ) ! !-- Define some global attributes of the dataset IF ( averaging_interval_pr /= 0.0 ) THEN WRITE (time_average_text,'('', '',F7.1,'' s average'')') & averaging_interval_pr nc_stat = NF90_PUT_ATT( id_set_pr, NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 218 ) WRITE ( time_average_text,'(F7.1,'' s avg'')' ) averaging_interval_pr nc_stat = NF90_PUT_ATT( id_set_pr, NF90_GLOBAL, 'time_avg', & TRIM( time_average_text ) ) ELSE nc_stat = NF90_PUT_ATT( id_set_pr, NF90_GLOBAL, 'title', & TRIM( run_description_header ) ) ENDIF CALL handle_netcdf_error( 'netcdf', 219 ) ! !-- Define time coordinate for profiles (unlimited dimension) nc_stat = NF90_DEF_DIM( id_set_pr, 'time', NF90_UNLIMITED, & id_dim_time_pr ) CALL handle_netcdf_error( 'netcdf', 220 ) nc_stat = NF90_DEF_VAR( id_set_pr, 'time', NF90_DOUBLE, & id_dim_time_pr, id_var_time_pr ) CALL handle_netcdf_error( 'netcdf', 221 ) nc_stat = NF90_PUT_ATT( id_set_pr, id_var_time_pr, 'units', 'seconds') CALL handle_netcdf_error( 'netcdf', 222 ) ! !-- Define the variables var_list = ';' DO i = 1, dopr_n IF ( statistic_regions == 0 ) THEN ! !-- Define the z-axes (each variable gets its own z-axis) nc_stat = NF90_DEF_DIM( id_set_pr, & 'z' // TRIM( data_output_pr(i) ), & nzt+2-nzb, id_dim_z_pr(i,0) ) CALL handle_netcdf_error( 'netcdf', 223 ) nc_stat = NF90_DEF_VAR( id_set_pr, & 'z' // TRIM( data_output_pr(i) ), & NF90_DOUBLE, id_dim_z_pr(i,0), & id_var_z_pr(i,0) ) CALL handle_netcdf_error( 'netcdf', 224 ) nc_stat = NF90_PUT_ATT( id_set_pr, id_var_z_pr(i,0), 'units', & 'meters' ) CALL handle_netcdf_error( 'netcdf', 225 ) ! !-- Define the variable nc_stat = NF90_DEF_VAR( id_set_pr, data_output_pr(i), & nc_precision(5), (/ id_dim_z_pr(i,0), & id_dim_time_pr /), id_var_dopr(i,0) ) CALL handle_netcdf_error( 'netcdf', 226 ) nc_stat = NF90_PUT_ATT( id_set_pr, id_var_dopr(i,0), & 'long_name', TRIM( data_output_pr(i) ) ) CALL handle_netcdf_error( 'netcdf', 227 ) nc_stat = NF90_PUT_ATT( id_set_pr, id_var_dopr(i,0), & 'units', TRIM( dopr_unit(i) ) ) CALL handle_netcdf_error( 'netcdf', 228 ) var_list = TRIM( var_list ) // TRIM( data_output_pr(i) ) // ';' ELSE ! !-- If statistic regions are defined, add suffix _SR+#SR to the !-- names DO j = 0, statistic_regions WRITE ( suffix, '(''_'',I1)' ) j ! !-- Define the z-axes (each variable gets it own z-axis) nc_stat = NF90_DEF_DIM( id_set_pr, & 'z'//TRIM(data_output_pr(i))//suffix, & nzt+2-nzb, id_dim_z_pr(i,j) ) CALL handle_netcdf_error( 'netcdf', 229 ) nc_stat = NF90_DEF_VAR( id_set_pr, & 'z'//TRIM(data_output_pr(i))//suffix, & nc_precision(5), id_dim_z_pr(i,j), & id_var_z_pr(i,j) ) CALL handle_netcdf_error( 'netcdf', 230 ) nc_stat = NF90_PUT_ATT( id_set_pr, id_var_z_pr(i,j), & 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 231 ) ! !-- Define the variable nc_stat = NF90_DEF_VAR( id_set_pr, & TRIM(data_output_pr(i)) // suffix, & nc_precision(5), & (/ id_dim_z_pr(i,j), & id_dim_time_pr /), id_var_dopr(i,j) ) CALL handle_netcdf_error( 'netcdf', 232 ) nc_stat = NF90_PUT_ATT( id_set_pr, id_var_dopr(i,j), & 'long_name', & TRIM( data_output_pr(i) ) // ' SR ' & // suffix(2:2) ) CALL handle_netcdf_error( 'netcdf', 233 ) nc_stat = NF90_PUT_ATT( id_set_pr, id_var_dopr(i,j), & 'units', TRIM( dopr_unit(i) ) ) CALL handle_netcdf_error( 'netcdf', 234 ) var_list = TRIM( var_list ) // TRIM( data_output_pr(i) ) // & suffix // ';' ENDDO ENDIF ENDDO ! !-- Write the list of variables as global attribute (this is used by !-- restart runs) nc_stat = NF90_PUT_ATT( id_set_pr, NF90_GLOBAL, 'VAR_LIST', var_list ) CALL handle_netcdf_error( 'netcdf', 235 ) ! !-- Define normalization variables (as time series) DO i = 1, dopr_norm_num nc_stat = NF90_DEF_VAR( id_set_pr, 'NORM_' // & TRIM( dopr_norm_names(i) ), & nc_precision(5), (/ id_dim_time_pr /), & id_var_norm_dopr(i) ) CALL handle_netcdf_error( 'netcdf', 236 ) nc_stat = NF90_PUT_ATT( id_set_pr, id_var_norm_dopr(i), & 'long_name', & TRIM( dopr_norm_longnames(i) ) ) CALL handle_netcdf_error( 'netcdf', 237 ) ENDDO ! !-- Leave NetCDF define mode nc_stat = NF90_ENDDEF( id_set_pr ) CALL handle_netcdf_error( 'netcdf', 238 ) ! !-- Write z-axes data DO i = 1, dopr_n DO j = 0, statistic_regions nc_stat = NF90_PUT_VAR( id_set_pr, id_var_z_pr(i,j), & hom(nzb:nzt+1,2,dopr_index(i),0), & start = (/ 1 /), & count = (/ nzt-nzb+2 /) ) CALL handle_netcdf_error( 'netcdf', 239 ) ENDDO ENDDO CASE ( 'pr_ext' ) ! !-- Get the list of variables and compare with the actual run. !-- First var_list_old has to be reset, since GET_ATT does not assign !-- trailing blanks. var_list_old = ' ' nc_stat = NF90_GET_ATT( id_set_pr, NF90_GLOBAL, 'VAR_LIST', & var_list_old ) CALL handle_netcdf_error( 'netcdf', 240 ) var_list = ';' DO i = 1, dopr_n IF ( statistic_regions == 0 ) THEN var_list = TRIM( var_list ) // TRIM( data_output_pr(i) ) // ';' ELSE DO j = 0, statistic_regions WRITE ( suffix, '(''_'',I1)' ) j var_list = TRIM( var_list ) // TRIM( data_output_pr(i) ) // & suffix // ';' ENDDO ENDIF ENDDO IF ( TRIM( var_list ) /= TRIM( var_list_old ) ) THEN message_string = 'NetCDF file for vertical profiles ' // & 'from previous run found,' // & '& but this file cannot be extended due to' // & ' variable mismatch.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0254', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Get the id of the time coordinate (unlimited coordinate) and its !-- last index on the file. The next time level is dopr..count+1. !-- The current time must be larger than the last output time !-- on the file. nc_stat = NF90_INQ_VARID( id_set_pr, 'time', id_var_time_pr ) CALL handle_netcdf_error( 'netcdf', 241 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_pr, id_var_time_pr, & dimids = id_dim_time_old ) CALL handle_netcdf_error( 'netcdf', 242 ) id_dim_time_pr = id_dim_time_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_pr, id_dim_time_pr, & len = dopr_time_count ) CALL handle_netcdf_error( 'netcdf', 243 ) nc_stat = NF90_GET_VAR( id_set_pr, id_var_time_pr, & last_time_coordinate, & start = (/ dopr_time_count /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'netcdf', 244 ) IF ( last_time_coordinate(1) >= simulated_time ) THEN message_string = 'NetCDF file for vertical profiles ' // & 'from previous run found,' // & '&but this file cannot be extended becaus' // & 'e the current output time' // & '&is less or equal than the last output t' // & 'ime on this file.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0255', 0, 1, 0, 6, 0 ) dopr_time_count = 0 extend = .FALSE. RETURN ENDIF ! !-- Dataset seems to be extendable. !-- Now get the variable ids. i = 1 DO i = 1, dopr_n IF ( statistic_regions == 0 ) THEN nc_stat = NF90_INQ_VARID( id_set_pr, data_output_pr(i), & id_var_dopr(i,0) ) CALL handle_netcdf_error( 'netcdf', 245 ) ELSE DO j = 0, statistic_regions WRITE ( suffix, '(''_'',I1)' ) j netcdf_var_name = TRIM( data_output_pr(i) ) // suffix nc_stat = NF90_INQ_VARID( id_set_pr, netcdf_var_name, & id_var_dopr(i,j) ) CALL handle_netcdf_error( 'netcdf', 246 ) ENDDO ENDIF ENDDO ! !-- Get ids of the normalization variables DO i = 1, dopr_norm_num nc_stat = NF90_INQ_VARID( id_set_pr, & 'NORM_' // TRIM( dopr_norm_names(i) ), & id_var_norm_dopr(i) ) CALL handle_netcdf_error( 'netcdf', 247 ) ENDDO ! !-- Update the title attribute on file !-- In order to avoid 'data mode' errors if updated attributes are larger !-- than their original size, NF90_PUT_ATT is called in 'define mode' !-- enclosed by NF90_REDEF and NF90_ENDDEF calls. This implies a possible !-- performance loss due to data copying; an alternative strategy would be !-- to ensure equal attribute size in a job chain. Maybe revise later. IF ( averaging_interval_pr == 0.0 ) THEN time_average_text = ' ' ELSE WRITE (time_average_text, '('', '',F7.1,'' s average'')') & averaging_interval_pr ENDIF nc_stat = NF90_REDEF( id_set_pr ) CALL handle_netcdf_error( 'netcdf', 437 ) nc_stat = NF90_PUT_ATT( id_set_pr, NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 248 ) nc_stat = NF90_ENDDEF( id_set_pr ) CALL handle_netcdf_error( 'netcdf', 438 ) message_string = 'NetCDF file for vertical profiles ' // & 'from previous run found.' // & '&This file will be extended.' CALL message( 'define_netcdf_header', 'PA0256', 0, 0, 0, 6, 0 ) CASE ( 'ts_new' ) ! !-- Define some global attributes of the dataset nc_stat = NF90_PUT_ATT( id_set_ts, NF90_GLOBAL, 'title', & TRIM( run_description_header ) ) CALL handle_netcdf_error( 'netcdf', 249 ) ! !-- Define time coordinate for time series (unlimited dimension) nc_stat = NF90_DEF_DIM( id_set_ts, 'time', NF90_UNLIMITED, & id_dim_time_ts ) CALL handle_netcdf_error( 'netcdf', 250 ) nc_stat = NF90_DEF_VAR( id_set_ts, 'time', NF90_DOUBLE, & id_dim_time_ts, id_var_time_ts ) CALL handle_netcdf_error( 'netcdf', 251 ) nc_stat = NF90_PUT_ATT( id_set_ts, id_var_time_ts, 'units', 'seconds') CALL handle_netcdf_error( 'netcdf', 252 ) ! !-- Define the variables var_list = ';' DO i = 1, dots_num IF ( statistic_regions == 0 ) THEN nc_stat = NF90_DEF_VAR( id_set_ts, dots_label(i), & nc_precision(6), (/ id_dim_time_ts /), & id_var_dots(i,0) ) CALL handle_netcdf_error( 'netcdf', 253 ) nc_stat = NF90_PUT_ATT( id_set_ts, id_var_dots(i,0), & 'long_name', TRIM( dots_label(i) ) ) CALL handle_netcdf_error( 'netcdf', 254 ) nc_stat = NF90_PUT_ATT( id_set_ts, id_var_dots(i,0), & 'units', TRIM( dots_unit(i) ) ) CALL handle_netcdf_error( 'netcdf', 255 ) var_list = TRIM( var_list ) // TRIM( dots_label(i) ) // ';' ELSE ! !-- If statistic regions are defined, add suffix _SR+#SR to the !-- names DO j = 0, statistic_regions WRITE ( suffix, '(''_'',I1)' ) j nc_stat = NF90_DEF_VAR( id_set_ts, & TRIM( dots_label(i) ) // suffix, & nc_precision(6), & (/ id_dim_time_ts /), & id_var_dots(i,j) ) CALL handle_netcdf_error( 'netcdf', 256 ) nc_stat = NF90_PUT_ATT( id_set_ts, id_var_dots(i,j), & 'long_name', & TRIM( dots_label(i) ) // ' SR ' // & suffix(2:2) ) CALL handle_netcdf_error( 'netcdf', 257 ) nc_stat = NF90_PUT_ATT( id_set_ts, id_var_dots(i,j), & 'units', TRIM( dots_unit(i) ) ) CALL handle_netcdf_error( 'netcdf', 347 ) var_list = TRIM( var_list ) // TRIM( dots_label(i) ) // & suffix // ';' ENDDO ENDIF ENDDO ! !-- Write the list of variables as global attribute (this is used by !-- restart runs) nc_stat = NF90_PUT_ATT( id_set_ts, NF90_GLOBAL, 'VAR_LIST', var_list ) CALL handle_netcdf_error( 'netcdf', 258 ) ! !-- Leave NetCDF define mode nc_stat = NF90_ENDDEF( id_set_ts ) CALL handle_netcdf_error( 'netcdf', 259 ) CASE ( 'ts_ext' ) ! !-- Get the list of variables and compare with the actual run. !-- First var_list_old has to be reset, since GET_ATT does not assign !-- trailing blanks. var_list_old = ' ' nc_stat = NF90_GET_ATT( id_set_ts, NF90_GLOBAL, 'VAR_LIST', & var_list_old ) CALL handle_netcdf_error( 'netcdf', 260 ) var_list = ';' i = 1 DO i = 1, dots_num IF ( statistic_regions == 0 ) THEN var_list = TRIM( var_list ) // TRIM( dots_label(i) ) // ';' ELSE DO j = 0, statistic_regions WRITE ( suffix, '(''_'',I1)' ) j var_list = TRIM( var_list ) // TRIM( dots_label(i) ) // & suffix // ';' ENDDO ENDIF ENDDO IF ( TRIM( var_list ) /= TRIM( var_list_old ) ) THEN message_string = 'NetCDF file for time series ' // & 'from previous run found,' // & '& but this file cannot be extended due to' // & ' variable mismatch.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0257', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Get the id of the time coordinate (unlimited coordinate) and its !-- last index on the file. The next time level is dots..count+1. !-- The current time must be larger than the last output time !-- on the file. nc_stat = NF90_INQ_VARID( id_set_ts, 'time', id_var_time_ts ) CALL handle_netcdf_error( 'netcdf', 261 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_ts, id_var_time_ts, & dimids = id_dim_time_old ) CALL handle_netcdf_error( 'netcdf', 262 ) id_dim_time_ts = id_dim_time_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_ts, id_dim_time_ts, & len = dots_time_count ) CALL handle_netcdf_error( 'netcdf', 263 ) nc_stat = NF90_GET_VAR( id_set_ts, id_var_time_ts, & last_time_coordinate, & start = (/ dots_time_count /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'netcdf', 264 ) IF ( last_time_coordinate(1) >= simulated_time ) THEN message_string = 'NetCDF file for time series ' // & 'from previous run found,' // & '&but this file cannot be extended becaus' // & 'e the current output time' // & '&is less or equal than the last output t' // & 'ime on this file.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0258', 0, 1, 0, 6, 0 ) dots_time_count = 0 extend = .FALSE. RETURN ENDIF ! !-- Dataset seems to be extendable. !-- Now get the variable ids i = 1 DO i = 1, dots_num IF ( statistic_regions == 0 ) THEN nc_stat = NF90_INQ_VARID( id_set_ts, dots_label(i), & id_var_dots(i,0) ) CALL handle_netcdf_error( 'netcdf', 265 ) ELSE DO j = 0, statistic_regions WRITE ( suffix, '(''_'',I1)' ) j netcdf_var_name = TRIM( dots_label(i) ) // suffix nc_stat = NF90_INQ_VARID( id_set_ts, netcdf_var_name, & id_var_dots(i,j) ) CALL handle_netcdf_error( 'netcdf', 266 ) ENDDO ENDIF ENDDO ! !-- Update the title attribute on file !-- In order to avoid 'data mode' errors if updated attributes are larger !-- than their original size, NF90_PUT_ATT is called in 'define mode' !-- enclosed by NF90_REDEF and NF90_ENDDEF calls. This implies a possible !-- performance loss due to data copying; an alternative strategy would be !-- to ensure equal attribute size in a job chain. Maybe revise later. nc_stat = NF90_REDEF( id_set_ts ) CALL handle_netcdf_error( 'netcdf', 439 ) nc_stat = NF90_PUT_ATT( id_set_ts, NF90_GLOBAL, 'title', & TRIM( run_description_header ) ) CALL handle_netcdf_error( 'netcdf', 267 ) nc_stat = NF90_ENDDEF( id_set_ts ) CALL handle_netcdf_error( 'netcdf', 440 ) message_string = 'NetCDF file for time series ' // & 'from previous run found.' // & '&This file will be extended.' CALL message( 'define_netcdf_header', 'PA0259', 0, 0, 0, 6, 0 ) CASE ( 'sp_new' ) ! !-- Define some global attributes of the dataset IF ( averaging_interval_sp /= 0.0 ) THEN WRITE (time_average_text,'('', '',F7.1,'' s average'')') & averaging_interval_sp nc_stat = NF90_PUT_ATT( id_set_sp, NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 268 ) WRITE ( time_average_text,'(F7.1,'' s avg'')' ) averaging_interval_sp nc_stat = NF90_PUT_ATT( id_set_sp, NF90_GLOBAL, 'time_avg', & TRIM( time_average_text ) ) ELSE nc_stat = NF90_PUT_ATT( id_set_sp, NF90_GLOBAL, 'title', & TRIM( run_description_header ) ) ENDIF CALL handle_netcdf_error( 'netcdf', 269 ) ! !-- Define time coordinate for spectra (unlimited dimension) nc_stat = NF90_DEF_DIM( id_set_sp, 'time', NF90_UNLIMITED, & id_dim_time_sp ) CALL handle_netcdf_error( 'netcdf', 270 ) nc_stat = NF90_DEF_VAR( id_set_sp, 'time', NF90_DOUBLE, & id_dim_time_sp, id_var_time_sp ) CALL handle_netcdf_error( 'netcdf', 271 ) nc_stat = NF90_PUT_ATT( id_set_sp, id_var_time_sp, 'units', 'seconds') CALL handle_netcdf_error( 'netcdf', 272 ) ! !-- Define the spatial dimensions and coordinates for spectra. !-- First, determine the number of vertical levels for which spectra !-- are to be output. ns = 1 DO WHILE ( comp_spectra_level(ns) /= 999999 .AND. ns <= 100 ) ns = ns + 1 ENDDO ns = ns - 1 ! !-- Define vertical coordinate grid (zu grid) nc_stat = NF90_DEF_DIM( id_set_sp, 'zu_sp', ns, id_dim_zu_sp ) CALL handle_netcdf_error( 'netcdf', 273 ) nc_stat = NF90_DEF_VAR( id_set_sp, 'zu_sp', NF90_DOUBLE, & id_dim_zu_sp, id_var_zu_sp ) CALL handle_netcdf_error( 'netcdf', 274 ) nc_stat = NF90_PUT_ATT( id_set_sp, id_var_zu_sp, 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 275 ) ! !-- Define vertical coordinate grid (zw grid) nc_stat = NF90_DEF_DIM( id_set_sp, 'zw_sp', ns, id_dim_zw_sp ) CALL handle_netcdf_error( 'netcdf', 276 ) nc_stat = NF90_DEF_VAR( id_set_sp, 'zw_sp', NF90_DOUBLE, & id_dim_zw_sp, id_var_zw_sp ) CALL handle_netcdf_error( 'netcdf', 277 ) nc_stat = NF90_PUT_ATT( id_set_sp, id_var_zw_sp, 'units', 'meters' ) CALL handle_netcdf_error( 'netcdf', 278 ) ! !-- Define x-axis nc_stat = NF90_DEF_DIM( id_set_sp, 'k_x', nx/2, id_dim_x_sp ) CALL handle_netcdf_error( 'netcdf', 279 ) nc_stat = NF90_DEF_VAR( id_set_sp, 'k_x', NF90_DOUBLE, id_dim_x_sp, & id_var_x_sp ) CALL handle_netcdf_error( 'netcdf', 280 ) nc_stat = NF90_PUT_ATT( id_set_sp, id_var_x_sp, 'units', 'm-1' ) CALL handle_netcdf_error( 'netcdf', 281 ) ! !-- Define y-axis nc_stat = NF90_DEF_DIM( id_set_sp, 'k_y', ny/2, id_dim_y_sp ) CALL handle_netcdf_error( 'netcdf', 282 ) nc_stat = NF90_DEF_VAR( id_set_sp, 'k_y', NF90_DOUBLE, id_dim_y_sp, & id_var_y_sp ) CALL handle_netcdf_error( 'netcdf', 283 ) nc_stat = NF90_PUT_ATT( id_set_sp, id_var_y_sp, 'units', 'm-1' ) CALL handle_netcdf_error( 'netcdf', 284 ) ! !-- Define the variables var_list = ';' i = 1 DO WHILE ( data_output_sp(i) /= ' ' .AND. i <= 10 ) ! !-- First check for the vertical grid found = .TRUE. SELECT CASE ( data_output_sp(i) ) ! !-- Most variables are defined on the zu levels CASE ( 'e', 'p', 'pc', 'pr', 'pt', 'q', 'ql', 'ql_c', 'ql_v', & 'ql_vp', 'qv', 'rho', 's', 'sa', 'u', 'v', 'vpt' ) grid_z = 'zu' ! !-- zw levels CASE ( 'w' ) grid_z = 'zw' CASE DEFAULT ! !-- Check for user-defined quantities (found, grid_x and grid_y !-- are dummies) CALL user_define_netcdf_grid( data_output_sp(i), found, & grid_x, grid_y, grid_z ) END SELECT IF ( INDEX( spectra_direction(i), 'x' ) /= 0 ) THEN ! !-- Define the variable netcdf_var_name = TRIM( data_output_sp(i) ) // '_x' IF ( TRIM( grid_z ) == 'zw' ) THEN nc_stat = NF90_DEF_VAR( id_set_sp, netcdf_var_name, & nc_precision(7), (/ id_dim_x_sp, & id_dim_zw_sp, id_dim_time_sp /), & id_var_dospx(i) ) ELSE nc_stat = NF90_DEF_VAR( id_set_sp, netcdf_var_name, & nc_precision(7), (/ id_dim_x_sp, & id_dim_zu_sp, id_dim_time_sp /), & id_var_dospx(i) ) ENDIF CALL handle_netcdf_error( 'netcdf', 285 ) nc_stat = NF90_PUT_ATT( id_set_sp, id_var_dospx(i), & 'long_name', netcdf_var_name ) CALL handle_netcdf_error( 'netcdf', 286 ) nc_stat = NF90_PUT_ATT( id_set_sp, id_var_dospx(i), & 'units', 'unknown' ) CALL handle_netcdf_error( 'netcdf', 287 ) var_list = TRIM( var_list ) // TRIM( netcdf_var_name ) // ';' ENDIF IF ( INDEX( spectra_direction(i), 'y' ) /= 0 ) THEN ! !-- Define the variable netcdf_var_name = TRIM( data_output_sp(i) ) // '_y' IF ( TRIM( grid_z ) == 'zw' ) THEN nc_stat = NF90_DEF_VAR( id_set_sp, netcdf_var_name, & nc_precision(7), (/ id_dim_y_sp, & id_dim_zw_sp, id_dim_time_sp /), & id_var_dospy(i) ) ELSE nc_stat = NF90_DEF_VAR( id_set_sp, netcdf_var_name, & nc_precision(7), (/ id_dim_y_sp, & id_dim_zu_sp, id_dim_time_sp /), & id_var_dospy(i) ) ENDIF CALL handle_netcdf_error( 'netcdf', 288 ) nc_stat = NF90_PUT_ATT( id_set_sp, id_var_dospy(i), & 'long_name', netcdf_var_name ) CALL handle_netcdf_error( 'netcdf', 289 ) nc_stat = NF90_PUT_ATT( id_set_sp, id_var_dospy(i), & 'units', 'unknown' ) CALL handle_netcdf_error( 'netcdf', 290 ) var_list = TRIM( var_list ) // TRIM( netcdf_var_name ) // ';' ENDIF i = i + 1 ENDDO ! !-- Write the list of variables as global attribute (this is used by !-- restart runs) nc_stat = NF90_PUT_ATT( id_set_sp, NF90_GLOBAL, 'VAR_LIST', var_list ) CALL handle_netcdf_error( 'netcdf', 291 ) ! !-- Leave NetCDF define mode nc_stat = NF90_ENDDEF( id_set_sp ) CALL handle_netcdf_error( 'netcdf', 292 ) ! !-- Write axis data: zu_sp, zw_sp, k_x, k_y ALLOCATE( netcdf_data(1:ns) ) ! !-- Write zu data netcdf_data(1:ns) = zu( comp_spectra_level(1:ns) ) nc_stat = NF90_PUT_VAR( id_set_sp, id_var_zu_sp, netcdf_data, & start = (/ 1 /), count = (/ ns /) ) CALL handle_netcdf_error( 'netcdf', 293 ) ! !-- Write zw data netcdf_data(1:ns) = zw( comp_spectra_level(1:ns) ) nc_stat = NF90_PUT_VAR( id_set_sp, id_var_zw_sp, netcdf_data, & start = (/ 1 /), count = (/ ns /) ) CALL handle_netcdf_error( 'netcdf', 294 ) DEALLOCATE( netcdf_data ) ! !-- Write data for x and y axis (wavenumbers) ALLOCATE( netcdf_data(nx/2) ) DO i = 1, nx/2 netcdf_data(i) = 2.0 * pi * i / ( dx * ( nx + 1 ) ) ENDDO nc_stat = NF90_PUT_VAR( id_set_sp, id_var_x_sp, netcdf_data, & start = (/ 1 /), count = (/ nx/2 /) ) CALL handle_netcdf_error( 'netcdf', 295 ) DEALLOCATE( netcdf_data ) ALLOCATE( netcdf_data(ny/2) ) DO i = 1, ny/2 netcdf_data(i) = 2.0 * pi * i / ( dy * ( ny + 1 ) ) ENDDO nc_stat = NF90_PUT_VAR( id_set_sp, id_var_y_sp, netcdf_data, & start = (/ 1 /), count = (/ ny/2 /) ) CALL handle_netcdf_error( 'netcdf', 296 ) DEALLOCATE( netcdf_data ) CASE ( 'sp_ext' ) ! !-- Get the list of variables and compare with the actual run. !-- First var_list_old has to be reset, since GET_ATT does not assign !-- trailing blanks. var_list_old = ' ' nc_stat = NF90_GET_ATT( id_set_sp, NF90_GLOBAL, 'VAR_LIST', & var_list_old ) CALL handle_netcdf_error( 'netcdf', 297 ) var_list = ';' i = 1 DO WHILE ( data_output_sp(i) /= ' ' .AND. i <= 10 ) IF ( INDEX( spectra_direction(i), 'x' ) /= 0 ) THEN netcdf_var_name = TRIM( data_output_sp(i) ) // '_x' var_list = TRIM( var_list ) // TRIM( netcdf_var_name ) // ';' ENDIF IF ( INDEX( spectra_direction(i), 'y' ) /= 0 ) THEN netcdf_var_name = TRIM( data_output_sp(i) ) // '_y' var_list = TRIM( var_list ) // TRIM( netcdf_var_name ) // ';' ENDIF i = i + 1 ENDDO IF ( TRIM( var_list ) /= TRIM( var_list_old ) ) THEN message_string = 'NetCDF file for spectra ' // & 'from previous run found,' // & '& but this file cannot be extended due to' // & ' variable mismatch.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0260', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Determine the number of current vertical levels for which spectra !-- shall be output ns = 1 DO WHILE ( comp_spectra_level(ns) /= 999999 .AND. ns <= 100 ) ns = ns + 1 ENDDO ns = ns - 1 ! !-- Get and compare the number of vertical levels nc_stat = NF90_INQ_VARID( id_set_sp, 'zu_sp', id_var_zu_sp ) CALL handle_netcdf_error( 'netcdf', 298 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_sp, id_var_zu_sp, & dimids = id_dim_zu_sp_old ) CALL handle_netcdf_error( 'netcdf', 299 ) id_dim_zu_sp = id_dim_zu_sp_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_sp, id_dim_zu_sp, & len = ns_old ) CALL handle_netcdf_error( 'netcdf', 300 ) IF ( ns /= ns_old ) THEN message_string = 'NetCDF file for spectra ' // & ' from previous run found,' // & '&but this file cannot be extended due to' // & ' mismatch in number of' // & '&vertical levels.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0261', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Get and compare the heights of the cross sections ALLOCATE( netcdf_data(1:ns_old) ) nc_stat = NF90_GET_VAR( id_set_sp, id_var_zu_sp, netcdf_data ) CALL handle_netcdf_error( 'netcdf', 301 ) DO i = 1, ns IF ( zu(comp_spectra_level(i)) /= netcdf_data(i) ) THEN message_string = 'NetCDF file for spectra ' // & ' from previous run found,' // & '&but this file cannot be extended due to' // & ' mismatch in heights of' // & '&vertical levels.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0262', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ENDDO DEALLOCATE( netcdf_data ) ! !-- Get the id of the time coordinate (unlimited coordinate) and its !-- last index on the file. The next time level is plsp..count+1. !-- The current time must be larger than the last output time !-- on the file. nc_stat = NF90_INQ_VARID( id_set_sp, 'time', id_var_time_sp ) CALL handle_netcdf_error( 'netcdf', 302 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_sp, id_var_time_sp, & dimids = id_dim_time_old ) CALL handle_netcdf_error( 'netcdf', 303 ) id_dim_time_sp = id_dim_time_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_sp, id_dim_time_sp, & len = dosp_time_count ) CALL handle_netcdf_error( 'netcdf', 304 ) nc_stat = NF90_GET_VAR( id_set_sp, id_var_time_sp, & last_time_coordinate, & start = (/ dosp_time_count /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'netcdf', 305 ) IF ( last_time_coordinate(1) >= simulated_time ) THEN message_string = 'NetCDF file for spectra ' // & 'from previous run found,' // & '&but this file cannot be extended becaus' // & 'e the current output time' // & '&is less or equal than the last output t' // & 'ime on this file.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0263', 0, 1, 0, 6, 0 ) dosp_time_count = 0 extend = .FALSE. RETURN ENDIF ! !-- Dataset seems to be extendable. !-- Now get the variable ids. i = 1 DO WHILE ( data_output_sp(i) /= ' ' .AND. i <= 10 ) IF ( INDEX( spectra_direction(i), 'x' ) /= 0 ) THEN netcdf_var_name = TRIM( data_output_sp(i) ) // '_x' nc_stat = NF90_INQ_VARID( id_set_sp, netcdf_var_name, & id_var_dospx(i) ) CALL handle_netcdf_error( 'netcdf', 306 ) ENDIF IF ( INDEX( spectra_direction(i), 'y' ) /= 0 ) THEN netcdf_var_name = TRIM( data_output_sp(i) ) // '_y' nc_stat = NF90_INQ_VARID( id_set_sp, netcdf_var_name, & id_var_dospy(i) ) CALL handle_netcdf_error( 'netcdf', 307 ) ENDIF i = i + 1 ENDDO ! !-- Update the title attribute on file !-- In order to avoid 'data mode' errors if updated attributes are larger !-- than their original size, NF90_PUT_ATT is called in 'define mode' !-- enclosed by NF90_REDEF and NF90_ENDDEF calls. This implies a possible !-- performance loss due to data copying; an alternative strategy would be !-- to ensure equal attribute size in a job chain. Maybe revise later. nc_stat = NF90_REDEF( id_set_sp ) CALL handle_netcdf_error( 'netcdf', 441 ) IF ( averaging_interval_sp /= 0.0 ) THEN WRITE (time_average_text,'('', '',F7.1,'' s average'')') & averaging_interval_sp nc_stat = NF90_PUT_ATT( id_set_sp, NF90_GLOBAL, 'title', & TRIM( run_description_header ) // & TRIM( time_average_text ) ) CALL handle_netcdf_error( 'netcdf', 308 ) WRITE ( time_average_text,'(F7.1,'' s avg'')' ) averaging_interval_sp nc_stat = NF90_PUT_ATT( id_set_sp, NF90_GLOBAL, 'time_avg', & TRIM( time_average_text ) ) ELSE nc_stat = NF90_PUT_ATT( id_set_sp, NF90_GLOBAL, 'title', & TRIM( run_description_header ) ) ENDIF CALL handle_netcdf_error( 'netcdf', 309 ) nc_stat = NF90_ENDDEF( id_set_sp ) CALL handle_netcdf_error( 'netcdf', 442 ) message_string = 'NetCDF file for spectra ' // & 'from previous run found.' // & '&This file will be extended.' CALL message( 'define_netcdf_header', 'PA0264', 0, 0, 0, 6, 0 ) CASE ( 'pt_new' ) ! !-- Define some global attributes of the dataset nc_stat = NF90_PUT_ATT( id_set_prt, NF90_GLOBAL, 'title', & TRIM( run_description_header ) ) CALL handle_netcdf_error( 'netcdf', 310 ) ! !-- Define time coordinate for particles (unlimited dimension) nc_stat = NF90_DEF_DIM( id_set_prt, 'time', NF90_UNLIMITED, & id_dim_time_prt ) CALL handle_netcdf_error( 'netcdf', 311 ) nc_stat = NF90_DEF_VAR( id_set_prt, 'time', NF90_DOUBLE, & id_dim_time_prt, id_var_time_prt ) CALL handle_netcdf_error( 'netcdf', 312 ) nc_stat = NF90_PUT_ATT( id_set_prt, id_var_time_prt, 'units', & 'seconds' ) CALL handle_netcdf_error( 'netcdf', 313 ) ! !-- Define particle coordinate (maximum particle number) IF ( netcdf_data_format < 3 ) THEN nc_stat = NF90_DEF_DIM( id_set_prt, 'particle_number', & maximum_number_of_particles, id_dim_prtnum) ELSE ! !-- NetCDF4 allows more than one unlimited dimension nc_stat = NF90_DEF_DIM( id_set_prt, 'particle_number', & NF90_UNLIMITED, id_dim_prtnum) ENDIF CALL handle_netcdf_error( 'netcdf', 314 ) nc_stat = NF90_DEF_VAR( id_set_prt, 'particle_number', NF90_DOUBLE, & id_dim_prtnum, id_var_prtnum ) CALL handle_netcdf_error( 'netcdf', 315 ) nc_stat = NF90_PUT_ATT( id_set_prt, id_var_prtnum, 'units', & 'particle number' ) CALL handle_netcdf_error( 'netcdf', 316 ) ! !-- Define variable which contains the real number of particles in use nc_stat = NF90_DEF_VAR( id_set_prt, 'real_num_of_prt', NF90_INT, & id_dim_time_prt, id_var_rnop_prt ) CALL handle_netcdf_error( 'netcdf', 317 ) nc_stat = NF90_PUT_ATT( id_set_prt, id_var_rnop_prt, 'units', & 'particle number' ) CALL handle_netcdf_error( 'netcdf', 318 ) ! !-- Define the variables DO i = 1, 17 nc_stat = NF90_DEF_VAR( id_set_prt, prt_var_names(i), & nc_precision(8), & (/ id_dim_prtnum, id_dim_time_prt /), & id_var_prt(i) ) CALL handle_netcdf_error( 'netcdf', 319 ) nc_stat = NF90_PUT_ATT( id_set_prt, id_var_prt(i), & 'long_name', TRIM( prt_var_names(i) ) ) CALL handle_netcdf_error( 'netcdf', 320 ) nc_stat = NF90_PUT_ATT( id_set_prt, id_var_prt(i), & 'units', TRIM( prt_var_units(i) ) ) CALL handle_netcdf_error( 'netcdf', 321 ) ENDDO ! !-- Leave NetCDF define mode nc_stat = NF90_ENDDEF( id_set_prt ) CALL handle_netcdf_error( 'netcdf', 322 ) CASE ( 'pt_ext' ) ! !-- Get the id of the time coordinate (unlimited coordinate) and its !-- last index on the file. The next time level is prt..count+1. !-- The current time must be larger than the last output time !-- on the file. nc_stat = NF90_INQ_VARID( id_set_prt, 'time', id_var_time_prt ) CALL handle_netcdf_error( 'netcdf', 323 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_prt, id_var_time_prt, & dimids = id_dim_time_old ) CALL handle_netcdf_error( 'netcdf', 324 ) id_dim_time_prt = id_dim_time_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_prt, id_dim_time_prt, & len = prt_time_count ) CALL handle_netcdf_error( 'netcdf', 325 ) nc_stat = NF90_GET_VAR( id_set_prt, id_var_time_prt, & last_time_coordinate, & start = (/ prt_time_count /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'netcdf', 326 ) IF ( last_time_coordinate(1) >= simulated_time ) THEN message_string = 'NetCDF file for particles ' // & 'from previous run found,' // & '&but this file cannot be extended becaus' // & 'e the current output time' // & '&is less or equal than the last output t' // & 'ime on this file.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0265', 0, 1, 0, 6, 0 ) prt_time_count = 0 extend = .FALSE. RETURN ENDIF ! !-- Dataset seems to be extendable. !-- Now get the variable ids. nc_stat = NF90_INQ_VARID( id_set_prt, 'real_num_of_prt', & id_var_rnop_prt ) CALL handle_netcdf_error( 'netcdf', 327 ) DO i = 1, 17 nc_stat = NF90_INQ_VARID( id_set_prt, prt_var_names(i), & id_var_prt(i) ) CALL handle_netcdf_error( 'netcdf', 328 ) ENDDO message_string = 'NetCDF file for particles ' // & 'from previous run found.' // & '&This file will be extended.' CALL message( 'define_netcdf_header', 'PA0266', 0, 0, 0, 6, 0 ) CASE ( 'ps_new' ) ! !-- Define some global attributes of the dataset nc_stat = NF90_PUT_ATT( id_set_pts, NF90_GLOBAL, 'title', & TRIM( run_description_header ) ) CALL handle_netcdf_error( 'netcdf', 396 ) ! !-- Define time coordinate for particle time series (unlimited dimension) nc_stat = NF90_DEF_DIM( id_set_pts, 'time', NF90_UNLIMITED, & id_dim_time_pts ) CALL handle_netcdf_error( 'netcdf', 397 ) nc_stat = NF90_DEF_VAR( id_set_pts, 'time', NF90_DOUBLE, & id_dim_time_pts, id_var_time_pts ) CALL handle_netcdf_error( 'netcdf', 398 ) nc_stat = NF90_PUT_ATT( id_set_pts, id_var_time_pts, 'units', & 'seconds') CALL handle_netcdf_error( 'netcdf', 399 ) ! !-- Define the variables. If more than one particle group is defined, !-- define seperate variables for each group var_list = ';' DO i = 1, dopts_num DO j = 0, number_of_particle_groups IF ( j == 0 ) THEN suffix1 = '' ELSE WRITE ( suffix1, '(''_'',I2.2)' ) j ENDIF nc_stat = NF90_DEF_VAR( id_set_pts, & TRIM( dopts_label(i) ) // suffix1, & nc_precision(6), & (/ id_dim_time_pts /), & id_var_dopts(i,j) ) CALL handle_netcdf_error( 'netcdf', 400 ) IF ( j == 0 ) THEN nc_stat = NF90_PUT_ATT( id_set_pts, id_var_dopts(i,j), & 'long_name', & TRIM( dopts_label(i) ) ) ELSE nc_stat = NF90_PUT_ATT( id_set_pts, id_var_dopts(i,j), & 'long_name', & TRIM( dopts_label(i) ) // ' PG ' // & suffix1(2:3) ) ENDIF CALL handle_netcdf_error( 'netcdf', 401 ) nc_stat = NF90_PUT_ATT( id_set_pts, id_var_dopts(i,j), & 'units', TRIM( dopts_unit(i) ) ) CALL handle_netcdf_error( 'netcdf', 402 ) var_list = TRIM( var_list ) // TRIM( dopts_label(i) ) // & suffix1 // ';' IF ( number_of_particle_groups == 1 ) EXIT ENDDO ENDDO ! !-- Write the list of variables as global attribute (this is used by !-- restart runs) nc_stat = NF90_PUT_ATT( id_set_pts, NF90_GLOBAL, 'VAR_LIST', & var_list ) CALL handle_netcdf_error( 'netcdf', 403 ) ! !-- Leave NetCDF define mode nc_stat = NF90_ENDDEF( id_set_pts ) CALL handle_netcdf_error( 'netcdf', 404 ) CASE ( 'ps_ext' ) ! !-- Get the list of variables and compare with the actual run. !-- First var_list_old has to be reset, since GET_ATT does not assign !-- trailing blanks. var_list_old = ' ' nc_stat = NF90_GET_ATT( id_set_pts, NF90_GLOBAL, 'VAR_LIST', & var_list_old ) CALL handle_netcdf_error( 'netcdf', 405 ) var_list = ';' i = 1 DO i = 1, dopts_num DO j = 0, number_of_particle_groups IF ( j == 0 ) THEN suffix1 = '' ELSE WRITE ( suffix1, '(''_'',I2.2)' ) j ENDIF var_list = TRIM( var_list ) // TRIM( dopts_label(i) ) // & suffix1 // ';' IF ( number_of_particle_groups == 1 ) EXIT ENDDO ENDDO IF ( TRIM( var_list ) /= TRIM( var_list_old ) ) THEN message_string = 'NetCDF file for particle time series ' // & 'from previous run found,' // & '& but this file cannot be extended due to' // & ' variable mismatch.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0267', 0, 1, 0, 6, 0 ) extend = .FALSE. RETURN ENDIF ! !-- Get the id of the time coordinate (unlimited coordinate) and its !-- last index on the file. The next time level is dots..count+1. !-- The current time must be larger than the last output time !-- on the file. nc_stat = NF90_INQ_VARID( id_set_pts, 'time', id_var_time_pts ) CALL handle_netcdf_error( 'netcdf', 406 ) nc_stat = NF90_INQUIRE_VARIABLE( id_set_pts, id_var_time_pts, & dimids = id_dim_time_old ) CALL handle_netcdf_error( 'netcdf', 407 ) id_dim_time_pts = id_dim_time_old(1) nc_stat = NF90_INQUIRE_DIMENSION( id_set_pts, id_dim_time_pts, & len = dopts_time_count ) CALL handle_netcdf_error( 'netcdf', 408 ) nc_stat = NF90_GET_VAR( id_set_pts, id_var_time_pts, & last_time_coordinate, & start = (/ dopts_time_count /), & count = (/ 1 /) ) CALL handle_netcdf_error( 'netcdf', 409 ) IF ( last_time_coordinate(1) >= simulated_time ) THEN message_string = 'NetCDF file for particle time series ' // & 'from previous run found,' // & '&but this file cannot be extended becaus' // & 'e the current output time' // & '&is less or equal than the last output t' // & 'ime on this file.' // & '&New file is created instead.' CALL message( 'define_netcdf_header', 'PA0268', 0, 1, 0, 6, 0 ) dopts_time_count = 0 extend = .FALSE. RETURN ENDIF ! !-- Dataset seems to be extendable. !-- Now get the variable ids i = 1 DO i = 1, dopts_num DO j = 0, number_of_particle_groups IF ( j == 0 ) THEN suffix1 = '' ELSE WRITE ( suffix1, '(''_'',I2.2)' ) j ENDIF netcdf_var_name = TRIM( dopts_label(i) ) // suffix1 nc_stat = NF90_INQ_VARID( id_set_pts, netcdf_var_name, & id_var_dopts(i,j) ) CALL handle_netcdf_error( 'netcdf', 410 ) IF ( number_of_particle_groups == 1 ) EXIT ENDDO ENDDO ! !-- Update the title attribute on file !-- In order to avoid 'data mode' errors if updated attributes are larger !-- than their original size, NF90_PUT_ATT is called in 'define mode' !-- enclosed by NF90_REDEF and NF90_ENDDEF calls. This implies a possible !-- performance loss due to data copying; an alternative strategy would be !-- to ensure equal attribute size in a job chain. Maybe revise later. nc_stat = NF90_REDEF( id_set_pts ) CALL handle_netcdf_error( 'netcdf', 443 ) nc_stat = NF90_PUT_ATT( id_set_pts, NF90_GLOBAL, 'title', & TRIM( run_description_header ) ) CALL handle_netcdf_error( 'netcdf', 411 ) nc_stat = NF90_ENDDEF( id_set_pts ) CALL handle_netcdf_error( 'netcdf', 444 ) message_string = 'NetCDF file for particle time series ' // & 'from previous run found.' // & '&This file will be extended.' CALL message( 'define_netcdf_header', 'PA0269', 0, 0, 0, 6, 0 ) CASE DEFAULT message_string = 'mode "' // TRIM( mode) // '" not supported' CALL message( 'define_netcdf_header', 'PA0270', 0, 0, 0, 6, 0 ) END SELECT #endif END SUBROUTINE define_netcdf_header SUBROUTINE handle_netcdf_error( routine_name, errno ) #if defined( __netcdf ) !------------------------------------------------------------------------------! ! ! Description: ! ------------ ! Prints out a text message corresponding to the current status. !------------------------------------------------------------------------------! USE control_parameters USE netcdf USE netcdf_control USE pegrid IMPLICIT NONE CHARACTER(LEN=6) :: message_identifier CHARACTER(LEN=*) :: routine_name INTEGER :: errno IF ( nc_stat /= NF90_NOERR ) THEN WRITE( message_identifier, '(''NC'',I4.4)' ) errno message_string = TRIM( NF90_STRERROR( nc_stat ) ) CALL message( routine_name, message_identifier, 2, 2, 0, 6, 1 ) ENDIF #endif END SUBROUTINE handle_netcdf_error