Ignore:
Timestamp:
Sep 10, 2019 5:03:24 PM (5 years ago)
Author:
suehring
Message:

Offline nesting: data input modularized; time variable is defined relative to time_utc_init, so that input data is correctly mapped to actual model time; checks rephrased and new checks for the time dimension added; Netcdf input: routine to retrieve dimension length renamed

File:
1 edited

Legend:

Unmodified
Added
Removed
  • palm/trunk/SOURCE/nesting_offl_mod.f90

    r4182 r4226  
    2525! -----------------
    2626! $Id$
     27! - Data input moved into nesting_offl_mod
     28! - check rephrased
     29! - time variable is now relative to time_utc_init
     30! - Define module specific data type for offline nesting in nesting_offl_mod
     31!
     32! 4182 2019-08-22 15:20:23Z scharf
    2733! Corrected "Former revisions" section
    2834!
     
    99105
    100106    USE arrays_3d,                                                             &
    101         ONLY:  dzw, e, diss, pt, pt_init, q, q_init, s, u, u_init, ug, v,      &
    102                v_init, vg, w, zu, zw
     107        ONLY:  dzw,                                                            &
     108               e,                                                              &
     109               diss,                                                           &
     110               pt,                                                             &
     111               pt_init,                                                        &
     112               q,                                                              &
     113               q_init,                                                         &
     114               rdf,                                                            &
     115               rdf_sc,                                                         &
     116               s,                                                              &
     117               u,                                                              &
     118               u_init,                                                         &
     119               ug,                                                             &
     120               v,                                                              &
     121               v_init,                                                         &
     122               vg,                                                             &
     123               w,                                                              &
     124               zu,                                                             &
     125               zw
     126
     127    USE basic_constants_and_equations_mod,                                     &
     128           ONLY:  g,                                                           &
     129                  pi
    103130                 
    104131    USE chem_modules,                                                          &
     
    111138               bc_dirichlet_r,                                                 &
    112139               bc_dirichlet_s,                                                 &
     140               coupling_char,                                                  &
    113141               dt_3d,                                                          &
    114142               dz,                                                             &
     
    132160               
    133161    USE cpulog,                                                                &
    134         ONLY:  cpu_log, log_point
     162        ONLY:  cpu_log,                                                        &
     163               log_point,                                                      &
     164               log_point_s
     165
     166    USE date_and_time_mod,                                                     &
     167        ONLY:  time_utc_init
    135168
    136169    USE grid_variables
     
    145178
    146179    USE netcdf_data_input_mod,                                                 &
    147         ONLY:  nest_offl
    148        
     180        ONLY:  check_existence,                                                &
     181               close_input_file,                                               &
     182               get_dimension_length,                                           &
     183               get_variable,                                                   &
     184               get_variable_pr,                                                &
     185               input_pids_dynamic,                                             &
     186               inquire_num_variables,                                          &
     187               inquire_variable_names,                                         &
     188               input_file_dynamic,                                             &
     189               num_var_pids,                                                   &
     190               open_read_file,                                                 &
     191               pids_id
     192
    149193    USE pegrid
    150194
     195    IMPLICIT NONE
     196
     197!
     198!-- Define data type for nesting in larger-scale models like COSMO.
     199!-- Data type comprises u, v, w, pt, and q at lateral and top boundaries.
     200    TYPE nest_offl_type
     201
     202       CHARACTER(LEN=16) ::  char_l = 'ls_forcing_left_'  !< leading substring for variables at left boundary
     203       CHARACTER(LEN=17) ::  char_n = 'ls_forcing_north_' !< leading substring for variables at north boundary 
     204       CHARACTER(LEN=17) ::  char_r = 'ls_forcing_right_' !< leading substring for variables at right boundary 
     205       CHARACTER(LEN=17) ::  char_s = 'ls_forcing_south_' !< leading substring for variables at south boundary
     206       CHARACTER(LEN=15) ::  char_t = 'ls_forcing_top_'   !< leading substring for variables at top boundary
     207
     208       CHARACTER(LEN=100), DIMENSION(:), ALLOCATABLE ::  var_names         !< list of variable in dynamic input file
     209       CHARACTER(LEN=100), DIMENSION(:), ALLOCATABLE ::  var_names_chem_l  !< names of mesoscale nested chemistry variables at left boundary
     210       CHARACTER(LEN=100), DIMENSION(:), ALLOCATABLE ::  var_names_chem_n  !< names of mesoscale nested chemistry variables at north boundary
     211       CHARACTER(LEN=100), DIMENSION(:), ALLOCATABLE ::  var_names_chem_r  !< names of mesoscale nested chemistry variables at right boundary
     212       CHARACTER(LEN=100), DIMENSION(:), ALLOCATABLE ::  var_names_chem_s  !< names of mesoscale nested chemistry variables at south boundary
     213       CHARACTER(LEN=100), DIMENSION(:), ALLOCATABLE ::  var_names_chem_t  !< names of mesoscale nested chemistry variables at top boundary
     214
     215       INTEGER(iwp) ::  nt     !< number of time levels in dynamic input file
     216       INTEGER(iwp) ::  nzu    !< number of vertical levels on scalar grid in dynamic input file
     217       INTEGER(iwp) ::  nzw    !< number of vertical levels on w grid in dynamic input file
     218       INTEGER(iwp) ::  tind   !< time index for reference time in mesoscale-offline nesting
     219       INTEGER(iwp) ::  tind_p !< time index for following time in mesoscale-offline nesting
     220
     221       LOGICAL      ::  init         = .FALSE. !< flag indicating that offline nesting is already initialized
     222
     223       LOGICAL, DIMENSION(:), ALLOCATABLE ::  chem_from_file_l !< flags inidicating whether left boundary data for chemistry is in dynamic input file 
     224       LOGICAL, DIMENSION(:), ALLOCATABLE ::  chem_from_file_n !< flags inidicating whether north boundary data for chemistry is in dynamic input file
     225       LOGICAL, DIMENSION(:), ALLOCATABLE ::  chem_from_file_r !< flags inidicating whether right boundary data for chemistry is in dynamic input file
     226       LOGICAL, DIMENSION(:), ALLOCATABLE ::  chem_from_file_s !< flags inidicating whether south boundary data for chemistry is in dynamic input file
     227       LOGICAL, DIMENSION(:), ALLOCATABLE ::  chem_from_file_t !< flags inidicating whether top boundary data for chemistry is in dynamic input file
     228
     229       REAL(wp), DIMENSION(:), ALLOCATABLE ::  surface_pressure !< time dependent surface pressure
     230       REAL(wp), DIMENSION(:), ALLOCATABLE ::  time             !< time levels in dynamic input file
     231       REAL(wp), DIMENSION(:), ALLOCATABLE ::  zu_atmos         !< vertical levels at scalar grid in dynamic input file
     232       REAL(wp), DIMENSION(:), ALLOCATABLE ::  zw_atmos         !< vertical levels at w grid in dynamic input file
     233
     234       REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ug         !< domain-averaged geostrophic component
     235       REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  vg         !< domain-averaged geostrophic component
     236
     237       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  u_left   !< u-component at left boundary
     238       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  v_left   !< v-component at left boundary
     239       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  w_left   !< w-component at left boundary
     240       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  q_left   !< mixing ratio at left boundary
     241       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  pt_left  !< potentital temperautre at left boundary
     242
     243       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  u_north  !< u-component at north boundary
     244       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  v_north  !< v-component at north boundary
     245       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  w_north  !< w-component at north boundary
     246       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  q_north  !< mixing ratio at north boundary
     247       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  pt_north !< potentital temperautre at north boundary
     248
     249       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  u_right  !< u-component at right boundary
     250       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  v_right  !< v-component at right boundary
     251       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  w_right  !< w-component at right boundary
     252       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  q_right  !< mixing ratio at right boundary
     253       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  pt_right !< potentital temperautre at right boundary
     254
     255       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  u_south  !< u-component at south boundary
     256       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  v_south  !< v-component at south boundary
     257       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  w_south  !< w-component at south boundary
     258       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  q_south  !< mixing ratio at south boundary
     259       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  pt_south !< potentital temperautre at south boundary
     260
     261       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  u_top    !< u-component at top boundary
     262       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  v_top    !< v-component at top boundary
     263       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  w_top    !< w-component at top boundary
     264       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  q_top    !< mixing ratio at top boundary
     265       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  pt_top   !< potentital temperautre at top boundary
     266       
     267       REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE ::  chem_left   !< chemical species at left boundary
     268       REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE ::  chem_north  !< chemical species at left boundary
     269       REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE ::  chem_right  !< chemical species at left boundary
     270       REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE ::  chem_south  !< chemical species at left boundary
     271       REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE ::  chem_top    !< chemical species at left boundary
     272
     273    END TYPE nest_offl_type
     274
     275    REAL(wp) ::  fac_dt              !< interpolation factor
    151276    REAL(wp) ::  zi_ribulk = 0.0_wp  !< boundary-layer depth according to bulk Richardson criterion, i.e. the height where Ri_bulk exceeds the critical
    152                                      !< bulk Richardson number of 0.25
     277                                     !< bulk Richardson number of 0.2
     278
     279    TYPE(nest_offl_type) ::  nest_offl  !< data structure for data input at lateral and top boundaries (provided by Inifor)
    153280   
    154281    SAVE
     
    159286           nesting_offl_calc_zi,                                               &
    160287           nesting_offl_check_parameters,                                      &
     288           nesting_offl_geostrophic_wind,                                      &
    161289           nesting_offl_header,                                                &
    162290           nesting_offl_init,                                                  &
     291           nesting_offl_input,                                                 &
     292           nesting_offl_interpolation_factor,                                  &
    163293           nesting_offl_mass_conservation,                                     &
    164294           nesting_offl_parin
     
    178308       MODULE PROCEDURE nesting_offl_check_parameters
    179309    END INTERFACE nesting_offl_check_parameters
     310
     311    INTERFACE nesting_offl_geostrophic_wind
     312       MODULE PROCEDURE nesting_offl_geostrophic_wind
     313    END INTERFACE nesting_offl_geostrophic_wind
    180314   
    181315    INTERFACE nesting_offl_header
     
    186320       MODULE PROCEDURE nesting_offl_init
    187321    END INTERFACE nesting_offl_init
     322
     323    INTERFACE nesting_offl_input
     324       MODULE PROCEDURE nesting_offl_input
     325    END INTERFACE nesting_offl_input
     326
     327    INTERFACE nesting_offl_interpolation_factor
     328       MODULE PROCEDURE nesting_offl_interpolation_factor
     329    END INTERFACE nesting_offl_interpolation_factor
    188330           
    189331    INTERFACE nesting_offl_mass_conservation
     
    197339 CONTAINS
    198340
     341!------------------------------------------------------------------------------!
     342! Description:
     343! ------------
     344!> Reads data at lateral and top boundaries derived from larger-scale model.
     345!------------------------------------------------------------------------------!
     346    SUBROUTINE nesting_offl_input
     347
     348       INTEGER(iwp) ::  n   !< running index for chemistry variables
     349       INTEGER(iwp) ::  t   !< running index time dimension
     350
     351!
     352!--    Initialize INIFOR forcing in first call.
     353       IF ( .NOT. nest_offl%init )  THEN
     354#if defined ( __netcdf )
     355!
     356!--       Open file in read-only mode
     357          CALL open_read_file( TRIM( input_file_dynamic ) //                   &
     358                               TRIM( coupling_char ), pids_id )
     359!
     360!--       At first, inquire all variable names.
     361          CALL inquire_num_variables( pids_id, num_var_pids )
     362!
     363!--       Allocate memory to store variable names.
     364          ALLOCATE( nest_offl%var_names(1:num_var_pids) )
     365          CALL inquire_variable_names( pids_id, nest_offl%var_names )
     366!
     367!--       Read time dimension, allocate memory and finally read time array
     368          CALL get_dimension_length( pids_id, nest_offl%nt, 'time' )
     369
     370          IF ( check_existence( nest_offl%var_names, 'time' ) )  THEN
     371             ALLOCATE( nest_offl%time(0:nest_offl%nt-1) )
     372             CALL get_variable( pids_id, 'time', nest_offl%time )
     373          ENDIF
     374!
     375!--       Read vertical dimension of scalar und w grid
     376          CALL get_dimension_length( pids_id, nest_offl%nzu, 'z' )
     377          CALL get_dimension_length( pids_id, nest_offl%nzw, 'zw' )
     378
     379          IF ( check_existence( nest_offl%var_names, 'z' ) )  THEN
     380             ALLOCATE( nest_offl%zu_atmos(1:nest_offl%nzu) )
     381             CALL get_variable( pids_id, 'z', nest_offl%zu_atmos )
     382          ENDIF
     383          IF ( check_existence( nest_offl%var_names, 'zw' ) )  THEN
     384             ALLOCATE( nest_offl%zw_atmos(1:nest_offl%nzw) )
     385             CALL get_variable( pids_id, 'zw', nest_offl%zw_atmos )
     386          ENDIF
     387!
     388!--       Read surface pressure
     389          IF ( check_existence( nest_offl%var_names,                           &
     390                                'surface_forcing_surface_pressure' ) )  THEN
     391             ALLOCATE( nest_offl%surface_pressure(0:nest_offl%nt-1) )
     392             CALL get_variable( pids_id,                                       &
     393                                'surface_forcing_surface_pressure',            &
     394                                nest_offl%surface_pressure )
     395          ENDIF
     396!
     397!--       Close input file
     398          CALL close_input_file( pids_id )
     399#endif
     400       ENDIF
     401!
     402!--    Check if dynamic driver data input is required.
     403       IF ( nest_offl%time(nest_offl%tind_p) <=                                &
     404            MAX( time_since_reference_point, 0.0_wp) + time_utc_init  .OR.                   &
     405            .NOT.  nest_offl%init )  THEN
     406          CONTINUE
     407!
     408!--    Return otherwise
     409       ELSE
     410          RETURN
     411       ENDIF
     412!
     413!--    CPU measurement
     414       CALL cpu_log( log_point_s(86), 'NetCDF input forcing', 'start' )
     415
     416!
     417!--    Obtain time index for current point in time. Note, the time coordinate
     418!--    in the input file is relative to time_utc_init. Since time_since_...
     419!--    is negativ when spinup is used, use MAX function to obtain correct
     420!--    time at the beginning.
     421       nest_offl%tind = MINLOC( ABS( nest_offl%time - (                        &
     422                                     time_utc_init +                           &
     423                                     MAX( time_since_reference_point, 0.0_wp) )&
     424                                   ), DIM = 1 ) - 1
     425       nest_offl%tind_p = nest_offl%tind + 1
     426!
     427!--    Open file in read-only mode
     428#if defined ( __netcdf )
     429       CALL open_read_file( TRIM( input_file_dynamic ) //                      &
     430                            TRIM( coupling_char ), pids_id )
     431!
     432!--    Read geostrophic wind components
     433       DO  t = nest_offl%tind, nest_offl%tind_p
     434          CALL get_variable_pr( pids_id, 'ls_forcing_ug', t+1,                 &
     435                                nest_offl%ug(t-nest_offl%tind,nzb+1:nzt) )
     436          CALL get_variable_pr( pids_id, 'ls_forcing_vg', t+1,                 &
     437                                nest_offl%vg(t-nest_offl%tind,nzb+1:nzt) )
     438       ENDDO
     439!
     440!--    Read data at lateral and top boundaries. Please note, at left and
     441!--    right domain boundary, yz-layers are read for u, v, w, pt and q.
     442!--    For the v-component, the data starts at nysv, while for the other
     443!--    quantities the data starts at nys. This is equivalent at the north
     444!--    and south domain boundary for the u-component.
     445!--    Note, lateral data is also accessed by parallel IO, which is the reason
     446!--    why different arguments are passed depending on the boundary control
     447!--    flags. Cores that do not belong to the respective boundary just make
     448!--    a dummy read with count = 0, just in order to participate the collective
     449!--    operation.
     450!--    Read data for western boundary   
     451       CALL get_variable( pids_id, 'ls_forcing_left_u',                        &
     452                          nest_offl%u_left,                                    & ! array to be read
     453                          MERGE( nys+1, 1, bc_dirichlet_l),                    & ! start index y direction
     454                          MERGE( nzb+1, 1, bc_dirichlet_l),                    & ! start index z direction
     455                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_l),         & ! start index time dimension
     456                          MERGE( nyn-nys+1, 0, bc_dirichlet_l),                & ! number of elements along y
     457                          MERGE( nest_offl%nzu, 0, bc_dirichlet_l),            & ! number of elements alogn z
     458                          MERGE( 2, 0, bc_dirichlet_l),                        & ! number of time steps (2 or 0)
     459                          .TRUE. )                                               ! parallel IO when compiled accordingly
     460     
     461       CALL get_variable( pids_id, 'ls_forcing_left_v',                        &
     462                          nest_offl%v_left,                                    &
     463                          MERGE( nysv, 1, bc_dirichlet_l),                     &
     464                          MERGE( nzb+1, 1, bc_dirichlet_l),                    &
     465                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_l),         &
     466                          MERGE( nyn-nysv+1, 0, bc_dirichlet_l),               &
     467                          MERGE( nest_offl%nzu, 0, bc_dirichlet_l),            &
     468                          MERGE( 2, 0, bc_dirichlet_l),                        &
     469                          .TRUE. )                                       
     470
     471       CALL get_variable( pids_id, 'ls_forcing_left_w',                        &
     472                          nest_offl%w_left,                                    &
     473                          MERGE( nys+1, 1, bc_dirichlet_l),                    &
     474                          MERGE( nzb+1, 1, bc_dirichlet_l),                    &
     475                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_l),         &
     476                          MERGE( nyn-nys+1, 0, bc_dirichlet_l),                &
     477                          MERGE( nest_offl%nzw, 0, bc_dirichlet_l),            &
     478                          MERGE( 2, 0, bc_dirichlet_l),                        &
     479                          .TRUE. )   
     480
     481       IF ( .NOT. neutral )  THEN
     482          CALL get_variable( pids_id, 'ls_forcing_left_pt',                    &
     483                             nest_offl%pt_left,                                &
     484                             MERGE( nys+1, 1, bc_dirichlet_l),                 &
     485                             MERGE( nzb+1, 1, bc_dirichlet_l),                 &
     486                             MERGE( nest_offl%tind+1, 1, bc_dirichlet_l),      &
     487                             MERGE( nyn-nys+1, 0, bc_dirichlet_l),             &
     488                             MERGE( nest_offl%nzu, 0, bc_dirichlet_l),         &
     489                             MERGE( 2, 0, bc_dirichlet_l),                     &
     490                             .TRUE. )
     491       ENDIF
     492
     493       IF ( humidity )  THEN
     494          CALL get_variable( pids_id, 'ls_forcing_left_qv',                    &
     495                             nest_offl%q_left,                                 &
     496                             MERGE( nys+1, 1, bc_dirichlet_l),                 &
     497                             MERGE( nzb+1, 1, bc_dirichlet_l),                 &
     498                             MERGE( nest_offl%tind+1, 1, bc_dirichlet_l),      &
     499                             MERGE( nyn-nys+1, 0, bc_dirichlet_l),             &
     500                             MERGE( nest_offl%nzu, 0, bc_dirichlet_l),         &
     501                             MERGE( 2, 0, bc_dirichlet_l),                     &
     502                             .TRUE. )
     503       ENDIF
     504       
     505       IF ( air_chemistry )  THEN
     506          DO  n = 1, UBOUND(nest_offl%var_names_chem_l, 1)
     507             IF ( check_existence( nest_offl%var_names,                        &
     508                                   nest_offl%var_names_chem_l(n) ) )  THEN 
     509                CALL get_variable( pids_id,                                    &
     510                           TRIM( nest_offl%var_names_chem_l(n) ),              &
     511                           nest_offl%chem_left(:,:,:,n),                       &
     512                           MERGE( nys+1, 1, bc_dirichlet_l),                   &
     513                           MERGE( nzb+1, 1, bc_dirichlet_l),                   &
     514                           MERGE( nest_offl%tind+1, 1, bc_dirichlet_l),        &
     515                           MERGE( nyn-nys+1, 0, bc_dirichlet_l),               &
     516                           MERGE( nest_offl%nzu, 0, bc_dirichlet_l),           &
     517                           MERGE( 2, 0, bc_dirichlet_l),                       &
     518                           .TRUE. )
     519                nest_offl%chem_from_file_l(n) = .TRUE.
     520             ENDIF
     521          ENDDO
     522       ENDIF
     523!
     524!--    Read data for eastern boundary   
     525       CALL get_variable( pids_id, 'ls_forcing_right_u',                       &
     526                          nest_offl%u_right,                                   &
     527                          MERGE( nys+1, 1, bc_dirichlet_r),                    &
     528                          MERGE( nzb+1, 1, bc_dirichlet_r),                    &
     529                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_r),         &
     530                          MERGE( nyn-nys+1, 0, bc_dirichlet_r),                &
     531                          MERGE( nest_offl%nzu, 0, bc_dirichlet_r),            &
     532                          MERGE( 2, 0, bc_dirichlet_r),                        &
     533                          .TRUE. )                                             
     534     
     535       CALL get_variable( pids_id, 'ls_forcing_right_v',                       &
     536                          nest_offl%v_right,                                   &
     537                          MERGE( nysv, 1, bc_dirichlet_r),                     &
     538                          MERGE( nzb+1, 1, bc_dirichlet_r),                    &
     539                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_r),         &
     540                          MERGE( nyn-nysv+1, 0, bc_dirichlet_r),               &
     541                          MERGE( nest_offl%nzu, 0, bc_dirichlet_r),            &
     542                          MERGE( 2, 0, bc_dirichlet_r),                        &
     543                          .TRUE. )                                             
     544
     545       CALL get_variable( pids_id, 'ls_forcing_right_w',                       &
     546                          nest_offl%w_right,                                   &
     547                          MERGE( nys+1, 1, bc_dirichlet_r),                    &
     548                          MERGE( nzb+1, 1, bc_dirichlet_r),                    &
     549                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_r),         &
     550                          MERGE( nyn-nys+1, 0, bc_dirichlet_r),                &
     551                          MERGE( nest_offl%nzw, 0, bc_dirichlet_r),            &
     552                          MERGE( 2, 0, bc_dirichlet_r),                        &
     553                          .TRUE. )   
     554
     555       IF ( .NOT. neutral )  THEN
     556          CALL get_variable( pids_id, 'ls_forcing_right_pt',                   &
     557                             nest_offl%pt_right,                               &
     558                             MERGE( nys+1, 1, bc_dirichlet_r),                 &
     559                             MERGE( nzb+1, 1, bc_dirichlet_r),                 &
     560                             MERGE( nest_offl%tind+1, 1, bc_dirichlet_r),      &
     561                             MERGE( nyn-nys+1, 0, bc_dirichlet_r),             &
     562                             MERGE( nest_offl%nzu, 0, bc_dirichlet_r),         &
     563                             MERGE( 2, 0, bc_dirichlet_r),                     &
     564                             .TRUE. )
     565       ENDIF
     566
     567       IF ( humidity )  THEN
     568          CALL get_variable( pids_id, 'ls_forcing_right_qv',                   &
     569                             nest_offl%q_right,                                &
     570                             MERGE( nys+1, 1, bc_dirichlet_r),                 &
     571                             MERGE( nzb+1, 1, bc_dirichlet_r),                 &
     572                             MERGE( nest_offl%tind+1, 1, bc_dirichlet_r),      &
     573                             MERGE( nyn-nys+1, 0, bc_dirichlet_r),             &
     574                             MERGE( nest_offl%nzu, 0, bc_dirichlet_r),         &
     575                             MERGE( 2, 0, bc_dirichlet_r),                     &
     576                             .TRUE. )
     577       ENDIF
     578       
     579       IF ( air_chemistry )  THEN
     580          DO  n = 1, UBOUND(nest_offl%var_names_chem_r, 1)
     581             IF ( check_existence( nest_offl%var_names,                        &
     582                                   nest_offl%var_names_chem_r(n) ) )  THEN     
     583                CALL get_variable( pids_id,                                    &
     584                           TRIM( nest_offl%var_names_chem_r(n) ),              &
     585                           nest_offl%chem_right(:,:,:,n),                      &
     586                           MERGE( nys+1, 1, bc_dirichlet_r),                   &
     587                           MERGE( nzb+1, 1, bc_dirichlet_r),                   &
     588                           MERGE( nest_offl%tind+1, 1, bc_dirichlet_r),        &
     589                           MERGE( nyn-nys+1, 0, bc_dirichlet_r),               &
     590                           MERGE( nest_offl%nzu, 0, bc_dirichlet_r),           &
     591                           MERGE( 2, 0, bc_dirichlet_r),                       &
     592                           .TRUE. )
     593                nest_offl%chem_from_file_r(n) = .TRUE.
     594             ENDIF
     595          ENDDO
     596       ENDIF
     597!
     598!--    Read data for northern boundary
     599       CALL get_variable( pids_id, 'ls_forcing_north_u',                       & ! array to be read
     600                          nest_offl%u_north,                                   & ! start index x direction
     601                          MERGE( nxlu, 1, bc_dirichlet_n ),                    & ! start index z direction
     602                          MERGE( nzb+1, 1, bc_dirichlet_n ),                   & ! start index time dimension
     603                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_n ),        & ! number of elements along x
     604                          MERGE( nxr-nxlu+1, 0, bc_dirichlet_n ),              & ! number of elements alogn z
     605                          MERGE( nest_offl%nzu, 0, bc_dirichlet_n ),           & ! number of time steps (2 or 0)
     606                          MERGE( 2, 0, bc_dirichlet_n ),                       & ! parallel IO when compiled accordingly
     607                          .TRUE. )                                             
     608                                                                               
     609       CALL get_variable( pids_id, 'ls_forcing_north_v',                       & ! array to be read
     610                          nest_offl%v_north,                                   & ! start index x direction
     611                          MERGE( nxl+1, 1, bc_dirichlet_n ),                   & ! start index z direction
     612                          MERGE( nzb+1, 1, bc_dirichlet_n ),                   & ! start index time dimension
     613                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_n ),        & ! number of elements along x
     614                          MERGE( nxr-nxl+1, 0, bc_dirichlet_n ),               & ! number of elements alogn z
     615                          MERGE( nest_offl%nzu, 0, bc_dirichlet_n ),           & ! number of time steps (2 or 0)
     616                          MERGE( 2, 0, bc_dirichlet_n ),                       & ! parallel IO when compiled accordingly
     617                          .TRUE. )                                             
     618                                                                               
     619       CALL get_variable( pids_id, 'ls_forcing_north_w',                       & ! array to be read
     620                          nest_offl%w_north,                                   & ! start index x direction
     621                          MERGE( nxl+1, 1, bc_dirichlet_n ),                   & ! start index z direction
     622                          MERGE( nzb+1, 1, bc_dirichlet_n ),                   & ! start index time dimension
     623                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_n ),        & ! number of elements along x
     624                          MERGE( nxr-nxl+1, 0, bc_dirichlet_n ),               & ! number of elements alogn z
     625                          MERGE( nest_offl%nzw, 0, bc_dirichlet_n ),           & ! number of time steps (2 or 0)
     626                          MERGE( 2, 0, bc_dirichlet_n ),                       & ! parallel IO when compiled accordingly
     627                          .TRUE. )                                             
     628                                                                               
     629       IF ( .NOT. neutral )  THEN                                             
     630          CALL get_variable( pids_id, 'ls_forcing_north_pt',                   & ! array to be read
     631                             nest_offl%pt_north,                               & ! start index x direction
     632                             MERGE( nxl+1, 1, bc_dirichlet_n ),                & ! start index z direction
     633                             MERGE( nzb+1, 1, bc_dirichlet_n ),                & ! start index time dimension
     634                             MERGE( nest_offl%tind+1, 1, bc_dirichlet_n ),     & ! number of elements along x
     635                             MERGE( nxr-nxl+1, 0, bc_dirichlet_n ),            & ! number of elements alogn z
     636                             MERGE( nest_offl%nzu, 0, bc_dirichlet_n ),        & ! number of time steps (2 or 0)
     637                             MERGE( 2, 0, bc_dirichlet_n ),                    & ! parallel IO when compiled accordingly
     638                             .TRUE. )                                             
     639       ENDIF                                                                   
     640       IF ( humidity )  THEN                                                   
     641          CALL get_variable( pids_id, 'ls_forcing_north_qv',                   & ! array to be read
     642                             nest_offl%q_north,                                & ! start index x direction
     643                             MERGE( nxl+1, 1, bc_dirichlet_n ),                & ! start index z direction
     644                             MERGE( nzb+1, 1, bc_dirichlet_n ),                & ! start index time dimension
     645                             MERGE( nest_offl%tind+1, 1, bc_dirichlet_n ),     & ! number of elements along x
     646                             MERGE( nxr-nxl+1, 0, bc_dirichlet_n ),            & ! number of elements alogn z
     647                             MERGE( nest_offl%nzu, 0, bc_dirichlet_n ),        & ! number of time steps (2 or 0)
     648                             MERGE( 2, 0, bc_dirichlet_n ),                    & ! parallel IO when compiled accordingly
     649                             .TRUE. )                                             
     650       ENDIF                                                                   
     651                                                                               
     652       IF ( air_chemistry )  THEN                                             
     653          DO  n = 1, UBOUND(nest_offl%var_names_chem_n, 1)                     
     654             IF ( check_existence( nest_offl%var_names,                        &
     655                                   nest_offl%var_names_chem_n(n) ) )  THEN     
     656                CALL get_variable( pids_id,                                    &
     657                           TRIM( nest_offl%var_names_chem_n(n) ),              &
     658                           nest_offl%chem_north(:,:,:,n),                      &
     659                           MERGE( nxl+1, 1, bc_dirichlet_n ),                  &
     660                           MERGE( nzb+1, 1, bc_dirichlet_n ),                  &
     661                           MERGE( nest_offl%tind+1, 1, bc_dirichlet_n ),       &
     662                           MERGE( nxr-nxl+1, 0, bc_dirichlet_n ),              &
     663                           MERGE( nest_offl%nzu, 0, bc_dirichlet_n ),          &
     664                           MERGE( 2, 0, bc_dirichlet_n ),                      &
     665                           .TRUE. )
     666                nest_offl%chem_from_file_n(n) = .TRUE.
     667             ENDIF
     668          ENDDO
     669       ENDIF
     670!
     671!--    Read data for southern boundary
     672       CALL get_variable( pids_id, 'ls_forcing_south_u',                       & ! array to be read
     673                          nest_offl%u_south,                                   & ! start index x direction
     674                          MERGE( nxlu, 1, bc_dirichlet_s ),                    & ! start index z direction
     675                          MERGE( nzb+1, 1, bc_dirichlet_s ),                   & ! start index time dimension
     676                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_s ),        & ! number of elements along x
     677                          MERGE( nxr-nxlu+1, 0, bc_dirichlet_s ),              & ! number of elements alogn z
     678                          MERGE( nest_offl%nzu, 0, bc_dirichlet_s ),           & ! number of time steps (2 or 0)
     679                          MERGE( 2, 0, bc_dirichlet_s ),                       & ! parallel IO when compiled accordingly
     680                          .TRUE. )                                             
     681                                                                               
     682       CALL get_variable( pids_id, 'ls_forcing_south_v',                       & ! array to be read
     683                          nest_offl%v_south,                                   & ! start index x direction
     684                          MERGE( nxl+1, 1, bc_dirichlet_s ),                   & ! start index z direction
     685                          MERGE( nzb+1, 1, bc_dirichlet_s ),                   & ! start index time dimension
     686                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_s ),        & ! number of elements along x
     687                          MERGE( nxr-nxl+1, 0, bc_dirichlet_s ),               & ! number of elements alogn z
     688                          MERGE( nest_offl%nzu, 0, bc_dirichlet_s ),           & ! number of time steps (2 or 0)
     689                          MERGE( 2, 0, bc_dirichlet_s ),                       & ! parallel IO when compiled accordingly
     690                          .TRUE. )                                             
     691                                                                               
     692       CALL get_variable( pids_id, 'ls_forcing_south_w',                       & ! array to be read
     693                          nest_offl%w_south,                                   & ! start index x direction
     694                          MERGE( nxl+1, 1, bc_dirichlet_s ),                   & ! start index z direction
     695                          MERGE( nzb+1, 1, bc_dirichlet_s ),                   & ! start index time dimension
     696                          MERGE( nest_offl%tind+1, 1, bc_dirichlet_s ),        & ! number of elements along x
     697                          MERGE( nxr-nxl+1, 0, bc_dirichlet_s ),               & ! number of elements alogn z
     698                          MERGE( nest_offl%nzw, 0, bc_dirichlet_s ),           & ! number of time steps (2 or 0)
     699                          MERGE( 2, 0, bc_dirichlet_s ),                       & ! parallel IO when compiled accordingly
     700                          .TRUE. )                                             
     701                                                                               
     702       IF ( .NOT. neutral )  THEN                                             
     703          CALL get_variable( pids_id, 'ls_forcing_south_pt',                   & ! array to be read
     704                             nest_offl%pt_south,                               & ! start index x direction
     705                             MERGE( nxl+1, 1, bc_dirichlet_s ),                & ! start index z direction
     706                             MERGE( nzb+1, 1, bc_dirichlet_s ),                & ! start index time dimension
     707                             MERGE( nest_offl%tind+1, 1, bc_dirichlet_s ),     & ! number of elements along x
     708                             MERGE( nxr-nxl+1, 0, bc_dirichlet_s ),            & ! number of elements alogn z
     709                             MERGE( nest_offl%nzu, 0, bc_dirichlet_s ),        & ! number of time steps (2 or 0)
     710                             MERGE( 2, 0, bc_dirichlet_s ),                    & ! parallel IO when compiled accordingly
     711                             .TRUE. )                                             
     712       ENDIF                                                                   
     713       IF ( humidity )  THEN                                                   
     714          CALL get_variable( pids_id, 'ls_forcing_south_qv',                   & ! array to be read
     715                             nest_offl%q_south,                                & ! start index x direction
     716                             MERGE( nxl+1, 1, bc_dirichlet_s ),                & ! start index z direction
     717                             MERGE( nzb+1, 1, bc_dirichlet_s ),                & ! start index time dimension
     718                             MERGE( nest_offl%tind+1, 1, bc_dirichlet_s ),     & ! number of elements along x
     719                             MERGE( nxr-nxl+1, 0, bc_dirichlet_s ),            & ! number of elements alogn z
     720                             MERGE( nest_offl%nzu, 0, bc_dirichlet_s ),        & ! number of time steps (2 or 0)
     721                             MERGE( 2, 0, bc_dirichlet_s ),                    & ! parallel IO when compiled accordingly
     722                             .TRUE. )                                             
     723       ENDIF                                                                   
     724                                                                               
     725       IF ( air_chemistry )  THEN                                             
     726          DO  n = 1, UBOUND(nest_offl%var_names_chem_s, 1)                     
     727             IF ( check_existence( nest_offl%var_names,                        &
     728                                   nest_offl%var_names_chem_s(n) ) )  THEN     
     729                CALL get_variable( pids_id,                                    &
     730                           TRIM( nest_offl%var_names_chem_s(n) ),              &
     731                           nest_offl%chem_south(:,:,:,n),                      &
     732                           MERGE( nxl+1, 1, bc_dirichlet_s ),                  &
     733                           MERGE( nzb+1, 1, bc_dirichlet_s ),                  &
     734                           MERGE( nest_offl%tind+1, 1, bc_dirichlet_s ),       &
     735                           MERGE( nxr-nxl+1, 0, bc_dirichlet_s ),              &
     736                           MERGE( nest_offl%nzu, 0, bc_dirichlet_s ),          &
     737                           MERGE( 2, 0, bc_dirichlet_s ),                      &
     738                           .TRUE. )
     739                nest_offl%chem_from_file_s(n) = .TRUE.
     740             ENDIF
     741          ENDDO
     742       ENDIF
     743!
     744!--    Top boundary
     745       CALL get_variable( pids_id, 'ls_forcing_top_u',                         &
     746                             nest_offl%u_top(0:1,nys:nyn,nxlu:nxr),            &
     747                             nxlu, nys+1, nest_offl%tind+1,                    &
     748                             nxr-nxlu+1, nyn-nys+1, 2, .TRUE. )
     749
     750       CALL get_variable( pids_id, 'ls_forcing_top_v',                         &
     751                             nest_offl%v_top(0:1,nysv:nyn,nxl:nxr),            &
     752                             nxl+1, nysv, nest_offl%tind+1,                    &
     753                             nxr-nxl+1, nyn-nysv+1, 2, .TRUE. )
     754                             
     755       CALL get_variable( pids_id, 'ls_forcing_top_w',                         &
     756                             nest_offl%w_top(0:1,nys:nyn,nxl:nxr),             &
     757                             nxl+1, nys+1, nest_offl%tind+1,                   &
     758                             nxr-nxl+1, nyn-nys+1, 2, .TRUE. )
     759                             
     760       IF ( .NOT. neutral )  THEN
     761          CALL get_variable( pids_id, 'ls_forcing_top_pt',                     &
     762                                nest_offl%pt_top(0:1,nys:nyn,nxl:nxr),         &
     763                                nxl+1, nys+1, nest_offl%tind+1,                &
     764                                nxr-nxl+1, nyn-nys+1, 2, .TRUE. )
     765       ENDIF
     766       IF ( humidity )  THEN
     767          CALL get_variable( pids_id, 'ls_forcing_top_qv',                     &
     768                                nest_offl%q_top(0:1,nys:nyn,nxl:nxr),          &
     769                                nxl+1, nys+1, nest_offl%tind+1,                &
     770                                nxr-nxl+1, nyn-nys+1, 2, .TRUE. )
     771       ENDIF
     772       
     773       IF ( air_chemistry )  THEN
     774          DO  n = 1, UBOUND(nest_offl%var_names_chem_t, 1)
     775             IF ( check_existence( nest_offl%var_names,                        &
     776                                   nest_offl%var_names_chem_t(n) ) )  THEN     
     777                CALL get_variable( pids_id,                                    &
     778                              TRIM( nest_offl%var_names_chem_t(n) ),           &
     779                              nest_offl%chem_top(0:1,nys:nyn,nxl:nxr,n),       &
     780                              nxl+1, nys+1, nest_offl%tind+1,                  &
     781                              nxr-nxl+1, nyn-nys+1, 2, .TRUE. )
     782                nest_offl%chem_from_file_t(n) = .TRUE.
     783             ENDIF
     784          ENDDO
     785       ENDIF
     786
     787!
     788!--    Close input file
     789       CALL close_input_file( pids_id )
     790#endif
     791!
     792!--    Set control flag to indicate that boundary data has been initially
     793!--    input.
     794       nest_offl%init = .TRUE.
     795!
     796!--    End of CPU measurement
     797       CALL cpu_log( log_point_s(86), 'NetCDF input forcing', 'stop' )
     798
     799    END SUBROUTINE nesting_offl_input
     800
    199801
    200802!------------------------------------------------------------------------------!
     
    210812    SUBROUTINE nesting_offl_mass_conservation
    211813
    212        IMPLICIT NONE
    213 
    214814       INTEGER(iwp) ::  i !< grid index in x-direction
    215815       INTEGER(iwp) ::  j !< grid index in y-direction
     
    281881#if defined( __parallel )
    282882       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
    283        CALL MPI_ALLREDUCE( volume_flow_l, volume_flow, 3, MPI_REAL, MPI_SUM,      &
     883       CALL MPI_ALLREDUCE( volume_flow_l, volume_flow, 3, MPI_REAL, MPI_SUM,   &
    284884                           comm2d, ierr )
    285885#else
     
    314914    SUBROUTINE nesting_offl_bc                     
    315915
    316        IMPLICIT NONE
    317 
    318916       INTEGER(iwp) ::  i !< running index x-direction
    319917       INTEGER(iwp) ::  j !< running index y-direction
    320918       INTEGER(iwp) ::  k !< running index z-direction
    321919       INTEGER(iwp) ::  n !< running index for chemical species
    322 
    323        REAL(wp) ::  fac_dt   !< interpolation factor
    324920       
    325921       REAL(wp), DIMENSION(nzb:nzt+1) ::  pt_ref   !< reference profile for potential temperature
     
    347943       u_ref_l  = 0.0_wp
    348944       v_ref_l  = 0.0_wp
    349 !
    350 !--    Determine interpolation factor and limit it to 1. This is because
    351 !--    t+dt can slightly exceed time(tind_p) before boundary data is updated
    352 !--    again.
    353        fac_dt = ( time_since_reference_point - nest_offl%time(nest_offl%tind)  &
    354                 + dt_3d ) /                                                    &
    355            ( nest_offl%time(nest_offl%tind_p) - nest_offl%time(nest_offl%tind) )
    356        fac_dt = MIN( 1.0_wp, fac_dt )
    357945!
    358946!--    Set boundary conditions of u-, v-, w-component, as well as q, and pt.
     
    8631451       CALL nesting_offl_calc_zi
    8641452       CALL adjust_sponge_layer
    865    
     1453       
     1454       CALL  cpu_log( log_point(58), 'offline nesting', 'stop' )
     1455
     1456       IF ( debug_output_timestep )  CALL debug_message( 'nesting_offl_bc', 'end' )
     1457
     1458
     1459    END SUBROUTINE nesting_offl_bc
     1460
     1461!------------------------------------------------------------------------------!
     1462! Description:
     1463!------------------------------------------------------------------------------!
     1464!>  Update of the geostrophic wind components.
     1465!>  @todo: update geostrophic wind also in the child domains (should be done
     1466!>         in the nesting.
     1467!------------------------------------------------------------------------------!
     1468    SUBROUTINE nesting_offl_geostrophic_wind
     1469
     1470       INTEGER(iwp) ::  k
    8661471!
    8671472!--    Update geostrophic wind components from dynamic input file.
     
    8741479       ug(nzt+1) = ug(nzt)
    8751480       vg(nzt+1) = vg(nzt)
    876    
    877        CALL  cpu_log( log_point(58), 'offline nesting', 'stop' )
    878 
    879        IF ( debug_output_timestep )  CALL debug_message( 'nesting_offl_bc', 'end' )
    880 
    881 
    882     END SUBROUTINE nesting_offl_bc
     1481
     1482    END SUBROUTINE nesting_offl_geostrophic_wind
     1483
     1484!------------------------------------------------------------------------------!
     1485! Description:
     1486!------------------------------------------------------------------------------!
     1487!>  Determine the interpolation constant for time interpolation. The
     1488!>  calculation is separated from the nesting_offl_bc and
     1489!>  nesting_offl_geostrophic_wind in order to be independent on the order
     1490!>  of calls.
     1491!------------------------------------------------------------------------------!
     1492    SUBROUTINE nesting_offl_interpolation_factor
     1493!
     1494!--    Determine interpolation factor and limit it to 1. This is because
     1495!--    t+dt can slightly exceed time(tind_p) before boundary data is updated
     1496!--    again.
     1497       fac_dt = ( time_utc_init + time_since_reference_point                   &
     1498                - nest_offl%time(nest_offl%tind) + dt_3d ) /                   &
     1499           ( nest_offl%time(nest_offl%tind_p) - nest_offl%time(nest_offl%tind) )
     1500
     1501       fac_dt = MIN( 1.0_wp, fac_dt )
     1502
     1503    END SUBROUTINE nesting_offl_interpolation_factor
    8831504
    8841505!------------------------------------------------------------------------------!
     
    8891510!------------------------------------------------------------------------------!
    8901511    SUBROUTINE nesting_offl_calc_zi
    891        
    892        USE basic_constants_and_equations_mod,                                  &
    893            ONLY:  g
    894        
    895        USE kinds
    896 
    897        IMPLICIT NONE
    8981512
    8991513       INTEGER(iwp) :: i                            !< loop index in x-direction
     
    10701684!------------------------------------------------------------------------------!
    10711685    SUBROUTINE adjust_sponge_layer
    1072        
    1073        USE arrays_3d,                                                          &
    1074            ONLY:  rdf, rdf_sc, zu
    1075        
    1076        USE basic_constants_and_equations_mod,                                  &
    1077            ONLY:  pi
    1078        
    1079        USE kinds
    1080 
    1081        IMPLICIT NONE
    10821686
    10831687       INTEGER(iwp) :: k   !< loop index in z-direction
     
    11161720!------------------------------------------------------------------------------!
    11171721    SUBROUTINE nesting_offl_check_parameters
    1118 
    1119        IMPLICIT NONE
    11201722!
    11211723!--    Check if offline nesting is applied in nested child domain.
     
    11231725          message_string = 'Offline nesting is only applicable in root model.'
    11241726          CALL message( 'offline_nesting_check_parameters', 'PA0622', 1, 2, 0, 6, 0 )       
    1125        ENDIF     
     1727       ENDIF
    11261728
    11271729    END SUBROUTINE nesting_offl_check_parameters
     
    11331735!------------------------------------------------------------------------------!
    11341736    SUBROUTINE nesting_offl_parin
    1135 
    1136        IMPLICIT NONE
    11371737       
    11381738       CHARACTER (LEN=80) ::  line   !< dummy string that contains the current line of the parameter file
     
    11741774    SUBROUTINE nesting_offl_header ( io )
    11751775
    1176        IMPLICIT NONE
    1177 
    11781776       INTEGER(iwp), INTENT(IN) ::  io !< Unit of the output file
    11791777
     
    11991797!------------------------------------------------------------------------------!
    12001798    SUBROUTINE nesting_offl_init
    1201    
    1202        USE netcdf_data_input_mod,                                              &
    1203            ONLY:  netcdf_data_input_offline_nesting 
    1204 
    1205        IMPLICIT NONE
    1206        
     1799           
    12071800       INTEGER(iwp) ::  n !< running index for chemical species
    12081801
     
    13401933       ENDIF
    13411934!
     1935!--    Before initial data input is initiated, check if dynamic input file is
     1936!--    present.
     1937       IF ( .NOT. input_pids_dynamic )  THEN
     1938          message_string = 'nesting_offline = .TRUE. requires dynamic '  //    &
     1939                            'input file ' //                                   &
     1940                            TRIM( input_file_dynamic ) // TRIM( coupling_char )
     1941          CALL message( 'nesting_offl_init', 'PA0546', 1, 2, 0, 6, 0 )
     1942       ENDIF
     1943!
    13421944!--    Read COSMO data at lateral and top boundaries
    1343        CALL netcdf_data_input_offline_nesting
     1945       CALL nesting_offl_input
    13441946!
    13451947!--    Check if sufficient time steps are provided to cover the entire
     
    13471949!--    not for the soil/wall spinup. However, as the spinup time is added
    13481950!--    to the end_time, this must be considered here.
    1349        IF ( end_time - spinup_time > nest_offl%time(nest_offl%nt-1) )  THEN
    1350           message_string = 'end_time > provided time in offline nesting.'
    1351           CALL message( 'offline_nesting_check_parameters', 'PA0183',          &
    1352                         1, 2, 0, 6, 0 )
     1951       IF ( end_time - spinup_time >                                           &
     1952            nest_offl%time(nest_offl%nt-1) - time_utc_init )  THEN
     1953          message_string = 'end_time of the simulation exceeds the ' //        &
     1954                           'time dimension in the dynamic input file.'
     1955          CALL message( 'nesting_offl_init', 'PA0183', 1, 2, 0, 6, 0 )
     1956       ENDIF
     1957
     1958       IF ( nest_offl%time(0) /= time_utc_init )  THEN
     1959          message_string = 'Offline nesting: time dimension must start at ' // &
     1960                           ' time_utc_init.'
     1961          CALL message( 'nesting_offl_init', 'PA0676', 1, 2, 0, 6, 0 )
    13531962       ENDIF
    13541963!
     
    14662075!------------------------------------------------------------------------------!
    14672076    FUNCTION interpolate_in_time( var_t1, var_t2, fac  )
    1468        
    1469        USE kinds
    1470 
    1471        IMPLICIT NONE
    14722077
    14732078       REAL(wp)            :: interpolate_in_time !< time-interpolated boundary value
Note: See TracChangeset for help on using the changeset viewer.