#if defined( __ibmy_special ) @PROCESS NOOPTimize #endif SUBROUTINE init_3d_model !------------------------------------------------------------------------------! ! Current revisions: ! ------------------ ! ! ! Former revisions: ! ----------------- ! $Id: init_3d_model.f90 760 2011-09-15 14:37:54Z suehring $ ! ! 759 2011-09-15 13:58:31Z raasch ! Splitting of parallel I/O in blocks of PEs ! Bugfix: No zero assignments to volume_flow_initial and volume_flow_area in ! case of normal restart runs. ! ! 713 2011-03-30 14:21:21Z suehring ! weight_substep and weight_pres are given as fractions. ! ! 709 2011-03-30 09:31:40Z raasch ! formatting adjustments ! ! 707 2011-03-29 11:39:40Z raasch ! p_sub renamed p_loc and allocated depending on the chosen pressure solver, ! initial assignments of zero to array p for iterative solvers only, ! bc_lr/ns replaced by bc_lr/ns_dirrad/raddir ! ! 680 2011-02-04 23:16:06Z gryschka ! bugfix: volume_flow_control ! ! 673 2011-01-18 16:19:48Z suehring ! weight_substep (moved from advec_ws) and weight_pres added. ! Allocate p_sub when using Multigrid or SOR solver. ! Call of ws_init moved behind the if requests. ! ! 667 2010-12-23 12:06:00Z suehring/gryschka ! nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng in loops and ! allocation of arrays. Calls of exchange_horiz are modified. ! Call ws_init to initialize arrays needed for calculating statisticas and for ! optimization when ws-scheme is used. ! Initial volume flow is now calculated by using the variable hom_sum. ! Therefore the correction of initial volume flow for non-flat topography ! removed (removed u_nzb_p1_for_vfc and v_nzb_p1_for_vfc) ! Changed surface boundary conditions for u and v in case of ibc_uv_b == 0 from ! mirror to Dirichlet boundary conditions (u=v=0), so that k=nzb is ! representative for the height z0. ! Bugfix: type conversion of '1' to 64bit for the MAX function (ngp_3d_inner) ! ! 622 2010-12-10 08:08:13Z raasch ! optional barriers included in order to speed up collective operations ! ! 560 2010-09-09 10:06:09Z weinreis ! bugfix: correction of calculating ngp_3d for 64 bit ! ! 485 2010-02-05 10:57:51Z raasch ! calculation of ngp_3d + ngp_3d_inner changed because they have now 64 bit ! ! 407 2009-12-01 15:01:15Z maronga ! var_ts is replaced by dots_max ! Enabled passive scalar/humidity wall fluxes for non-flat topography ! ! 388 2009-09-23 09:40:33Z raasch ! Initialization of prho added. ! bugfix: correction of initial volume flow for non-flat topography ! bugfix: zero initialization of arrays within buildings for 'cyclic_fill' ! bugfix: avoid that ngp_2dh_s_inner becomes zero ! initializing_actions='read_data_for_recycling' renamed to 'cyclic_fill', now ! independent of turbulent_inflow ! Output of messages replaced by message handling routine. ! Set the starting level and the vertical smoothing factor used for ! the external pressure gradient ! +conserve_volume_flow_mode: 'default', 'initial_profiles', 'inflow_profile' ! and 'bulk_velocity' ! If the inversion height calculated by the prerun is zero, ! inflow_damping_height must be explicitly specified. ! ! 181 2008-07-30 07:07:47Z raasch ! bugfix: zero assignments to tendency arrays in case of restarts, ! further extensions and modifications in the initialisation of the plant ! canopy model, ! allocation of hom_sum moved to parin, initialization of spectrum_x|y directly ! after allocating theses arrays, ! read data for recycling added as new initialization option, ! dummy allocation for diss ! ! 138 2007-11-28 10:03:58Z letzel ! New counter ngp_2dh_s_inner. ! Allow new case bc_uv_t = 'dirichlet_0' for channel flow. ! Corrected calculation of initial volume flow for 'set_1d-model_profiles' and ! 'set_constant_profiles' in case of buildings in the reference cross-sections. ! ! 108 2007-08-24 15:10:38Z letzel ! Flux initialization in case of coupled runs, +momentum fluxes at top boundary, ! +arrays for phase speed c_u, c_v, c_w, indices for u|v|w_m_l|r changed ! +qswst_remote in case of atmosphere model with humidity coupled to ocean ! Rayleigh damping for ocean, optionally calculate km and kh from initial ! TKE e_init ! ! 97 2007-06-21 08:23:15Z raasch ! Initialization of salinity, call of init_ocean ! ! 87 2007-05-22 15:46:47Z raasch ! var_hom and var_sum renamed pr_palm ! ! 75 2007-03-22 09:54:05Z raasch ! Arrays for radiation boundary conditions are allocated (u_m_l, u_m_r, etc.), ! bugfix for cases with the outflow damping layer extending over more than one ! subdomain, moisture renamed humidity, ! new initializing action "by_user" calls user_init_3d_model, ! precipitation_amount/rate, ts_value are allocated, +module netcdf_control, ! initial velocities at nzb+1 are regarded for volume ! flow control in case they have been set zero before (to avoid small timesteps) ! -uvmean_outflow, uxrp, vynp eliminated ! ! 19 2007-02-23 04:53:48Z raasch ! +handling of top fluxes ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.49 2006/08/22 15:59:07 raasch ! No optimization of this file on the ibmy (Yonsei Univ.) ! ! Revision 1.1 1998/03/09 16:22:22 raasch ! Initial revision ! ! ! Description: ! ------------ ! Allocation of arrays and initialization of the 3D model via ! a) pre-run the 1D model ! or ! b) pre-set constant linear profiles ! or ! c) read values of a previous run !------------------------------------------------------------------------------! USE advec_ws USE arrays_3d USE averaging USE cloud_parameters USE constants USE control_parameters USE cpulog USE indices USE interfaces USE model_1d USE netcdf_control USE particle_attributes USE pegrid USE profil_parameter USE random_function_mod USE statistics IMPLICIT NONE INTEGER :: i, ind_array(1), j, k, sr INTEGER, DIMENSION(:), ALLOCATABLE :: ngp_2dh_l INTEGER, DIMENSION(:,:), ALLOCATABLE :: ngp_2dh_outer_l, & ngp_2dh_s_inner_l REAL :: a, b REAL, DIMENSION(1:2) :: volume_flow_area_l, volume_flow_initial_l REAL, DIMENSION(:), ALLOCATABLE :: ngp_3d_inner_l, ngp_3d_inner_tmp ! !-- Allocate arrays ALLOCATE( ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions), & ngp_3d(0:statistic_regions), & ngp_3d_inner(0:statistic_regions), & ngp_3d_inner_l(0:statistic_regions), & ngp_3d_inner_tmp(0:statistic_regions), & sums_divnew_l(0:statistic_regions), & sums_divold_l(0:statistic_regions) ) ALLOCATE( dp_smooth_factor(nzb:nzt), rdf(nzb+1:nzt) ) ALLOCATE( ngp_2dh_outer(nzb:nzt+1,0:statistic_regions), & ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions), & ngp_2dh_s_inner(nzb:nzt+1,0:statistic_regions), & ngp_2dh_s_inner_l(nzb:nzt+1,0:statistic_regions), & rmask(nysg:nyng,nxlg:nxrg,0:statistic_regions), & sums(nzb:nzt+1,pr_palm+max_pr_user), & sums_l(nzb:nzt+1,pr_palm+max_pr_user,0:threads_per_task-1), & sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1), & sums_up_fraction_l(10,3,0:statistic_regions), & sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions), & ts_value(dots_max,0:statistic_regions) ) ALLOCATE( km_damp_x(nxlg:nxrg), km_damp_y(nysg:nyng) ) ALLOCATE( rif_1(nysg:nyng,nxlg:nxrg), shf_1(nysg:nyng,nxlg:nxrg), & ts(nysg:nyng,nxlg:nxrg), tswst_1(nysg:nyng,nxlg:nxrg), & us(nysg:nyng,nxlg:nxrg), usws_1(nysg:nyng,nxlg:nxrg), & uswst_1(nysg:nyng,nxlg:nxrg), & vsws_1(nysg:nyng,nxlg:nxrg), & vswst_1(nysg:nyng,nxlg:nxrg), z0(nysg:nyng,nxlg:nxrg) ) IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Leapfrog scheme needs two timelevels of diffusion quantities ALLOCATE( rif_2(nysg:nyng,nxlg:nxrg), & shf_2(nysg:nyng,nxlg:nxrg), & tswst_2(nysg:nyng,nxlg:nxrg), & usws_2(nysg:nyng,nxlg:nxrg), & uswst_2(nysg:nyng,nxlg:nxrg), & vswst_2(nysg:nyng,nxlg:nxrg), & vsws_2(nysg:nyng,nxlg:nxrg) ) ENDIF ALLOCATE( d(nzb+1:nzta,nys:nyna,nxl:nxra), & e_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & e_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & e_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & kh_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & km_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & p(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & pt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & pt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & pt_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & u_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & u_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & u_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & v_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & v_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & v_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & w_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & w_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & w_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ! !-- Following array is required for perturbation pressure within the iterative !-- pressure solvers. For the multistep schemes (Runge-Kutta), array p holds !-- the weighted average of the substeps and cannot be used in the Poisson !-- solver. IF ( psolver == 'sor' ) THEN ALLOCATE( p_loc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ELSEIF ( psolver == 'multigrid' ) THEN ! !-- For performance reasons, multigrid is using one ghost layer only ALLOCATE( p_loc(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) ) ENDIF IF ( timestep_scheme(1:5) /= 'runge' ) THEN ALLOCATE( kh_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & km_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF IF ( humidity .OR. passive_scalar ) THEN ! !-- 2D-humidity/scalar arrays ALLOCATE ( qs(nysg:nyng,nxlg:nxrg), & qsws_1(nysg:nyng,nxlg:nxrg), & qswst_1(nysg:nyng,nxlg:nxrg) ) IF ( timestep_scheme(1:5) /= 'runge' ) THEN ALLOCATE( qsws_2(nysg:nyng,nxlg:nxrg), & qswst_2(nysg:nyng,nxlg:nxrg) ) ENDIF ! !-- 3D-humidity/scalar arrays ALLOCATE( q_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & q_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & q_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ! !-- 3D-arrays needed for humidity only IF ( humidity ) THEN ALLOCATE( vpt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) IF ( timestep_scheme(1:5) /= 'runge' ) THEN ALLOCATE( vpt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF IF ( cloud_physics ) THEN ! !-- Liquid water content ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ! !-- Precipitation amount and rate (only needed if output is switched) ALLOCATE( precipitation_amount(nysg:nyng,nxlg:nxrg), & precipitation_rate(nysg:nyng,nxlg:nxrg) ) ENDIF IF ( cloud_droplets ) THEN ! !-- Liquid water content, change in liquid water content, !-- real volume of particles (with weighting), volume of particles ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & ql_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & ql_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & ql_vp(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF ENDIF ENDIF IF ( ocean ) THEN ALLOCATE( saswsb_1(nysg:nyng,nxlg:nxrg), & saswst_1(nysg:nyng,nxlg:nxrg) ) ALLOCATE( prho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & rho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & sa_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & sa_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & sa_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) prho => prho_1 rho => rho_1 ! routines calc_mean_profile and diffusion_e require ! density to be apointer IF ( humidity_remote ) THEN ALLOCATE( qswst_remote(nysg:nyng,nxlg:nxrg)) qswst_remote = 0.0 ENDIF ENDIF ! !-- 3D-array for storing the dissipation, needed for calculating the sgs !-- particle velocities IF ( use_sgs_for_particles ) THEN ALLOCATE ( diss(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ELSE ALLOCATE ( diss(2,2,2) ) ! required because diss is used as a ! formal parameter ENDIF IF ( dt_dosp /= 9999999.9 ) THEN ALLOCATE( spectrum_x( 1:nx/2, 1:10, 1:10 ), & spectrum_y( 1:ny/2, 1:10, 1:10 ) ) spectrum_x = 0.0 spectrum_y = 0.0 ENDIF ! !-- 3D-arrays for the leaf area density and the canopy drag coefficient IF ( plant_canopy ) THEN ALLOCATE ( lad_s(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & lad_u(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & lad_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & lad_w(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & cdc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) IF ( passive_scalar ) THEN ALLOCATE ( sls(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & sec(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF IF ( cthf /= 0.0 ) THEN ALLOCATE ( lai(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & canopy_heat_flux(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF ENDIF ! !-- 4D-array for storing the Rif-values at vertical walls IF ( topography /= 'flat' ) THEN ALLOCATE( rif_wall(nzb:nzt+1,nysg:nyng,nxlg:nxrg,1:4) ) rif_wall = 0.0 ENDIF ! !-- Arrays to store velocity data from t-dt and the phase speeds which !-- are needed for radiation boundary conditions IF ( outflow_l ) THEN ALLOCATE( u_m_l(nzb:nzt+1,nysg:nyng,1:2), & v_m_l(nzb:nzt+1,nysg:nyng,0:1), & w_m_l(nzb:nzt+1,nysg:nyng,0:1) ) ENDIF IF ( outflow_r ) THEN ALLOCATE( u_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx), & v_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx), & w_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx) ) ENDIF IF ( outflow_l .OR. outflow_r ) THEN ALLOCATE( c_u(nzb:nzt+1,nysg:nyng), c_v(nzb:nzt+1,nysg:nyng), & c_w(nzb:nzt+1,nysg:nyng) ) ENDIF IF ( outflow_s ) THEN ALLOCATE( u_m_s(nzb:nzt+1,0:1,nxlg:nxrg), & v_m_s(nzb:nzt+1,1:2,nxlg:nxrg), & w_m_s(nzb:nzt+1,0:1,nxlg:nxrg) ) ENDIF IF ( outflow_n ) THEN ALLOCATE( u_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg), & v_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg), & w_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg) ) ENDIF IF ( outflow_s .OR. outflow_n ) THEN ALLOCATE( c_u(nzb:nzt+1,nxlg:nxrg), c_v(nzb:nzt+1,nxlg:nxrg), & c_w(nzb:nzt+1,nxlg:nxrg) ) ENDIF ! !-- Initial assignment of the pointers IF ( timestep_scheme(1:5) /= 'runge' ) THEN rif_m => rif_1; rif => rif_2 shf_m => shf_1; shf => shf_2 tswst_m => tswst_1; tswst => tswst_2 usws_m => usws_1; usws => usws_2 uswst_m => uswst_1; uswst => uswst_2 vsws_m => vsws_1; vsws => vsws_2 vswst_m => vswst_1; vswst => vswst_2 e_m => e_1; e => e_2; e_p => e_3; te_m => e_3 kh_m => kh_1; kh => kh_2 km_m => km_1; km => km_2 pt_m => pt_1; pt => pt_2; pt_p => pt_3; tpt_m => pt_3 u_m => u_1; u => u_2; u_p => u_3; tu_m => u_3 v_m => v_1; v => v_2; v_p => v_3; tv_m => v_3 w_m => w_1; w => w_2; w_p => w_3; tw_m => w_3 IF ( humidity .OR. passive_scalar ) THEN qsws_m => qsws_1; qsws => qsws_2 qswst_m => qswst_1; qswst => qswst_2 q_m => q_1; q => q_2; q_p => q_3; tq_m => q_3 IF ( humidity ) vpt_m => vpt_1; vpt => vpt_2 IF ( cloud_physics ) ql => ql_1 IF ( cloud_droplets ) THEN ql => ql_1 ql_c => ql_2 ENDIF ENDIF ELSE rif => rif_1 shf => shf_1 tswst => tswst_1 usws => usws_1 uswst => uswst_1 vsws => vsws_1 vswst => vswst_1 e => e_1; e_p => e_2; te_m => e_3; e_m => e_3 kh => kh_1 km => km_1 pt => pt_1; pt_p => pt_2; tpt_m => pt_3; pt_m => pt_3 u => u_1; u_p => u_2; tu_m => u_3; u_m => u_3 v => v_1; v_p => v_2; tv_m => v_3; v_m => v_3 w => w_1; w_p => w_2; tw_m => w_3; w_m => w_3 IF ( humidity .OR. passive_scalar ) THEN qsws => qsws_1 qswst => qswst_1 q => q_1; q_p => q_2; tq_m => q_3; q_m => q_3 IF ( humidity ) vpt => vpt_1 IF ( cloud_physics ) ql => ql_1 IF ( cloud_droplets ) THEN ql => ql_1 ql_c => ql_2 ENDIF ENDIF IF ( ocean ) THEN saswsb => saswsb_1 saswst => saswst_1 sa => sa_1; sa_p => sa_2; tsa_m => sa_3 ENDIF ENDIF ! !-- Allocate arrays containing the RK coefficient for calculation of !-- perturbation pressure and turbulent fluxes. At this point values are !-- set for pressure calculation during initialization (where no timestep !-- is done). Further below the values needed within the timestep scheme !-- will be set. ALLOCATE( weight_substep(1:intermediate_timestep_count_max), & weight_pres(1:intermediate_timestep_count_max) ) weight_substep = 1.0 weight_pres = 1.0 intermediate_timestep_count = 1 ! needed when simulated_time = 0.0 ! !-- Initialize model variables IF ( TRIM( initializing_actions ) /= 'read_restart_data' .AND. & TRIM( initializing_actions ) /= 'cyclic_fill' ) THEN ! !-- First model run of a possible job queue. !-- Initial profiles of the variables must be computes. IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 ) THEN ! !-- Use solutions of the 1D model as initial profiles, !-- start 1D model CALL init_1d_model ! !-- Transfer initial profiles to the arrays of the 3D model DO i = nxlg, nxrg DO j = nysg, nyng e(:,j,i) = e1d kh(:,j,i) = kh1d km(:,j,i) = km1d pt(:,j,i) = pt_init u(:,j,i) = u1d v(:,j,i) = v1d ENDDO ENDDO IF ( humidity .OR. passive_scalar ) THEN DO i = nxlg, nxrg DO j = nysg, nyng q(:,j,i) = q_init ENDDO ENDDO ENDIF IF ( .NOT. constant_diffusion ) THEN DO i = nxlg, nxrg DO j = nysg, nyng e(:,j,i) = e1d ENDDO ENDDO ! !-- Store initial profiles for output purposes etc. hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 ) IF ( prandtl_layer ) THEN rif = rif1d(nzb+1) ts = 0.0 ! could actually be computed more accurately in the ! 1D model. Update when opportunity arises. us = us1d usws = usws1d vsws = vsws1d ELSE ts = 0.0 ! must be set, because used in rif = 0.0 ! flowste us = 0.0 usws = 0.0 vsws = 0.0 ENDIF ELSE e = 0.0 ! must be set, because used in rif = 0.0 ! flowste ts = 0.0 us = 0.0 usws = 0.0 vsws = 0.0 ENDIF uswst = top_momentumflux_u vswst = top_momentumflux_v ! !-- In every case qs = 0.0 (see also pt) !-- This could actually be computed more accurately in the 1D model. !-- Update when opportunity arises! IF ( humidity .OR. passive_scalar ) qs = 0.0 ! !-- inside buildings set velocities back to zero IF ( topography /= 'flat' ) THEN DO i = nxl-1, nxr+1 DO j = nys-1, nyn+1 u(nzb:nzb_u_inner(j,i),j,i) = 0.0 v(nzb:nzb_v_inner(j,i),j,i) = 0.0 ENDDO ENDDO ! !-- WARNING: The extra boundary conditions set after running the !-- ------- 1D model impose an error on the divergence one layer !-- below the topography; need to correct later !-- ATTENTION: Provisional correction for Piacsek & Williams !-- --------- advection scheme: keep u and v zero one layer below !-- the topography. ! !-- Following was removed, because mirror boundary condition are !-- replaced by dirichlet boundary conditions ! ! IF ( ibc_uv_b == 0 ) THEN !! !!-- Satisfying the Dirichlet condition with an extra layer below !!-- the surface where the u and v component change their sign. ! DO i = nxl-1, nxr+1 ! DO j = nys-1, nyn+1 ! IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = -u(1,j,i) ! IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = -v(1,j,i) ! ENDDO ! ENDDO ! ! ELSE IF ( ibc_uv_b == 1 ) THEN ! !-- Neumann condition DO i = nxl-1, nxr+1 DO j = nys-1, nyn+1 IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = u(1,j,i) IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = v(1,j,i) ENDDO ENDDO ENDIF ENDIF ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 ) & THEN ! !-- Use constructed initial profiles (velocity constant with height, !-- temperature profile with constant gradient) DO i = nxlg, nxrg DO j = nysg, nyng pt(:,j,i) = pt_init u(:,j,i) = u_init v(:,j,i) = v_init ENDDO ENDDO ! !-- Set initial horizontal velocities at the lowest computational grid !-- levels to zero in order to avoid too small time steps caused by the !-- diffusion limit in the initial phase of a run (at k=1, dz/2 occurs !-- in the limiting formula!). The original values are stored to be later !-- used for volume flow control. DO i = nxlg, nxrg DO j = nysg, nyng u(nzb:nzb_u_inner(j,i)+1,j,i) = 0.0 v(nzb:nzb_v_inner(j,i)+1,j,i) = 0.0 ENDDO ENDDO IF ( humidity .OR. passive_scalar ) THEN DO i = nxlg, nxrg DO j = nysg, nyng q(:,j,i) = q_init ENDDO ENDDO ENDIF IF ( ocean ) THEN DO i = nxlg, nxrg DO j = nysg, nyng sa(:,j,i) = sa_init ENDDO ENDDO ENDIF IF ( constant_diffusion ) THEN km = km_constant kh = km / prandtl_number e = 0.0 ELSEIF ( e_init > 0.0 ) THEN DO k = nzb+1, nzt km(k,:,:) = 0.1 * l_grid(k) * SQRT( e_init ) ENDDO km(nzb,:,:) = km(nzb+1,:,:) km(nzt+1,:,:) = km(nzt,:,:) kh = km / prandtl_number e = e_init ELSE IF ( .NOT. ocean ) THEN kh = 0.01 ! there must exist an initial diffusion, because km = 0.01 ! otherwise no TKE would be produced by the ! production terms, as long as not yet ! e = (u*/cm)**2 at k=nzb+1 ELSE kh = 0.00001 km = 0.00001 ENDIF e = 0.0 ENDIF rif = 0.0 ts = 0.0 us = 0.0 usws = 0.0 uswst = top_momentumflux_u vsws = 0.0 vswst = top_momentumflux_v IF ( humidity .OR. passive_scalar ) qs = 0.0 ! !-- Compute initial temperature field and other constants used in case !-- of a sloping surface IF ( sloping_surface ) CALL init_slope ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 ) & THEN ! !-- Initialization will completely be done by the user CALL user_init_3d_model ENDIF ! !-- Bottom boundary IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2 ) THEN u(nzb,:,:) = 0.0 v(nzb,:,:) = 0.0 ENDIF ! !-- Apply channel flow boundary condition IF ( TRIM( bc_uv_t ) == 'dirichlet_0' ) THEN u(nzt+1,:,:) = 0.0 v(nzt+1,:,:) = 0.0 !-- For the Dirichlet condition to be correctly applied at the top, set !-- ug and vg to zero there ug(nzt+1) = 0.0 vg(nzt+1) = 0.0 ENDIF ! !-- Calculate virtual potential temperature IF ( humidity ) vpt = pt * ( 1.0 + 0.61 * q ) ! !-- Store initial profiles for output purposes etc. hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 ) hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 ) IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2) THEN hom(nzb,1,5,:) = 0.0 hom(nzb,1,6,:) = 0.0 ENDIF hom(:,1,7,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 ) hom(:,1,23,:) = SPREAD( km(:,nys,nxl), 2, statistic_regions+1 ) hom(:,1,24,:) = SPREAD( kh(:,nys,nxl), 2, statistic_regions+1 ) IF ( ocean ) THEN ! !-- Store initial salinity profile hom(:,1,26,:) = SPREAD( sa(:,nys,nxl), 2, statistic_regions+1 ) ENDIF IF ( humidity ) THEN ! !-- Store initial profile of total water content, virtual potential !-- temperature hom(:,1,26,:) = SPREAD( q(:,nys,nxl), 2, statistic_regions+1 ) hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 ) IF ( cloud_physics .OR. cloud_droplets ) THEN ! !-- Store initial profile of specific humidity and potential !-- temperature hom(:,1,27,:) = SPREAD( q(:,nys,nxl), 2, statistic_regions+1 ) hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 ) ENDIF ENDIF IF ( passive_scalar ) THEN ! !-- Store initial scalar profile hom(:,1,26,:) = SPREAD( q(:,nys,nxl), 2, statistic_regions+1 ) ENDIF ! !-- Initialize fluxes at bottom surface IF ( use_surface_fluxes ) THEN IF ( constant_heatflux ) THEN ! !-- Heat flux is prescribed IF ( random_heatflux ) THEN CALL disturb_heatflux ELSE shf = surface_heatflux ! !-- Over topography surface_heatflux is replaced by wall_heatflux(0) IF ( TRIM( topography ) /= 'flat' ) THEN DO i = nxlg, nxrg DO j = nysg, nyng IF ( nzb_s_inner(j,i) /= 0 ) THEN shf(j,i) = wall_heatflux(0) ENDIF ENDDO ENDDO ENDIF ENDIF IF ( ASSOCIATED( shf_m ) ) shf_m = shf ENDIF ! !-- Determine the near-surface water flux IF ( humidity .OR. passive_scalar ) THEN IF ( constant_waterflux ) THEN qsws = surface_waterflux ! !-- Over topography surface_waterflux is replaced by !-- wall_humidityflux(0) IF ( TRIM( topography ) /= 'flat' ) THEN wall_qflux = wall_humidityflux DO i = nxlg, nxrg DO j = nysg, nyng IF ( nzb_s_inner(j,i) /= 0 ) THEN qsws(j,i) = wall_qflux(0) ENDIF ENDDO ENDDO ENDIF IF ( ASSOCIATED( qsws_m ) ) qsws_m = qsws ENDIF ENDIF ENDIF ! !-- Initialize fluxes at top surface !-- Currently, only the heatflux and salinity flux can be prescribed. !-- The latent flux is zero in this case! IF ( use_top_fluxes ) THEN IF ( constant_top_heatflux ) THEN ! !-- Heat flux is prescribed tswst = top_heatflux IF ( ASSOCIATED( tswst_m ) ) tswst_m = tswst IF ( humidity .OR. passive_scalar ) THEN qswst = 0.0 IF ( ASSOCIATED( qswst_m ) ) qswst_m = qswst ENDIF IF ( ocean ) THEN saswsb = bottom_salinityflux saswst = top_salinityflux ENDIF ENDIF ! !-- Initialization in case of a coupled model run IF ( coupling_mode == 'ocean_to_atmosphere' ) THEN tswst = 0.0 IF ( ASSOCIATED( tswst_m ) ) tswst_m = tswst ENDIF ENDIF ! !-- Initialize Prandtl layer quantities IF ( prandtl_layer ) THEN z0 = roughness_length IF ( .NOT. constant_heatflux ) THEN ! !-- Surface temperature is prescribed. Here the heat flux cannot be !-- simply estimated, because therefore rif, u* and theta* would have !-- to be computed by iteration. This is why the heat flux is assumed !-- to be zero before the first time step. It approaches its correct !-- value in the course of the first few time steps. shf = 0.0 IF ( ASSOCIATED( shf_m ) ) shf_m = 0.0 ENDIF IF ( humidity .OR. passive_scalar ) THEN IF ( .NOT. constant_waterflux ) THEN qsws = 0.0 IF ( ASSOCIATED( qsws_m ) ) qsws_m = 0.0 ENDIF ENDIF ENDIF ! !-- For the moment, vertical velocity is zero w = 0.0 ! !-- Initialize array sums (must be defined in first call of pres) sums = 0.0 ! !-- In case of iterative solvers, p must get an initial value IF ( psolver == 'multigrid' .OR. psolver == 'sor' ) p = 0.0 ! !-- Treating cloud physics, liquid water content and precipitation amount !-- are zero at beginning of the simulation IF ( cloud_physics ) THEN ql = 0.0 IF ( precipitation ) precipitation_amount = 0.0 ENDIF ! !-- Initialize quantities for special advections schemes CALL init_advec ! !-- Impose vortex with vertical axis on the initial velocity profile IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 ) THEN CALL init_rankine ENDIF ! !-- Impose temperature anomaly (advection test only) IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0 ) THEN CALL init_pt_anomaly ENDIF ! !-- If required, change the surface temperature at the start of the 3D run IF ( pt_surface_initial_change /= 0.0 ) THEN pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change ENDIF ! !-- If required, change the surface humidity/scalar at the start of the 3D !-- run IF ( ( humidity .OR. passive_scalar ) .AND. & q_surface_initial_change /= 0.0 ) THEN q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change ENDIF ! !-- Initialize the random number generator (from numerical recipes) CALL random_function_ini ! !-- Initialize old and new time levels. IF ( timestep_scheme(1:5) /= 'runge' ) THEN e_m = e; pt_m = pt; u_m = u; v_m = v; w_m = w; kh_m = kh; km_m = km ELSE te_m = 0.0; tpt_m = 0.0; tu_m = 0.0; tv_m = 0.0; tw_m = 0.0 ENDIF e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w IF ( humidity .OR. passive_scalar ) THEN IF ( ASSOCIATED( q_m ) ) q_m = q IF ( timestep_scheme(1:5) == 'runge' ) tq_m = 0.0 q_p = q IF ( humidity .AND. ASSOCIATED( vpt_m ) ) vpt_m = vpt ENDIF IF ( ocean ) THEN tsa_m = 0.0 sa_p = sa ENDIF ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data' .OR. & TRIM( initializing_actions ) == 'cyclic_fill' ) & THEN ! !-- When reading data for cyclic fill of 3D prerun data, first read !-- some of the global variables from restart file IF ( TRIM( initializing_actions ) == 'cyclic_fill' ) THEN DO i = 0, io_blocks-1 IF ( i == io_group ) THEN CALL read_parts_of_var_list CALL close_file( 13 ) ENDIF #if defined( __parallel ) CALL MPI_BARRIER( comm2d, ierr ) #endif ENDDO ! !-- Initialization of the turbulence recycling method IF ( turbulent_inflow ) THEN ! !-- Store temporally and horizontally averaged vertical profiles to be !-- used as mean inflow profiles ALLOCATE( mean_inflow_profiles(nzb:nzt+1,5) ) mean_inflow_profiles(:,1) = hom_sum(:,1,0) ! u mean_inflow_profiles(:,2) = hom_sum(:,2,0) ! v mean_inflow_profiles(:,4) = hom_sum(:,4,0) ! pt mean_inflow_profiles(:,5) = hom_sum(:,8,0) ! e ! !-- Use these mean profiles for the inflow (provided that Dirichlet !-- conditions are used) IF ( inflow_l ) THEN DO j = nysg, nyng DO k = nzb, nzt+1 u(k,j,nxlg:-1) = mean_inflow_profiles(k,1) v(k,j,nxlg:-1) = mean_inflow_profiles(k,2) w(k,j,nxlg:-1) = 0.0 pt(k,j,nxlg:-1) = mean_inflow_profiles(k,4) e(k,j,nxlg:-1) = mean_inflow_profiles(k,5) ENDDO ENDDO ENDIF ! !-- Calculate the damping factors to be used at the inflow. For a !-- turbulent inflow the turbulent fluctuations have to be limited !-- vertically because otherwise the turbulent inflow layer will grow !-- in time. IF ( inflow_damping_height == 9999999.9 ) THEN ! !-- Default: use the inversion height calculated by the prerun; if !-- this is zero, inflow_damping_height must be explicitly !-- specified. IF ( hom_sum(nzb+6,pr_palm,0) /= 0.0 ) THEN inflow_damping_height = hom_sum(nzb+6,pr_palm,0) ELSE WRITE( message_string, * ) 'inflow_damping_height must be ',& 'explicitly specified because&the inversion height ', & 'calculated by the prerun is zero.' CALL message( 'init_3d_model', 'PA0318', 1, 2, 0, 6, 0 ) ENDIF ENDIF IF ( inflow_damping_width == 9999999.9 ) THEN ! !-- Default for the transition range: one tenth of the undamped !-- layer inflow_damping_width = 0.1 * inflow_damping_height ENDIF ALLOCATE( inflow_damping_factor(nzb:nzt+1) ) DO k = nzb, nzt+1 IF ( zu(k) <= inflow_damping_height ) THEN inflow_damping_factor(k) = 1.0 ELSEIF ( zu(k) <= inflow_damping_height + & inflow_damping_width ) THEN inflow_damping_factor(k) = 1.0 - & ( zu(k) - inflow_damping_height ) / & inflow_damping_width ELSE inflow_damping_factor(k) = 0.0 ENDIF ENDDO ENDIF ENDIF ! !-- Read binary data from restart file DO i = 0, io_blocks-1 IF ( i == io_group ) THEN CALL read_3d_binary ENDIF #if defined( __parallel ) CALL MPI_BARRIER( comm2d, ierr ) #endif ENDDO ! !-- Inside buildings set velocities and TKE back to zero IF ( TRIM( initializing_actions ) == 'cyclic_fill' .AND. & topography /= 'flat' ) THEN ! !-- Inside buildings set velocities and TKE back to zero. !-- Other scalars (pt, q, s, km, kh, p, sa, ...) are ignored at present, !-- maybe revise later. IF ( timestep_scheme(1:5) == 'runge' ) THEN DO i = nxlg, nxrg DO j = nysg, nyng u (nzb:nzb_u_inner(j,i),j,i) = 0.0 v (nzb:nzb_v_inner(j,i),j,i) = 0.0 w (nzb:nzb_w_inner(j,i),j,i) = 0.0 e (nzb:nzb_w_inner(j,i),j,i) = 0.0 u_m(nzb:nzb_u_inner(j,i),j,i) = 0.0 v_m(nzb:nzb_v_inner(j,i),j,i) = 0.0 w_m(nzb:nzb_w_inner(j,i),j,i) = 0.0 e_m(nzb:nzb_w_inner(j,i),j,i) = 0.0 tu_m(nzb:nzb_u_inner(j,i),j,i) = 0.0 tv_m(nzb:nzb_v_inner(j,i),j,i) = 0.0 tw_m(nzb:nzb_w_inner(j,i),j,i) = 0.0 te_m(nzb:nzb_w_inner(j,i),j,i) = 0.0 tpt_m(nzb:nzb_w_inner(j,i),j,i) = 0.0 ENDDO ENDDO ELSE DO i = nxlg, nxrg DO j = nysg, nyng u (nzb:nzb_u_inner(j,i),j,i) = 0.0 v (nzb:nzb_v_inner(j,i),j,i) = 0.0 w (nzb:nzb_w_inner(j,i),j,i) = 0.0 e (nzb:nzb_w_inner(j,i),j,i) = 0.0 u_m(nzb:nzb_u_inner(j,i),j,i) = 0.0 v_m(nzb:nzb_v_inner(j,i),j,i) = 0.0 w_m(nzb:nzb_w_inner(j,i),j,i) = 0.0 e_m(nzb:nzb_w_inner(j,i),j,i) = 0.0 u_p(nzb:nzb_u_inner(j,i),j,i) = 0.0 v_p(nzb:nzb_v_inner(j,i),j,i) = 0.0 w_p(nzb:nzb_w_inner(j,i),j,i) = 0.0 e_p(nzb:nzb_w_inner(j,i),j,i) = 0.0 ENDDO ENDDO ENDIF ENDIF ! !-- Calculate initial temperature field and other constants used in case !-- of a sloping surface IF ( sloping_surface ) CALL init_slope ! !-- Initialize new time levels (only done in order to set boundary values !-- including ghost points) e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w IF ( humidity .OR. passive_scalar ) q_p = q IF ( ocean ) sa_p = sa ! !-- Allthough tendency arrays are set in prognostic_equations, they have !-- have to be predefined here because they are used (but multiplied with 0) !-- there before they are set. IF ( timestep_scheme(1:5) == 'runge' ) THEN te_m = 0.0; tpt_m = 0.0; tu_m = 0.0; tv_m = 0.0; tw_m = 0.0 IF ( humidity .OR. passive_scalar ) tq_m = 0.0 IF ( ocean ) tsa_m = 0.0 ENDIF ELSE ! !-- Actually this part of the programm should not be reached message_string = 'unknown initializing problem' CALL message( 'init_3d_model', 'PA0193', 1, 2, 0, 6, 0 ) ENDIF IF ( TRIM( initializing_actions ) /= 'read_restart_data' ) THEN ! !-- Initialize old timelevels needed for radiation boundary conditions IF ( outflow_l ) THEN u_m_l(:,:,:) = u(:,:,1:2) v_m_l(:,:,:) = v(:,:,0:1) w_m_l(:,:,:) = w(:,:,0:1) ENDIF IF ( outflow_r ) THEN u_m_r(:,:,:) = u(:,:,nx-1:nx) v_m_r(:,:,:) = v(:,:,nx-1:nx) w_m_r(:,:,:) = w(:,:,nx-1:nx) ENDIF IF ( outflow_s ) THEN u_m_s(:,:,:) = u(:,0:1,:) v_m_s(:,:,:) = v(:,1:2,:) w_m_s(:,:,:) = w(:,0:1,:) ENDIF IF ( outflow_n ) THEN u_m_n(:,:,:) = u(:,ny-1:ny,:) v_m_n(:,:,:) = v(:,ny-1:ny,:) w_m_n(:,:,:) = w(:,ny-1:ny,:) ENDIF ENDIF ! !-- Calculate the initial volume flow at the right and north boundary IF ( conserve_volume_flow ) THEN IF ( TRIM( initializing_actions ) == 'cyclic_fill' ) THEN volume_flow_initial_l = 0.0 volume_flow_area_l = 0.0 IF ( nxr == nx ) THEN DO j = nys, nyn DO k = nzb_2d(j,nx)+1, nzt volume_flow_initial_l(1) = volume_flow_initial_l(1) + & hom_sum(k,1,0) * dzw(k) volume_flow_area_l(1) = volume_flow_area_l(1) + dzw(k) ENDDO ENDDO ENDIF IF ( nyn == ny ) THEN DO i = nxl, nxr DO k = nzb_2d(ny,i)+1, nzt volume_flow_initial_l(2) = volume_flow_initial_l(2) + & hom_sum(k,2,0) * dzw(k) volume_flow_area_l(2) = volume_flow_area_l(2) + dzw(k) ENDDO ENDDO ENDIF #if defined( __parallel ) CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),& 2, MPI_REAL, MPI_SUM, comm2d, ierr ) CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1), & 2, MPI_REAL, MPI_SUM, comm2d, ierr ) #else volume_flow_initial = volume_flow_initial_l volume_flow_area = volume_flow_area_l #endif ELSEIF ( TRIM( initializing_actions ) /= 'read_restart_data' ) THEN volume_flow_initial_l = 0.0 volume_flow_area_l = 0.0 IF ( nxr == nx ) THEN DO j = nys, nyn DO k = nzb_2d(j,nx)+1, nzt volume_flow_initial_l(1) = volume_flow_initial_l(1) + & u(k,j,nx) * dzw(k) volume_flow_area_l(1) = volume_flow_area_l(1) + dzw(k) ENDDO ENDDO ENDIF IF ( nyn == ny ) THEN DO i = nxl, nxr DO k = nzb_2d(ny,i)+1, nzt volume_flow_initial_l(2) = volume_flow_initial_l(2) + & v(k,ny,i) * dzw(k) volume_flow_area_l(2) = volume_flow_area_l(2) + dzw(k) ENDDO ENDDO ENDIF #if defined( __parallel ) CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),& 2, MPI_REAL, MPI_SUM, comm2d, ierr ) CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1), & 2, MPI_REAL, MPI_SUM, comm2d, ierr ) #else volume_flow_initial = volume_flow_initial_l volume_flow_area = volume_flow_area_l #endif ENDIF ! !-- In case of 'bulk_velocity' mode, volume_flow_initial is calculated !-- from u|v_bulk instead IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' ) THEN volume_flow_initial(1) = u_bulk * volume_flow_area(1) volume_flow_initial(2) = v_bulk * volume_flow_area(2) ENDIF ENDIF ! !-- Impose random perturbation on the horizontal velocity field and then !-- remove the divergences from the velocity field at the initial stage IF ( create_disturbances .AND. & TRIM( initializing_actions ) /= 'read_restart_data' .AND. & TRIM( initializing_actions ) /= 'cyclic_fill' ) THEN CALL disturb_field( nzb_u_inner, tend, u ) CALL disturb_field( nzb_v_inner, tend, v ) n_sor = nsor_ini CALL pres n_sor = nsor ENDIF ! !-- Initialization of the leaf area density IF ( plant_canopy ) THEN SELECT CASE ( TRIM( canopy_mode ) ) CASE( 'block' ) DO i = nxlg, nxrg DO j = nysg, nyng lad_s(:,j,i) = lad(:) cdc(:,j,i) = drag_coefficient IF ( passive_scalar ) THEN sls(:,j,i) = leaf_surface_concentration sec(:,j,i) = scalar_exchange_coefficient ENDIF ENDDO ENDDO CASE DEFAULT ! !-- The DEFAULT case is reached either if the parameter !-- canopy mode contains a wrong character string or if the !-- user has coded a special case in the user interface. !-- There, the subroutine user_init_plant_canopy checks !-- which of these two conditions applies. CALL user_init_plant_canopy END SELECT CALL exchange_horiz( lad_s, nbgp ) CALL exchange_horiz( cdc, nbgp ) IF ( passive_scalar ) THEN CALL exchange_horiz( sls, nbgp ) CALL exchange_horiz( sec, nbgp ) ENDIF ! !-- Sharp boundaries of the plant canopy in horizontal directions !-- In vertical direction the interpolation is retained, as the leaf !-- area density is initialised by prescribing a vertical profile !-- consisting of piecewise linear segments. The upper boundary !-- of the plant canopy is now defined by lad_w(pch_index,:,:) = 0.0. DO i = nxl, nxr DO j = nys, nyn DO k = nzb, nzt+1 IF ( lad_s(k,j,i) > 0.0 ) THEN lad_u(k,j,i) = lad_s(k,j,i) lad_u(k,j,i+1) = lad_s(k,j,i) lad_v(k,j,i) = lad_s(k,j,i) lad_v(k,j+1,i) = lad_s(k,j,i) ENDIF ENDDO DO k = nzb, nzt lad_w(k,j,i) = 0.5 * ( lad_s(k+1,j,i) + lad_s(k,j,i) ) ENDDO ENDDO ENDDO lad_w(pch_index,:,:) = 0.0 lad_w(nzt+1,:,:) = lad_w(nzt,:,:) CALL exchange_horiz( lad_u, nbgp ) CALL exchange_horiz( lad_v, nbgp ) CALL exchange_horiz( lad_w, nbgp ) ! !-- Initialisation of the canopy heat source distribution IF ( cthf /= 0.0 ) THEN ! !-- Piecewise evaluation of the leaf area index by !-- integration of the leaf area density lai(:,:,:) = 0.0 DO i = nxlg, nxrg DO j = nysg, nyng DO k = pch_index-1, 0, -1 lai(k,j,i) = lai(k+1,j,i) + & ( 0.5 * ( lad_w(k+1,j,i) + & lad_s(k+1,j,i) ) * & ( zw(k+1) - zu(k+1) ) ) + & ( 0.5 * ( lad_w(k,j,i) + & lad_s(k+1,j,i) ) * & ( zu(k+1) - zw(k) ) ) ENDDO ENDDO ENDDO ! !-- Evaluation of the upward kinematic vertical heat flux within the !-- canopy DO i = nxlg, nxrg DO j = nysg, nyng DO k = 0, pch_index canopy_heat_flux(k,j,i) = cthf * & exp( -0.6 * lai(k,j,i) ) ENDDO ENDDO ENDDO ! !-- The near surface heat flux is derived from the heat flux !-- distribution within the canopy shf(:,:) = canopy_heat_flux(0,:,:) IF ( ASSOCIATED( shf_m ) ) shf_m = shf ENDIF ENDIF ! !-- If required, initialize dvrp-software IF ( dt_dvrp /= 9999999.9 ) CALL init_dvrp IF ( ocean ) THEN ! !-- Initialize quantities needed for the ocean model CALL init_ocean ELSE ! !-- Initialize quantities for handling cloud physics !-- This routine must be called before init_particles, because !-- otherwise, array pt_d_t, needed in data_output_dvrp (called by !-- init_particles) is not defined. CALL init_cloud_physics ENDIF ! !-- If required, initialize particles IF ( particle_advection ) CALL init_particles ! !-- Initialize the ws-scheme. IF ( ws_scheme_sca .OR. ws_scheme_mom ) CALL ws_init ! !-- Setting weighting factors for calculation of perturbation pressure !-- and turbulent quantities from the RK substeps IF ( TRIM(timestep_scheme) == 'runge-kutta-3' ) THEN ! for RK3-method weight_substep(1) = 1./6. weight_substep(2) = 3./10. weight_substep(3) = 8./15. weight_pres(1) = 1./3. weight_pres(2) = 5./12. weight_pres(3) = 1./4. ELSEIF ( TRIM(timestep_scheme) == 'runge-kutta-2' ) THEN ! for RK2-method weight_substep(1) = 1./2. weight_substep(2) = 1./2. weight_pres(1) = 1./2. weight_pres(2) = 1./2. ELSE ! for Euler- and leapfrog-method weight_substep(1) = 1.0 weight_pres(1) = 1.0 ENDIF ! !-- Initialize Rayleigh damping factors rdf = 0.0 IF ( rayleigh_damping_factor /= 0.0 ) THEN IF ( .NOT. ocean ) THEN DO k = nzb+1, nzt IF ( zu(k) >= rayleigh_damping_height ) THEN rdf(k) = rayleigh_damping_factor * & ( SIN( pi * 0.5 * ( zu(k) - rayleigh_damping_height ) & / ( zu(nzt) - rayleigh_damping_height ) )& )**2 ENDIF ENDDO ELSE DO k = nzt, nzb+1, -1 IF ( zu(k) <= rayleigh_damping_height ) THEN rdf(k) = rayleigh_damping_factor * & ( SIN( pi * 0.5 * ( rayleigh_damping_height - zu(k) ) & / ( rayleigh_damping_height - zu(nzb+1)))& )**2 ENDIF ENDDO ENDIF ENDIF ! !-- Initialize the starting level and the vertical smoothing factor used for !-- the external pressure gradient dp_smooth_factor = 1.0 IF ( dp_external ) THEN ! !-- Set the starting level dp_level_ind_b only if it has not been set before !-- (e.g. in init_grid). IF ( dp_level_ind_b == 0 ) THEN ind_array = MINLOC( ABS( dp_level_b - zu ) ) dp_level_ind_b = ind_array(1) - 1 + nzb ! MINLOC uses lower array bound 1 ENDIF IF ( dp_smooth ) THEN dp_smooth_factor(:dp_level_ind_b) = 0.0 DO k = dp_level_ind_b+1, nzt dp_smooth_factor(k) = 0.5 * ( 1.0 + SIN( pi * & ( REAL( k - dp_level_ind_b ) / & REAL( nzt - dp_level_ind_b ) - 0.5 ) ) ) ENDDO ENDIF ENDIF ! !-- Initialize diffusivities used within the outflow damping layer in case of !-- non-cyclic lateral boundaries. A linear increase is assumed over the first !-- half of the width of the damping layer IF ( bc_lr_dirrad ) THEN DO i = nxlg, nxrg IF ( i >= nx - outflow_damping_width ) THEN km_damp_x(i) = km_damp_max * MIN( 1.0, & ( i - ( nx - outflow_damping_width ) ) / & REAL( outflow_damping_width/2 ) & ) ELSE km_damp_x(i) = 0.0 ENDIF ENDDO ELSEIF ( bc_lr_raddir ) THEN DO i = nxlg, nxrg IF ( i <= outflow_damping_width ) THEN km_damp_x(i) = km_damp_max * MIN( 1.0, & ( outflow_damping_width - i ) / & REAL( outflow_damping_width/2 ) & ) ELSE km_damp_x(i) = 0.0 ENDIF ENDDO ENDIF IF ( bc_ns_dirrad ) THEN DO j = nysg, nyng IF ( j >= ny - outflow_damping_width ) THEN km_damp_y(j) = km_damp_max * MIN( 1.0, & ( j - ( ny - outflow_damping_width ) ) / & REAL( outflow_damping_width/2 ) & ) ELSE km_damp_y(j) = 0.0 ENDIF ENDDO ELSEIF ( bc_ns_raddir ) THEN DO j = nysg, nyng IF ( j <= outflow_damping_width ) THEN km_damp_y(j) = km_damp_max * MIN( 1.0, & ( outflow_damping_width - j ) / & REAL( outflow_damping_width/2 ) & ) ELSE km_damp_y(j) = 0.0 ENDIF ENDDO ENDIF ! !-- Initialize local summation arrays for routine flow_statistics. !-- This is necessary because they may not yet have been initialized when they !-- are called from flow_statistics (or - depending on the chosen model run - !-- are never initialized) sums_divnew_l = 0.0 sums_divold_l = 0.0 sums_l_l = 0.0 sums_up_fraction_l = 0.0 sums_wsts_bc_l = 0.0 ! !-- Pre-set masks for regional statistics. Default is the total model domain. rmask = 1.0 ! !-- User-defined initializing actions. Check afterwards, if maximum number !-- of allowed timeseries is exceeded CALL user_init IF ( dots_num > dots_max ) THEN WRITE( message_string, * ) 'number of time series quantities exceeds', & ' its maximum of dots_max = ', dots_max, & ' &Please increase dots_max in modules.f90.' CALL message( 'init_3d_model', 'PA0194', 1, 2, 0, 6, 0 ) ENDIF ! !-- Input binary data file is not needed anymore. This line must be placed !-- after call of user_init! CALL close_file( 13 ) ! !-- Compute total sum of active mask grid points !-- ngp_2dh: number of grid points of a horizontal cross section through the !-- total domain !-- ngp_3d: number of grid points of the total domain ngp_2dh_outer_l = 0 ngp_2dh_outer = 0 ngp_2dh_s_inner_l = 0 ngp_2dh_s_inner = 0 ngp_2dh_l = 0 ngp_2dh = 0 ngp_3d_inner_l = 0.0 ngp_3d_inner = 0 ngp_3d = 0 ngp_sums = ( nz + 2 ) * ( pr_palm + max_pr_user ) DO sr = 0, statistic_regions DO i = nxl, nxr DO j = nys, nyn IF ( rmask(j,i,sr) == 1.0 ) THEN ! !-- All xy-grid points ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1 ! !-- xy-grid points above topography DO k = nzb_s_outer(j,i), nz + 1 ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr) + 1 ENDDO DO k = nzb_s_inner(j,i), nz + 1 ngp_2dh_s_inner_l(k,sr) = ngp_2dh_s_inner_l(k,sr) + 1 ENDDO ! !-- All grid points of the total domain above topography ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr) + & ( nz - nzb_s_inner(j,i) + 2 ) ENDIF ENDDO ENDDO ENDDO sr = statistic_regions + 1 #if defined( __parallel ) IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM, & comm2d, ierr ) IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr, & MPI_INTEGER, MPI_SUM, comm2d, ierr ) IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( ngp_2dh_s_inner_l(0,0), ngp_2dh_s_inner(0,0), & (nz+2)*sr, MPI_INTEGER, MPI_SUM, comm2d, ierr ) IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner_tmp(0), sr, MPI_REAL, & MPI_SUM, comm2d, ierr ) ngp_3d_inner = INT( ngp_3d_inner_tmp, KIND = SELECTED_INT_KIND( 18 ) ) #else ngp_2dh = ngp_2dh_l ngp_2dh_outer = ngp_2dh_outer_l ngp_2dh_s_inner = ngp_2dh_s_inner_l ngp_3d_inner = INT( ngp_3d_inner_l, KIND = SELECTED_INT_KIND( 18 ) ) #endif ngp_3d = INT ( ngp_2dh, KIND = SELECTED_INT_KIND( 18 ) ) * & INT ( (nz + 2 ), KIND = SELECTED_INT_KIND( 18 ) ) ! !-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics, !-- buoyancy, etc. A zero value will occur for cases where all grid points of !-- the respective subdomain lie below the surface topography ngp_2dh_outer = MAX( 1, ngp_2dh_outer(:,:) ) ngp_3d_inner = MAX( INT(1, KIND = SELECTED_INT_KIND( 18 )), & ngp_3d_inner(:) ) ngp_2dh_s_inner = MAX( 1, ngp_2dh_s_inner(:,:) ) DEALLOCATE( ngp_2dh_l, ngp_2dh_outer_l, ngp_3d_inner_l, ngp_3d_inner_tmp ) END SUBROUTINE init_3d_model