[1] | 1 | SUBROUTINE init_3d_model |
---|
| 2 | |
---|
[1036] | 3 | !--------------------------------------------------------------------------------! |
---|
| 4 | ! This file is part of PALM. |
---|
| 5 | ! |
---|
| 6 | ! PALM is free software: you can redistribute it and/or modify it under the terms |
---|
| 7 | ! of the GNU General Public License as published by the Free Software Foundation, |
---|
| 8 | ! either version 3 of the License, or (at your option) any later version. |
---|
| 9 | ! |
---|
| 10 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
| 11 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
| 12 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
| 13 | ! |
---|
| 14 | ! You should have received a copy of the GNU General Public License along with |
---|
| 15 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
| 16 | ! |
---|
[1310] | 17 | ! Copyright 1997-2014 Leibniz Universitaet Hannover |
---|
[1036] | 18 | !--------------------------------------------------------------------------------! |
---|
| 19 | ! |
---|
[254] | 20 | ! Current revisions: |
---|
[732] | 21 | ! ------------------ |
---|
[1320] | 22 | ! ONLY-attribute added to USE-statements, |
---|
| 23 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
| 24 | ! kinds are defined in new module kinds, |
---|
| 25 | ! old module precision_kind is removed, |
---|
| 26 | ! revision history before 2012 removed, |
---|
| 27 | ! comment fields (!:) to be used for variable explanations added to |
---|
| 28 | ! all variable declaration statements |
---|
[1054] | 29 | ! |
---|
[1320] | 30 | ! module interfaces removed |
---|
[1319] | 31 | ! |
---|
[1054] | 32 | ! Former revisions: |
---|
| 33 | ! ----------------- |
---|
| 34 | ! $Id: init_3d_model.f90 1320 2014-03-20 08:40:49Z raasch $ |
---|
| 35 | ! |
---|
[1317] | 36 | ! 1316 2014-03-17 07:44:59Z heinze |
---|
| 37 | ! Bugfix: allocation of w_subs |
---|
| 38 | ! |
---|
[1300] | 39 | ! 1299 2014-03-06 13:15:21Z heinze |
---|
| 40 | ! Allocate w_subs due to extension of large scale subsidence in combination |
---|
| 41 | ! with large scale forcing data (LSF_DATA) |
---|
| 42 | ! |
---|
[1242] | 43 | ! 1241 2013-10-30 11:36:58Z heinze |
---|
| 44 | ! Overwrite initial profiles in case of nudging |
---|
| 45 | ! Inititialize shf and qsws in case of large_scale_forcing |
---|
| 46 | ! |
---|
[1222] | 47 | ! 1221 2013-09-10 08:59:13Z raasch |
---|
| 48 | ! +rflags_s_inner in copyin statement, use copyin for most arrays instead of |
---|
| 49 | ! copy |
---|
| 50 | ! |
---|
[1213] | 51 | ! 1212 2013-08-15 08:46:27Z raasch |
---|
| 52 | ! array tri is allocated and included in data copy statement |
---|
| 53 | ! |
---|
[1196] | 54 | ! 1195 2013-07-01 12:27:57Z heinze |
---|
| 55 | ! Bugfix: move allocation of ref_state to parin.f90 and read_var_list.f90 |
---|
| 56 | ! |
---|
[1182] | 57 | ! 1179 2013-06-14 05:57:58Z raasch |
---|
| 58 | ! allocate and set ref_state to be used in buoyancy terms |
---|
| 59 | ! |
---|
[1172] | 60 | ! 1171 2013-05-30 11:27:45Z raasch |
---|
| 61 | ! diss array is allocated with full size if accelerator boards are used |
---|
| 62 | ! |
---|
[1160] | 63 | ! 1159 2013-05-21 11:58:22Z fricke |
---|
| 64 | ! -bc_lr_dirneu, bc_lr_neudir, bc_ns_dirneu, bc_ns_neudir |
---|
| 65 | ! |
---|
[1154] | 66 | ! 1153 2013-05-10 14:33:08Z raasch |
---|
| 67 | ! diss array is allocated with dummy elements even if it is not needed |
---|
[1171] | 68 | ! (required by PGI 13.4 / CUDA 5.0) |
---|
[1154] | 69 | ! |
---|
[1116] | 70 | ! 1115 2013-03-26 18:16:16Z hoffmann |
---|
| 71 | ! unused variables removed |
---|
| 72 | ! |
---|
[1114] | 73 | ! 1113 2013-03-10 02:48:14Z raasch |
---|
| 74 | ! openACC directive modified |
---|
| 75 | ! |
---|
[1112] | 76 | ! 1111 2013-03-08 23:54:10Z raasch |
---|
| 77 | ! openACC directives added for pres |
---|
| 78 | ! array diss allocated only if required |
---|
| 79 | ! |
---|
[1093] | 80 | ! 1092 2013-02-02 11:24:22Z raasch |
---|
| 81 | ! unused variables removed |
---|
| 82 | ! |
---|
[1066] | 83 | ! 1065 2012-11-22 17:42:36Z hoffmann |
---|
| 84 | ! allocation of diss (dissipation rate) in case of turbulence = .TRUE. added |
---|
| 85 | ! |
---|
[1054] | 86 | ! 1053 2012-11-13 17:11:03Z hoffmann |
---|
[1053] | 87 | ! allocation and initialisation of necessary data arrays for the two-moment |
---|
| 88 | ! cloud physics scheme the two new prognostic equations (nr, qr): |
---|
| 89 | ! +dr, lambda_r, mu_r, sed_*, xr, *s, *sws, *swst, *, *_p, t*_m, *_1, *_2, *_3, |
---|
| 90 | ! +tend_*, prr |
---|
[979] | 91 | ! |
---|
[1037] | 92 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
| 93 | ! code put under GPL (PALM 3.9) |
---|
| 94 | ! |
---|
[1033] | 95 | ! 1032 2012-10-21 13:03:21Z letzel |
---|
| 96 | ! save memory by not allocating pt_2 in case of neutral = .T. |
---|
| 97 | ! |
---|
[1026] | 98 | ! 1025 2012-10-07 16:04:41Z letzel |
---|
| 99 | ! bugfix: swap indices of mask for ghost boundaries |
---|
| 100 | ! |
---|
[1017] | 101 | ! 1015 2012-09-27 09:23:24Z raasch |
---|
| 102 | ! mask is set to zero for ghost boundaries |
---|
| 103 | ! |
---|
[1011] | 104 | ! 1010 2012-09-20 07:59:54Z raasch |
---|
| 105 | ! cpp switch __nopointer added for pointer free version |
---|
| 106 | ! |
---|
[1004] | 107 | ! 1003 2012-09-14 14:35:53Z raasch |
---|
| 108 | ! nxra,nyna, nzta replaced ny nxr, nyn, nzt |
---|
| 109 | ! |
---|
[1002] | 110 | ! 1001 2012-09-13 14:08:46Z raasch |
---|
| 111 | ! all actions concerning leapfrog scheme removed |
---|
| 112 | ! |
---|
[997] | 113 | ! 996 2012-09-07 10:41:47Z raasch |
---|
| 114 | ! little reformatting |
---|
| 115 | ! |
---|
[979] | 116 | ! 978 2012-08-09 08:28:32Z fricke |
---|
[978] | 117 | ! outflow damping layer removed |
---|
| 118 | ! roughness length for scalar quantites z0h added |
---|
| 119 | ! damping zone for the potential temperatur in case of non-cyclic lateral |
---|
| 120 | ! boundaries added |
---|
| 121 | ! initialization of ptdf_x, ptdf_y |
---|
| 122 | ! initialization of c_u_m, c_u_m_l, c_v_m, c_v_m_l, c_w_m, c_w_m_l |
---|
[708] | 123 | ! |
---|
[850] | 124 | ! 849 2012-03-15 10:35:09Z raasch |
---|
| 125 | ! init_particles renamed lpm_init |
---|
| 126 | ! |
---|
[826] | 127 | ! 825 2012-02-19 03:03:44Z raasch |
---|
| 128 | ! wang_collision_kernel renamed wang_kernel |
---|
| 129 | ! |
---|
[1] | 130 | ! Revision 1.1 1998/03/09 16:22:22 raasch |
---|
| 131 | ! Initial revision |
---|
| 132 | ! |
---|
| 133 | ! |
---|
| 134 | ! Description: |
---|
| 135 | ! ------------ |
---|
| 136 | ! Allocation of arrays and initialization of the 3D model via |
---|
| 137 | ! a) pre-run the 1D model |
---|
| 138 | ! or |
---|
| 139 | ! b) pre-set constant linear profiles |
---|
| 140 | ! or |
---|
| 141 | ! c) read values of a previous run |
---|
| 142 | !------------------------------------------------------------------------------! |
---|
| 143 | |
---|
[667] | 144 | USE advec_ws |
---|
[1320] | 145 | |
---|
[1] | 146 | USE arrays_3d |
---|
[1320] | 147 | |
---|
| 148 | USE cloud_parameters, & |
---|
| 149 | ONLY: nc_const, precipitation_amount, precipitation_rate, prr |
---|
| 150 | |
---|
| 151 | USE constants, & |
---|
| 152 | ONLY: pi |
---|
| 153 | |
---|
[1] | 154 | USE control_parameters |
---|
[1320] | 155 | |
---|
| 156 | USE grid_variables, & |
---|
| 157 | ONLY: dx, dy |
---|
| 158 | |
---|
[1] | 159 | USE indices |
---|
[1320] | 160 | |
---|
| 161 | USE kinds |
---|
| 162 | |
---|
[1241] | 163 | USE ls_forcing_mod |
---|
[1320] | 164 | |
---|
| 165 | USE model_1d, & |
---|
| 166 | ONLY: e1d, kh1d, km1d, l1d, rif1d, u1d, us1d, usws1d, v1d, vsws1d |
---|
| 167 | |
---|
[51] | 168 | USE netcdf_control |
---|
[1320] | 169 | |
---|
| 170 | USE particle_attributes, & |
---|
| 171 | ONLY: particle_advection, use_sgs_for_particles, wang_kernel |
---|
| 172 | |
---|
[1] | 173 | USE pegrid |
---|
[1320] | 174 | |
---|
| 175 | USE random_function_mod |
---|
| 176 | |
---|
| 177 | USE statistics, & |
---|
| 178 | ONLY: hom, hom_sum, pr_palm, rmask, spectrum_x, spectrum_y, & |
---|
| 179 | statistic_regions, sums, sums_divnew_l, sums_divold_l, sums_l, & |
---|
| 180 | sums_l_l, sums_up_fraction_l, sums_wsts_bc_l, ts_value, & |
---|
| 181 | weight_pres, weight_substep |
---|
| 182 | |
---|
| 183 | USE transpose_indices |
---|
[1] | 184 | |
---|
| 185 | IMPLICIT NONE |
---|
| 186 | |
---|
[1320] | 187 | INTEGER(iwp) :: i !: |
---|
| 188 | INTEGER(iwp) :: ind_array(1) !: |
---|
| 189 | INTEGER(iwp) :: j !: |
---|
| 190 | INTEGER(iwp) :: k !: |
---|
| 191 | INTEGER(iwp) :: sr !: |
---|
[1] | 192 | |
---|
[1320] | 193 | INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: ngp_2dh_l !: |
---|
[1] | 194 | |
---|
[1320] | 195 | INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE :: ngp_2dh_outer_l !: |
---|
| 196 | INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE :: ngp_2dh_s_inner_l !: |
---|
[1] | 197 | |
---|
[1320] | 198 | REAL(wp), DIMENSION(1:2) :: volume_flow_area_l !: |
---|
| 199 | REAL(wp), DIMENSION(1:2) :: volume_flow_initial_l !: |
---|
[1] | 200 | |
---|
[1320] | 201 | REAL(wp), DIMENSION(:), ALLOCATABLE :: ngp_3d_inner_l !: |
---|
| 202 | REAL(wp), DIMENSION(:), ALLOCATABLE :: ngp_3d_inner_tmp !: |
---|
[1] | 203 | |
---|
[485] | 204 | |
---|
[1] | 205 | ! |
---|
| 206 | !-- Allocate arrays |
---|
| 207 | ALLOCATE( ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions), & |
---|
| 208 | ngp_3d(0:statistic_regions), & |
---|
| 209 | ngp_3d_inner(0:statistic_regions), & |
---|
| 210 | ngp_3d_inner_l(0:statistic_regions), & |
---|
[485] | 211 | ngp_3d_inner_tmp(0:statistic_regions), & |
---|
[1] | 212 | sums_divnew_l(0:statistic_regions), & |
---|
| 213 | sums_divold_l(0:statistic_regions) ) |
---|
[1195] | 214 | ALLOCATE( dp_smooth_factor(nzb:nzt), rdf(nzb+1:nzt), rdf_sc(nzb+1:nzt) ) |
---|
[143] | 215 | ALLOCATE( ngp_2dh_outer(nzb:nzt+1,0:statistic_regions), & |
---|
[1] | 216 | ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions), & |
---|
[132] | 217 | ngp_2dh_s_inner(nzb:nzt+1,0:statistic_regions), & |
---|
| 218 | ngp_2dh_s_inner_l(nzb:nzt+1,0:statistic_regions), & |
---|
[996] | 219 | rmask(nysg:nyng,nxlg:nxrg,0:statistic_regions), & |
---|
[87] | 220 | sums(nzb:nzt+1,pr_palm+max_pr_user), & |
---|
| 221 | sums_l(nzb:nzt+1,pr_palm+max_pr_user,0:threads_per_task-1), & |
---|
[1] | 222 | sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1), & |
---|
| 223 | sums_up_fraction_l(10,3,0:statistic_regions), & |
---|
[48] | 224 | sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions), & |
---|
[394] | 225 | ts_value(dots_max,0:statistic_regions) ) |
---|
[978] | 226 | ALLOCATE( ptdf_x(nxlg:nxrg), ptdf_y(nysg:nyng) ) |
---|
[1] | 227 | |
---|
[1001] | 228 | ALLOCATE( rif(nysg:nyng,nxlg:nxrg), shf(nysg:nyng,nxlg:nxrg), & |
---|
| 229 | ts(nysg:nyng,nxlg:nxrg), tswst(nysg:nyng,nxlg:nxrg), & |
---|
| 230 | us(nysg:nyng,nxlg:nxrg), usws(nysg:nyng,nxlg:nxrg), & |
---|
| 231 | uswst(nysg:nyng,nxlg:nxrg), vsws(nysg:nyng,nxlg:nxrg), & |
---|
| 232 | vswst(nysg:nyng,nxlg:nxrg), z0(nysg:nyng,nxlg:nxrg), & |
---|
[978] | 233 | z0h(nysg:nyng,nxlg:nxrg) ) |
---|
[1] | 234 | |
---|
[1010] | 235 | ALLOCATE( d(nzb+1:nzt,nys:nyn,nxl:nxr), & |
---|
| 236 | kh(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 237 | km(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 238 | p(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 239 | tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 240 | |
---|
| 241 | #if defined( __nopointer ) |
---|
| 242 | ALLOCATE( e(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 243 | e_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 244 | pt(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 245 | pt_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 246 | u(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 247 | u_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 248 | v(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 249 | v_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 250 | w(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 251 | w_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 252 | te_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 253 | tpt_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 254 | tu_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 255 | tv_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 256 | tw_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 257 | #else |
---|
| 258 | ALLOCATE( e_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 259 | e_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 260 | e_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 261 | pt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 262 | pt_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 263 | u_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 264 | u_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 265 | u_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 266 | v_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 267 | v_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 268 | v_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 269 | w_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 270 | w_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
[667] | 271 | w_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[1032] | 272 | IF ( .NOT. neutral ) THEN |
---|
| 273 | ALLOCATE( pt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 274 | ENDIF |
---|
[1010] | 275 | #endif |
---|
| 276 | |
---|
[673] | 277 | ! |
---|
[707] | 278 | !-- Following array is required for perturbation pressure within the iterative |
---|
| 279 | !-- pressure solvers. For the multistep schemes (Runge-Kutta), array p holds |
---|
| 280 | !-- the weighted average of the substeps and cannot be used in the Poisson |
---|
| 281 | !-- solver. |
---|
| 282 | IF ( psolver == 'sor' ) THEN |
---|
| 283 | ALLOCATE( p_loc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 284 | ELSEIF ( psolver == 'multigrid' ) THEN |
---|
| 285 | ! |
---|
| 286 | !-- For performance reasons, multigrid is using one ghost layer only |
---|
| 287 | ALLOCATE( p_loc(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) ) |
---|
[673] | 288 | ENDIF |
---|
[1] | 289 | |
---|
[1111] | 290 | ! |
---|
| 291 | !-- Array for storing constant coeffficients of the tridiagonal solver |
---|
| 292 | IF ( psolver == 'poisfft' ) THEN |
---|
[1212] | 293 | ALLOCATE( tri(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1,2) ) |
---|
[1111] | 294 | ALLOCATE( tric(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1) ) |
---|
| 295 | ENDIF |
---|
| 296 | |
---|
[75] | 297 | IF ( humidity .OR. passive_scalar ) THEN |
---|
[1] | 298 | ! |
---|
[75] | 299 | !-- 2D-humidity/scalar arrays |
---|
[1001] | 300 | ALLOCATE ( qs(nysg:nyng,nxlg:nxrg), & |
---|
| 301 | qsws(nysg:nyng,nxlg:nxrg), & |
---|
| 302 | qswst(nysg:nyng,nxlg:nxrg) ) |
---|
[1] | 303 | |
---|
| 304 | ! |
---|
[75] | 305 | !-- 3D-humidity/scalar arrays |
---|
[1010] | 306 | #if defined( __nopointer ) |
---|
| 307 | ALLOCATE( q(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 308 | q_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 309 | tq_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 310 | #else |
---|
[667] | 311 | ALLOCATE( q_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 312 | q_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 313 | q_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[1010] | 314 | #endif |
---|
[1] | 315 | |
---|
| 316 | ! |
---|
[75] | 317 | !-- 3D-arrays needed for humidity only |
---|
| 318 | IF ( humidity ) THEN |
---|
[1010] | 319 | #if defined( __nopointer ) |
---|
| 320 | ALLOCATE( vpt(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 321 | #else |
---|
[667] | 322 | ALLOCATE( vpt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[1010] | 323 | #endif |
---|
[1] | 324 | |
---|
| 325 | IF ( cloud_physics ) THEN |
---|
[1053] | 326 | |
---|
[1] | 327 | ! |
---|
| 328 | !-- Liquid water content |
---|
[1010] | 329 | #if defined( __nopointer ) |
---|
| 330 | ALLOCATE ( ql(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 331 | #else |
---|
[667] | 332 | ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[1010] | 333 | #endif |
---|
[72] | 334 | ! |
---|
| 335 | !-- Precipitation amount and rate (only needed if output is switched) |
---|
[667] | 336 | ALLOCATE( precipitation_amount(nysg:nyng,nxlg:nxrg), & |
---|
| 337 | precipitation_rate(nysg:nyng,nxlg:nxrg) ) |
---|
[1053] | 338 | |
---|
| 339 | IF ( icloud_scheme == 0 ) THEN |
---|
| 340 | ! |
---|
[1115] | 341 | !-- 1D-arrays |
---|
| 342 | ALLOCATE ( nc_1d(nzb:nzt+1), pt_1d(nzb:nzt+1), & |
---|
| 343 | q_1d(nzb:nzt+1), qc_1d(nzb:nzt+1) ) |
---|
[1053] | 344 | ! |
---|
[1115] | 345 | !-- 3D-cloud water content |
---|
[1053] | 346 | #if defined( __nopointer ) |
---|
[1115] | 347 | ALLOCATE( qc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[1053] | 348 | #else |
---|
[1115] | 349 | ALLOCATE( qc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[1053] | 350 | #endif |
---|
| 351 | ! |
---|
| 352 | !-- 3D-tendency arrays |
---|
[1115] | 353 | ALLOCATE( tend_pt(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 354 | tend_q(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 355 | |
---|
| 356 | IF ( precipitation ) THEN |
---|
[1053] | 357 | ! |
---|
[1115] | 358 | !-- 1D-arrays |
---|
| 359 | ALLOCATE ( nr_1d(nzb:nzt+1), qr_1d(nzb:nzt+1) ) |
---|
| 360 | ! |
---|
| 361 | ! |
---|
| 362 | !-- 3D-tendency arrays |
---|
| 363 | ALLOCATE( tend_nr(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 364 | tend_qr(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 365 | ! |
---|
| 366 | !-- 2D-rain water content and rain drop concentration arrays |
---|
| 367 | ALLOCATE ( qrs(nysg:nyng,nxlg:nxrg), & |
---|
| 368 | qrsws(nysg:nyng,nxlg:nxrg), & |
---|
| 369 | qrswst(nysg:nyng,nxlg:nxrg), & |
---|
| 370 | nrs(nysg:nyng,nxlg:nxrg), & |
---|
| 371 | nrsws(nysg:nyng,nxlg:nxrg), & |
---|
| 372 | nrswst(nysg:nyng,nxlg:nxrg) ) |
---|
| 373 | ! |
---|
| 374 | !-- 3D-rain water content, rain drop concentration arrays |
---|
| 375 | #if defined( __nopointer ) |
---|
| 376 | ALLOCATE( nr(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 377 | nr_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 378 | qr(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 379 | qr_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 380 | tnr_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 381 | tqr_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 382 | #else |
---|
| 383 | ALLOCATE( nr_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 384 | nr_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 385 | nr_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 386 | qr_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 387 | qr_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 388 | qr_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 389 | #endif |
---|
| 390 | ! |
---|
| 391 | !-- 3d-precipitation rate |
---|
[1053] | 392 | ALLOCATE( prr(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 393 | ENDIF |
---|
| 394 | |
---|
| 395 | ENDIF |
---|
[1] | 396 | ENDIF |
---|
| 397 | |
---|
| 398 | IF ( cloud_droplets ) THEN |
---|
| 399 | ! |
---|
[1010] | 400 | !-- Liquid water content, change in liquid water content |
---|
| 401 | #if defined( __nopointer ) |
---|
| 402 | ALLOCATE ( ql(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 403 | ql_c(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 404 | #else |
---|
[667] | 405 | ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
[1010] | 406 | ql_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 407 | #endif |
---|
| 408 | ! |
---|
| 409 | !-- Real volume of particles (with weighting), volume of particles |
---|
| 410 | ALLOCATE ( ql_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
[667] | 411 | ql_vp(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[1] | 412 | ENDIF |
---|
| 413 | |
---|
| 414 | ENDIF |
---|
| 415 | |
---|
| 416 | ENDIF |
---|
| 417 | |
---|
[94] | 418 | IF ( ocean ) THEN |
---|
[1001] | 419 | ALLOCATE( saswsb(nysg:nyng,nxlg:nxrg), & |
---|
| 420 | saswst(nysg:nyng,nxlg:nxrg) ) |
---|
[1010] | 421 | #if defined( __nopointer ) |
---|
| 422 | ALLOCATE( prho(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 423 | rho(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 424 | sa(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 425 | sa_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 426 | tsa_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 427 | #else |
---|
[667] | 428 | ALLOCATE( prho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 429 | rho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 430 | sa_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 431 | sa_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 432 | sa_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[388] | 433 | prho => prho_1 |
---|
| 434 | rho => rho_1 ! routines calc_mean_profile and diffusion_e require |
---|
| 435 | ! density to be apointer |
---|
[1010] | 436 | #endif |
---|
[108] | 437 | IF ( humidity_remote ) THEN |
---|
[667] | 438 | ALLOCATE( qswst_remote(nysg:nyng,nxlg:nxrg)) |
---|
[108] | 439 | qswst_remote = 0.0 |
---|
| 440 | ENDIF |
---|
[94] | 441 | ENDIF |
---|
| 442 | |
---|
[1] | 443 | ! |
---|
| 444 | !-- 3D-array for storing the dissipation, needed for calculating the sgs |
---|
| 445 | !-- particle velocities |
---|
[1171] | 446 | IF ( use_sgs_for_particles .OR. wang_kernel .OR. turbulence .OR. & |
---|
| 447 | num_acc_per_node > 0 ) THEN |
---|
[1153] | 448 | ALLOCATE( diss(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[1] | 449 | ENDIF |
---|
| 450 | |
---|
| 451 | IF ( dt_dosp /= 9999999.9 ) THEN |
---|
| 452 | ALLOCATE( spectrum_x( 1:nx/2, 1:10, 1:10 ), & |
---|
| 453 | spectrum_y( 1:ny/2, 1:10, 1:10 ) ) |
---|
[146] | 454 | spectrum_x = 0.0 |
---|
| 455 | spectrum_y = 0.0 |
---|
[1] | 456 | ENDIF |
---|
| 457 | |
---|
| 458 | ! |
---|
[1299] | 459 | !-- 1D-array for large scale subsidence velocity |
---|
[1316] | 460 | ALLOCATE ( w_subs(nzb:nzt+1) ) |
---|
| 461 | w_subs = 0.0 |
---|
[1299] | 462 | |
---|
[1316] | 463 | |
---|
[1299] | 464 | ! |
---|
[138] | 465 | !-- 3D-arrays for the leaf area density and the canopy drag coefficient |
---|
| 466 | IF ( plant_canopy ) THEN |
---|
[667] | 467 | ALLOCATE ( lad_s(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 468 | lad_u(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 469 | lad_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 470 | lad_w(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 471 | cdc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[153] | 472 | |
---|
| 473 | IF ( passive_scalar ) THEN |
---|
[996] | 474 | ALLOCATE ( sls(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
[667] | 475 | sec(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[153] | 476 | ENDIF |
---|
| 477 | |
---|
| 478 | IF ( cthf /= 0.0 ) THEN |
---|
[996] | 479 | ALLOCATE ( lai(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
[667] | 480 | canopy_heat_flux(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
[153] | 481 | ENDIF |
---|
| 482 | |
---|
[138] | 483 | ENDIF |
---|
| 484 | |
---|
| 485 | ! |
---|
[51] | 486 | !-- 4D-array for storing the Rif-values at vertical walls |
---|
| 487 | IF ( topography /= 'flat' ) THEN |
---|
[667] | 488 | ALLOCATE( rif_wall(nzb:nzt+1,nysg:nyng,nxlg:nxrg,1:4) ) |
---|
[51] | 489 | rif_wall = 0.0 |
---|
| 490 | ENDIF |
---|
| 491 | |
---|
| 492 | ! |
---|
[106] | 493 | !-- Arrays to store velocity data from t-dt and the phase speeds which |
---|
| 494 | !-- are needed for radiation boundary conditions |
---|
[73] | 495 | IF ( outflow_l ) THEN |
---|
[667] | 496 | ALLOCATE( u_m_l(nzb:nzt+1,nysg:nyng,1:2), & |
---|
| 497 | v_m_l(nzb:nzt+1,nysg:nyng,0:1), & |
---|
| 498 | w_m_l(nzb:nzt+1,nysg:nyng,0:1) ) |
---|
[73] | 499 | ENDIF |
---|
| 500 | IF ( outflow_r ) THEN |
---|
[667] | 501 | ALLOCATE( u_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx), & |
---|
| 502 | v_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx), & |
---|
| 503 | w_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx) ) |
---|
[73] | 504 | ENDIF |
---|
[106] | 505 | IF ( outflow_l .OR. outflow_r ) THEN |
---|
[667] | 506 | ALLOCATE( c_u(nzb:nzt+1,nysg:nyng), c_v(nzb:nzt+1,nysg:nyng), & |
---|
| 507 | c_w(nzb:nzt+1,nysg:nyng) ) |
---|
[106] | 508 | ENDIF |
---|
[73] | 509 | IF ( outflow_s ) THEN |
---|
[667] | 510 | ALLOCATE( u_m_s(nzb:nzt+1,0:1,nxlg:nxrg), & |
---|
| 511 | v_m_s(nzb:nzt+1,1:2,nxlg:nxrg), & |
---|
| 512 | w_m_s(nzb:nzt+1,0:1,nxlg:nxrg) ) |
---|
[73] | 513 | ENDIF |
---|
| 514 | IF ( outflow_n ) THEN |
---|
[667] | 515 | ALLOCATE( u_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg), & |
---|
| 516 | v_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg), & |
---|
| 517 | w_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg) ) |
---|
[73] | 518 | ENDIF |
---|
[106] | 519 | IF ( outflow_s .OR. outflow_n ) THEN |
---|
[667] | 520 | ALLOCATE( c_u(nzb:nzt+1,nxlg:nxrg), c_v(nzb:nzt+1,nxlg:nxrg), & |
---|
| 521 | c_w(nzb:nzt+1,nxlg:nxrg) ) |
---|
[106] | 522 | ENDIF |
---|
[996] | 523 | IF ( outflow_l .OR. outflow_r .OR. outflow_s .OR. outflow_n ) THEN |
---|
[978] | 524 | ALLOCATE( c_u_m_l(nzb:nzt+1), c_v_m_l(nzb:nzt+1), c_w_m_l(nzb:nzt+1) ) |
---|
| 525 | ALLOCATE( c_u_m(nzb:nzt+1), c_v_m(nzb:nzt+1), c_w_m(nzb:nzt+1) ) |
---|
| 526 | ENDIF |
---|
[73] | 527 | |
---|
[978] | 528 | |
---|
[1010] | 529 | #if ! defined( __nopointer ) |
---|
[73] | 530 | ! |
---|
[1] | 531 | !-- Initial assignment of the pointers |
---|
[1001] | 532 | e => e_1; e_p => e_2; te_m => e_3 |
---|
[1032] | 533 | IF ( .NOT. neutral ) THEN |
---|
| 534 | pt => pt_1; pt_p => pt_2; tpt_m => pt_3 |
---|
| 535 | ELSE |
---|
| 536 | pt => pt_1; pt_p => pt_1; tpt_m => pt_3 |
---|
| 537 | ENDIF |
---|
[1001] | 538 | u => u_1; u_p => u_2; tu_m => u_3 |
---|
| 539 | v => v_1; v_p => v_2; tv_m => v_3 |
---|
| 540 | w => w_1; w_p => w_2; tw_m => w_3 |
---|
[1] | 541 | |
---|
[1001] | 542 | IF ( humidity .OR. passive_scalar ) THEN |
---|
| 543 | q => q_1; q_p => q_2; tq_m => q_3 |
---|
[1053] | 544 | IF ( humidity ) THEN |
---|
| 545 | vpt => vpt_1 |
---|
| 546 | IF ( cloud_physics ) THEN |
---|
| 547 | ql => ql_1 |
---|
| 548 | IF ( icloud_scheme == 0 ) THEN |
---|
[1115] | 549 | qc => qc_1 |
---|
| 550 | IF ( precipitation ) THEN |
---|
| 551 | qr => qr_1; qr_p => qr_2; tqr_m => qr_3 |
---|
| 552 | nr => nr_1; nr_p => nr_2; tnr_m => nr_3 |
---|
| 553 | ENDIF |
---|
[1053] | 554 | ENDIF |
---|
| 555 | ENDIF |
---|
| 556 | ENDIF |
---|
[1001] | 557 | IF ( cloud_droplets ) THEN |
---|
| 558 | ql => ql_1 |
---|
| 559 | ql_c => ql_2 |
---|
[1] | 560 | ENDIF |
---|
[1001] | 561 | ENDIF |
---|
[1] | 562 | |
---|
[1001] | 563 | IF ( ocean ) THEN |
---|
| 564 | sa => sa_1; sa_p => sa_2; tsa_m => sa_3 |
---|
| 565 | ENDIF |
---|
[1010] | 566 | #endif |
---|
[1] | 567 | |
---|
| 568 | ! |
---|
[709] | 569 | !-- Allocate arrays containing the RK coefficient for calculation of |
---|
| 570 | !-- perturbation pressure and turbulent fluxes. At this point values are |
---|
| 571 | !-- set for pressure calculation during initialization (where no timestep |
---|
| 572 | !-- is done). Further below the values needed within the timestep scheme |
---|
| 573 | !-- will be set. |
---|
| 574 | ALLOCATE( weight_substep(1:intermediate_timestep_count_max), & |
---|
[673] | 575 | weight_pres(1:intermediate_timestep_count_max) ) |
---|
[709] | 576 | weight_substep = 1.0 |
---|
| 577 | weight_pres = 1.0 |
---|
| 578 | intermediate_timestep_count = 1 ! needed when simulated_time = 0.0 |
---|
[673] | 579 | |
---|
| 580 | ! |
---|
[1] | 581 | !-- Initialize model variables |
---|
[147] | 582 | IF ( TRIM( initializing_actions ) /= 'read_restart_data' .AND. & |
---|
[328] | 583 | TRIM( initializing_actions ) /= 'cyclic_fill' ) THEN |
---|
[1] | 584 | ! |
---|
| 585 | !-- First model run of a possible job queue. |
---|
| 586 | !-- Initial profiles of the variables must be computes. |
---|
| 587 | IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 ) THEN |
---|
| 588 | ! |
---|
| 589 | !-- Use solutions of the 1D model as initial profiles, |
---|
| 590 | !-- start 1D model |
---|
| 591 | CALL init_1d_model |
---|
| 592 | ! |
---|
| 593 | !-- Transfer initial profiles to the arrays of the 3D model |
---|
[667] | 594 | DO i = nxlg, nxrg |
---|
| 595 | DO j = nysg, nyng |
---|
[1] | 596 | e(:,j,i) = e1d |
---|
| 597 | kh(:,j,i) = kh1d |
---|
| 598 | km(:,j,i) = km1d |
---|
| 599 | pt(:,j,i) = pt_init |
---|
| 600 | u(:,j,i) = u1d |
---|
| 601 | v(:,j,i) = v1d |
---|
| 602 | ENDDO |
---|
| 603 | ENDDO |
---|
| 604 | |
---|
[75] | 605 | IF ( humidity .OR. passive_scalar ) THEN |
---|
[667] | 606 | DO i = nxlg, nxrg |
---|
| 607 | DO j = nysg, nyng |
---|
[1] | 608 | q(:,j,i) = q_init |
---|
| 609 | ENDDO |
---|
| 610 | ENDDO |
---|
[1115] | 611 | IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & |
---|
| 612 | precipitation ) THEN |
---|
[1053] | 613 | DO i = nxlg, nxrg |
---|
| 614 | DO j = nysg, nyng |
---|
[1115] | 615 | qr(:,j,i) = 0.0 |
---|
| 616 | nr(:,j,i) = 0.0 |
---|
[1053] | 617 | ENDDO |
---|
| 618 | ENDDO |
---|
[1115] | 619 | ! |
---|
| 620 | !-- Initialze nc_1d with default value |
---|
| 621 | nc_1d(:) = nc_const |
---|
| 622 | |
---|
[1053] | 623 | ENDIF |
---|
[1] | 624 | ENDIF |
---|
| 625 | |
---|
| 626 | IF ( .NOT. constant_diffusion ) THEN |
---|
[667] | 627 | DO i = nxlg, nxrg |
---|
| 628 | DO j = nysg, nyng |
---|
[1] | 629 | e(:,j,i) = e1d |
---|
| 630 | ENDDO |
---|
| 631 | ENDDO |
---|
| 632 | ! |
---|
| 633 | !-- Store initial profiles for output purposes etc. |
---|
| 634 | hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 ) |
---|
| 635 | |
---|
| 636 | IF ( prandtl_layer ) THEN |
---|
| 637 | rif = rif1d(nzb+1) |
---|
| 638 | ts = 0.0 ! could actually be computed more accurately in the |
---|
| 639 | ! 1D model. Update when opportunity arises. |
---|
| 640 | us = us1d |
---|
| 641 | usws = usws1d |
---|
| 642 | vsws = vsws1d |
---|
| 643 | ELSE |
---|
| 644 | ts = 0.0 ! must be set, because used in |
---|
| 645 | rif = 0.0 ! flowste |
---|
| 646 | us = 0.0 |
---|
| 647 | usws = 0.0 |
---|
| 648 | vsws = 0.0 |
---|
| 649 | ENDIF |
---|
| 650 | |
---|
| 651 | ELSE |
---|
| 652 | e = 0.0 ! must be set, because used in |
---|
| 653 | rif = 0.0 ! flowste |
---|
| 654 | ts = 0.0 |
---|
| 655 | us = 0.0 |
---|
| 656 | usws = 0.0 |
---|
| 657 | vsws = 0.0 |
---|
| 658 | ENDIF |
---|
[102] | 659 | uswst = top_momentumflux_u |
---|
| 660 | vswst = top_momentumflux_v |
---|
[1] | 661 | |
---|
| 662 | ! |
---|
| 663 | !-- In every case qs = 0.0 (see also pt) |
---|
| 664 | !-- This could actually be computed more accurately in the 1D model. |
---|
| 665 | !-- Update when opportunity arises! |
---|
[1053] | 666 | IF ( humidity .OR. passive_scalar ) THEN |
---|
| 667 | qs = 0.0 |
---|
[1115] | 668 | IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & |
---|
| 669 | precipitation ) THEN |
---|
[1053] | 670 | qrs = 0.0 |
---|
| 671 | nrs = 0.0 |
---|
| 672 | ENDIF |
---|
| 673 | ENDIF |
---|
[1] | 674 | |
---|
| 675 | ! |
---|
| 676 | !-- inside buildings set velocities back to zero |
---|
| 677 | IF ( topography /= 'flat' ) THEN |
---|
| 678 | DO i = nxl-1, nxr+1 |
---|
| 679 | DO j = nys-1, nyn+1 |
---|
| 680 | u(nzb:nzb_u_inner(j,i),j,i) = 0.0 |
---|
| 681 | v(nzb:nzb_v_inner(j,i),j,i) = 0.0 |
---|
| 682 | ENDDO |
---|
| 683 | ENDDO |
---|
[667] | 684 | |
---|
[1] | 685 | ! |
---|
| 686 | !-- WARNING: The extra boundary conditions set after running the |
---|
| 687 | !-- ------- 1D model impose an error on the divergence one layer |
---|
| 688 | !-- below the topography; need to correct later |
---|
| 689 | !-- ATTENTION: Provisional correction for Piacsek & Williams |
---|
| 690 | !-- --------- advection scheme: keep u and v zero one layer below |
---|
| 691 | !-- the topography. |
---|
[667] | 692 | IF ( ibc_uv_b == 1 ) THEN |
---|
| 693 | ! |
---|
[1] | 694 | !-- Neumann condition |
---|
| 695 | DO i = nxl-1, nxr+1 |
---|
| 696 | DO j = nys-1, nyn+1 |
---|
| 697 | IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = u(1,j,i) |
---|
| 698 | IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = v(1,j,i) |
---|
| 699 | ENDDO |
---|
| 700 | ENDDO |
---|
| 701 | |
---|
| 702 | ENDIF |
---|
| 703 | |
---|
| 704 | ENDIF |
---|
| 705 | |
---|
| 706 | ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 ) & |
---|
| 707 | THEN |
---|
[1241] | 708 | |
---|
[1] | 709 | ! |
---|
[1241] | 710 | !-- Overwrite initial profiles in case of nudging |
---|
| 711 | IF ( nudging ) THEN |
---|
| 712 | pt_init = ptnudge(:,1) |
---|
| 713 | u_init = unudge(:,1) |
---|
| 714 | v_init = vnudge(:,1) |
---|
| 715 | IF ( humidity .OR. passive_scalar ) THEN |
---|
| 716 | q_init = qnudge(:,1) |
---|
| 717 | ENDIF |
---|
| 718 | |
---|
| 719 | WRITE( message_string, * ) 'Initial profiles of u, v and ', & |
---|
| 720 | 'scalars from NUDGING_DATA are used.' |
---|
| 721 | CALL message( 'init_3d_model', 'PA0370', 0, 0, 0, 6, 0 ) |
---|
| 722 | ENDIF |
---|
| 723 | |
---|
| 724 | ! |
---|
[1] | 725 | !-- Use constructed initial profiles (velocity constant with height, |
---|
| 726 | !-- temperature profile with constant gradient) |
---|
[667] | 727 | DO i = nxlg, nxrg |
---|
| 728 | DO j = nysg, nyng |
---|
[1] | 729 | pt(:,j,i) = pt_init |
---|
| 730 | u(:,j,i) = u_init |
---|
| 731 | v(:,j,i) = v_init |
---|
| 732 | ENDDO |
---|
| 733 | ENDDO |
---|
[75] | 734 | |
---|
[1] | 735 | ! |
---|
[292] | 736 | !-- Set initial horizontal velocities at the lowest computational grid |
---|
| 737 | !-- levels to zero in order to avoid too small time steps caused by the |
---|
| 738 | !-- diffusion limit in the initial phase of a run (at k=1, dz/2 occurs |
---|
| 739 | !-- in the limiting formula!). The original values are stored to be later |
---|
| 740 | !-- used for volume flow control. |
---|
[667] | 741 | DO i = nxlg, nxrg |
---|
| 742 | DO j = nysg, nyng |
---|
[1] | 743 | u(nzb:nzb_u_inner(j,i)+1,j,i) = 0.0 |
---|
| 744 | v(nzb:nzb_v_inner(j,i)+1,j,i) = 0.0 |
---|
| 745 | ENDDO |
---|
| 746 | ENDDO |
---|
| 747 | |
---|
[75] | 748 | IF ( humidity .OR. passive_scalar ) THEN |
---|
[667] | 749 | DO i = nxlg, nxrg |
---|
| 750 | DO j = nysg, nyng |
---|
[1] | 751 | q(:,j,i) = q_init |
---|
| 752 | ENDDO |
---|
| 753 | ENDDO |
---|
[1115] | 754 | IF ( cloud_physics .AND. icloud_scheme == 0 ) THEN |
---|
| 755 | ! |
---|
| 756 | !-- Initialze nc_1d with default value |
---|
| 757 | nc_1d(:) = nc_const |
---|
| 758 | |
---|
| 759 | IF ( precipitation ) THEN |
---|
| 760 | DO i = nxlg, nxrg |
---|
| 761 | DO j = nysg, nyng |
---|
| 762 | qr(:,j,i) = 0.0 |
---|
| 763 | nr(:,j,i) = 0.0 |
---|
| 764 | ENDDO |
---|
[1053] | 765 | ENDDO |
---|
[1115] | 766 | ENDIF |
---|
| 767 | |
---|
[1053] | 768 | ENDIF |
---|
[1] | 769 | ENDIF |
---|
| 770 | |
---|
[94] | 771 | IF ( ocean ) THEN |
---|
[667] | 772 | DO i = nxlg, nxrg |
---|
| 773 | DO j = nysg, nyng |
---|
[94] | 774 | sa(:,j,i) = sa_init |
---|
| 775 | ENDDO |
---|
| 776 | ENDDO |
---|
| 777 | ENDIF |
---|
[1] | 778 | |
---|
| 779 | IF ( constant_diffusion ) THEN |
---|
| 780 | km = km_constant |
---|
| 781 | kh = km / prandtl_number |
---|
[108] | 782 | e = 0.0 |
---|
| 783 | ELSEIF ( e_init > 0.0 ) THEN |
---|
| 784 | DO k = nzb+1, nzt |
---|
| 785 | km(k,:,:) = 0.1 * l_grid(k) * SQRT( e_init ) |
---|
| 786 | ENDDO |
---|
| 787 | km(nzb,:,:) = km(nzb+1,:,:) |
---|
| 788 | km(nzt+1,:,:) = km(nzt,:,:) |
---|
| 789 | kh = km / prandtl_number |
---|
| 790 | e = e_init |
---|
[1] | 791 | ELSE |
---|
[108] | 792 | IF ( .NOT. ocean ) THEN |
---|
| 793 | kh = 0.01 ! there must exist an initial diffusion, because |
---|
| 794 | km = 0.01 ! otherwise no TKE would be produced by the |
---|
| 795 | ! production terms, as long as not yet |
---|
| 796 | ! e = (u*/cm)**2 at k=nzb+1 |
---|
| 797 | ELSE |
---|
| 798 | kh = 0.00001 |
---|
| 799 | km = 0.00001 |
---|
| 800 | ENDIF |
---|
| 801 | e = 0.0 |
---|
[1] | 802 | ENDIF |
---|
[102] | 803 | rif = 0.0 |
---|
| 804 | ts = 0.0 |
---|
| 805 | us = 0.0 |
---|
| 806 | usws = 0.0 |
---|
| 807 | uswst = top_momentumflux_u |
---|
| 808 | vsws = 0.0 |
---|
| 809 | vswst = top_momentumflux_v |
---|
[75] | 810 | IF ( humidity .OR. passive_scalar ) qs = 0.0 |
---|
[1] | 811 | |
---|
| 812 | ! |
---|
| 813 | !-- Compute initial temperature field and other constants used in case |
---|
| 814 | !-- of a sloping surface |
---|
| 815 | IF ( sloping_surface ) CALL init_slope |
---|
| 816 | |
---|
[46] | 817 | ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 ) & |
---|
| 818 | THEN |
---|
| 819 | ! |
---|
| 820 | !-- Initialization will completely be done by the user |
---|
| 821 | CALL user_init_3d_model |
---|
| 822 | |
---|
[1] | 823 | ENDIF |
---|
[667] | 824 | ! |
---|
| 825 | !-- Bottom boundary |
---|
| 826 | IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2 ) THEN |
---|
| 827 | u(nzb,:,:) = 0.0 |
---|
| 828 | v(nzb,:,:) = 0.0 |
---|
| 829 | ENDIF |
---|
[1] | 830 | |
---|
| 831 | ! |
---|
[151] | 832 | !-- Apply channel flow boundary condition |
---|
[132] | 833 | IF ( TRIM( bc_uv_t ) == 'dirichlet_0' ) THEN |
---|
| 834 | u(nzt+1,:,:) = 0.0 |
---|
| 835 | v(nzt+1,:,:) = 0.0 |
---|
| 836 | ENDIF |
---|
| 837 | |
---|
| 838 | ! |
---|
[1] | 839 | !-- Calculate virtual potential temperature |
---|
[75] | 840 | IF ( humidity ) vpt = pt * ( 1.0 + 0.61 * q ) |
---|
[1] | 841 | |
---|
| 842 | ! |
---|
| 843 | !-- Store initial profiles for output purposes etc. |
---|
| 844 | hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 ) |
---|
| 845 | hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 ) |
---|
[667] | 846 | IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2) THEN |
---|
| 847 | hom(nzb,1,5,:) = 0.0 |
---|
| 848 | hom(nzb,1,6,:) = 0.0 |
---|
[1] | 849 | ENDIF |
---|
| 850 | hom(:,1,7,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 ) |
---|
| 851 | hom(:,1,23,:) = SPREAD( km(:,nys,nxl), 2, statistic_regions+1 ) |
---|
| 852 | hom(:,1,24,:) = SPREAD( kh(:,nys,nxl), 2, statistic_regions+1 ) |
---|
| 853 | |
---|
[97] | 854 | IF ( ocean ) THEN |
---|
| 855 | ! |
---|
| 856 | !-- Store initial salinity profile |
---|
| 857 | hom(:,1,26,:) = SPREAD( sa(:,nys,nxl), 2, statistic_regions+1 ) |
---|
| 858 | ENDIF |
---|
[1] | 859 | |
---|
[75] | 860 | IF ( humidity ) THEN |
---|
[1] | 861 | ! |
---|
| 862 | !-- Store initial profile of total water content, virtual potential |
---|
| 863 | !-- temperature |
---|
| 864 | hom(:,1,26,:) = SPREAD( q(:,nys,nxl), 2, statistic_regions+1 ) |
---|
| 865 | hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 ) |
---|
| 866 | IF ( cloud_physics .OR. cloud_droplets ) THEN |
---|
| 867 | ! |
---|
| 868 | !-- Store initial profile of specific humidity and potential |
---|
| 869 | !-- temperature |
---|
| 870 | hom(:,1,27,:) = SPREAD( q(:,nys,nxl), 2, statistic_regions+1 ) |
---|
| 871 | hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 ) |
---|
| 872 | ENDIF |
---|
| 873 | ENDIF |
---|
| 874 | |
---|
| 875 | IF ( passive_scalar ) THEN |
---|
| 876 | ! |
---|
| 877 | !-- Store initial scalar profile |
---|
| 878 | hom(:,1,26,:) = SPREAD( q(:,nys,nxl), 2, statistic_regions+1 ) |
---|
| 879 | ENDIF |
---|
| 880 | |
---|
| 881 | ! |
---|
[19] | 882 | !-- Initialize fluxes at bottom surface |
---|
[1] | 883 | IF ( use_surface_fluxes ) THEN |
---|
| 884 | |
---|
| 885 | IF ( constant_heatflux ) THEN |
---|
| 886 | ! |
---|
| 887 | !-- Heat flux is prescribed |
---|
| 888 | IF ( random_heatflux ) THEN |
---|
| 889 | CALL disturb_heatflux |
---|
| 890 | ELSE |
---|
| 891 | shf = surface_heatflux |
---|
| 892 | ! |
---|
[1241] | 893 | !-- Initialize shf with data from external file LSF_DATA |
---|
| 894 | IF ( large_scale_forcing .AND. lsf_surf ) THEN |
---|
| 895 | CALL ls_forcing_surf ( simulated_time ) |
---|
| 896 | ENDIF |
---|
| 897 | |
---|
| 898 | ! |
---|
[1] | 899 | !-- Over topography surface_heatflux is replaced by wall_heatflux(0) |
---|
| 900 | IF ( TRIM( topography ) /= 'flat' ) THEN |
---|
[667] | 901 | DO i = nxlg, nxrg |
---|
| 902 | DO j = nysg, nyng |
---|
[1] | 903 | IF ( nzb_s_inner(j,i) /= 0 ) THEN |
---|
| 904 | shf(j,i) = wall_heatflux(0) |
---|
| 905 | ENDIF |
---|
| 906 | ENDDO |
---|
| 907 | ENDDO |
---|
| 908 | ENDIF |
---|
| 909 | ENDIF |
---|
| 910 | ENDIF |
---|
| 911 | |
---|
| 912 | ! |
---|
| 913 | !-- Determine the near-surface water flux |
---|
[75] | 914 | IF ( humidity .OR. passive_scalar ) THEN |
---|
[1115] | 915 | IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & |
---|
| 916 | precipitation ) THEN |
---|
| 917 | qrsws = 0.0 |
---|
| 918 | nrsws = 0.0 |
---|
[1053] | 919 | ENDIF |
---|
[1] | 920 | IF ( constant_waterflux ) THEN |
---|
| 921 | qsws = surface_waterflux |
---|
[407] | 922 | ! |
---|
| 923 | !-- Over topography surface_waterflux is replaced by |
---|
| 924 | !-- wall_humidityflux(0) |
---|
| 925 | IF ( TRIM( topography ) /= 'flat' ) THEN |
---|
| 926 | wall_qflux = wall_humidityflux |
---|
[667] | 927 | DO i = nxlg, nxrg |
---|
| 928 | DO j = nysg, nyng |
---|
[407] | 929 | IF ( nzb_s_inner(j,i) /= 0 ) THEN |
---|
| 930 | qsws(j,i) = wall_qflux(0) |
---|
| 931 | ENDIF |
---|
| 932 | ENDDO |
---|
| 933 | ENDDO |
---|
| 934 | ENDIF |
---|
[1] | 935 | ENDIF |
---|
| 936 | ENDIF |
---|
| 937 | |
---|
| 938 | ENDIF |
---|
| 939 | |
---|
| 940 | ! |
---|
[19] | 941 | !-- Initialize fluxes at top surface |
---|
[94] | 942 | !-- Currently, only the heatflux and salinity flux can be prescribed. |
---|
| 943 | !-- The latent flux is zero in this case! |
---|
[19] | 944 | IF ( use_top_fluxes ) THEN |
---|
| 945 | |
---|
| 946 | IF ( constant_top_heatflux ) THEN |
---|
| 947 | ! |
---|
| 948 | !-- Heat flux is prescribed |
---|
| 949 | tswst = top_heatflux |
---|
| 950 | |
---|
[1053] | 951 | IF ( humidity .OR. passive_scalar ) THEN |
---|
| 952 | qswst = 0.0 |
---|
[1115] | 953 | IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & |
---|
| 954 | precipitation ) THEN |
---|
[1053] | 955 | nrswst = 0.0 |
---|
| 956 | qrswst = 0.0 |
---|
| 957 | ENDIF |
---|
| 958 | ENDIF |
---|
[94] | 959 | |
---|
| 960 | IF ( ocean ) THEN |
---|
[95] | 961 | saswsb = bottom_salinityflux |
---|
[94] | 962 | saswst = top_salinityflux |
---|
| 963 | ENDIF |
---|
[102] | 964 | ENDIF |
---|
[19] | 965 | |
---|
[102] | 966 | ! |
---|
| 967 | !-- Initialization in case of a coupled model run |
---|
| 968 | IF ( coupling_mode == 'ocean_to_atmosphere' ) THEN |
---|
| 969 | tswst = 0.0 |
---|
| 970 | ENDIF |
---|
| 971 | |
---|
[19] | 972 | ENDIF |
---|
| 973 | |
---|
| 974 | ! |
---|
[1] | 975 | !-- Initialize Prandtl layer quantities |
---|
| 976 | IF ( prandtl_layer ) THEN |
---|
| 977 | |
---|
| 978 | z0 = roughness_length |
---|
[978] | 979 | z0h = z0h_factor * z0 |
---|
[1] | 980 | |
---|
| 981 | IF ( .NOT. constant_heatflux ) THEN |
---|
| 982 | ! |
---|
| 983 | !-- Surface temperature is prescribed. Here the heat flux cannot be |
---|
| 984 | !-- simply estimated, because therefore rif, u* and theta* would have |
---|
| 985 | !-- to be computed by iteration. This is why the heat flux is assumed |
---|
| 986 | !-- to be zero before the first time step. It approaches its correct |
---|
| 987 | !-- value in the course of the first few time steps. |
---|
| 988 | shf = 0.0 |
---|
| 989 | ENDIF |
---|
| 990 | |
---|
[75] | 991 | IF ( humidity .OR. passive_scalar ) THEN |
---|
[1001] | 992 | IF ( .NOT. constant_waterflux ) qsws = 0.0 |
---|
[1115] | 993 | IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & |
---|
| 994 | precipitation ) THEN |
---|
| 995 | qrsws = 0.0 |
---|
| 996 | nrsws = 0.0 |
---|
[1053] | 997 | ENDIF |
---|
[1] | 998 | ENDIF |
---|
| 999 | |
---|
| 1000 | ENDIF |
---|
| 1001 | |
---|
[1179] | 1002 | ! |
---|
| 1003 | !-- Set the reference state to be used in the buoyancy terms (for ocean runs |
---|
| 1004 | !-- the reference state will be set (overwritten) in init_ocean) |
---|
| 1005 | IF ( use_single_reference_value ) THEN |
---|
| 1006 | IF ( .NOT. humidity ) THEN |
---|
| 1007 | ref_state(:) = pt_reference |
---|
| 1008 | ELSE |
---|
| 1009 | ref_state(:) = vpt_reference |
---|
| 1010 | ENDIF |
---|
| 1011 | ELSE |
---|
| 1012 | IF ( .NOT. humidity ) THEN |
---|
| 1013 | ref_state(:) = pt_init(:) |
---|
| 1014 | ELSE |
---|
| 1015 | ref_state(:) = vpt(:,nys,nxl) |
---|
| 1016 | ENDIF |
---|
| 1017 | ENDIF |
---|
[152] | 1018 | |
---|
| 1019 | ! |
---|
[707] | 1020 | !-- For the moment, vertical velocity is zero |
---|
| 1021 | w = 0.0 |
---|
[1] | 1022 | |
---|
| 1023 | ! |
---|
| 1024 | !-- Initialize array sums (must be defined in first call of pres) |
---|
| 1025 | sums = 0.0 |
---|
| 1026 | |
---|
| 1027 | ! |
---|
[707] | 1028 | !-- In case of iterative solvers, p must get an initial value |
---|
| 1029 | IF ( psolver == 'multigrid' .OR. psolver == 'sor' ) p = 0.0 |
---|
| 1030 | |
---|
| 1031 | ! |
---|
[72] | 1032 | !-- Treating cloud physics, liquid water content and precipitation amount |
---|
| 1033 | !-- are zero at beginning of the simulation |
---|
| 1034 | IF ( cloud_physics ) THEN |
---|
| 1035 | ql = 0.0 |
---|
| 1036 | IF ( precipitation ) precipitation_amount = 0.0 |
---|
[1115] | 1037 | IF ( icloud_scheme == 0 ) THEN |
---|
| 1038 | qc = 0.0 |
---|
| 1039 | nc_1d = nc_const |
---|
| 1040 | ENDIF |
---|
[72] | 1041 | ENDIF |
---|
[673] | 1042 | ! |
---|
[1] | 1043 | !-- Impose vortex with vertical axis on the initial velocity profile |
---|
| 1044 | IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 ) THEN |
---|
| 1045 | CALL init_rankine |
---|
| 1046 | ENDIF |
---|
| 1047 | |
---|
| 1048 | ! |
---|
| 1049 | !-- Impose temperature anomaly (advection test only) |
---|
| 1050 | IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0 ) THEN |
---|
| 1051 | CALL init_pt_anomaly |
---|
| 1052 | ENDIF |
---|
| 1053 | |
---|
| 1054 | ! |
---|
| 1055 | !-- If required, change the surface temperature at the start of the 3D run |
---|
| 1056 | IF ( pt_surface_initial_change /= 0.0 ) THEN |
---|
| 1057 | pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change |
---|
| 1058 | ENDIF |
---|
| 1059 | |
---|
| 1060 | ! |
---|
| 1061 | !-- If required, change the surface humidity/scalar at the start of the 3D |
---|
| 1062 | !-- run |
---|
[75] | 1063 | IF ( ( humidity .OR. passive_scalar ) .AND. & |
---|
[1] | 1064 | q_surface_initial_change /= 0.0 ) THEN |
---|
| 1065 | q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change |
---|
| 1066 | ENDIF |
---|
| 1067 | ! |
---|
| 1068 | !-- Initialize the random number generator (from numerical recipes) |
---|
| 1069 | CALL random_function_ini |
---|
| 1070 | |
---|
| 1071 | ! |
---|
| 1072 | !-- Initialize old and new time levels. |
---|
[1001] | 1073 | te_m = 0.0; tpt_m = 0.0; tu_m = 0.0; tv_m = 0.0; tw_m = 0.0 |
---|
[1] | 1074 | e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w |
---|
| 1075 | |
---|
[75] | 1076 | IF ( humidity .OR. passive_scalar ) THEN |
---|
[1001] | 1077 | tq_m = 0.0 |
---|
[1] | 1078 | q_p = q |
---|
[1115] | 1079 | IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & |
---|
| 1080 | precipitation ) THEN |
---|
[1053] | 1081 | tqr_m = 0.0 |
---|
| 1082 | qr_p = qr |
---|
| 1083 | tnr_m = 0.0 |
---|
| 1084 | nr_p = nr |
---|
| 1085 | ENDIF |
---|
[1] | 1086 | ENDIF |
---|
| 1087 | |
---|
[94] | 1088 | IF ( ocean ) THEN |
---|
| 1089 | tsa_m = 0.0 |
---|
| 1090 | sa_p = sa |
---|
| 1091 | ENDIF |
---|
[667] | 1092 | |
---|
[94] | 1093 | |
---|
[147] | 1094 | ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data' .OR. & |
---|
[667] | 1095 | TRIM( initializing_actions ) == 'cyclic_fill' ) & |
---|
[1] | 1096 | THEN |
---|
| 1097 | ! |
---|
[767] | 1098 | !-- When reading data for cyclic fill of 3D prerun data files, read |
---|
| 1099 | !-- some of the global variables from the restart file which are required |
---|
| 1100 | !-- for initializing the inflow |
---|
[328] | 1101 | IF ( TRIM( initializing_actions ) == 'cyclic_fill' ) THEN |
---|
[559] | 1102 | |
---|
[759] | 1103 | DO i = 0, io_blocks-1 |
---|
| 1104 | IF ( i == io_group ) THEN |
---|
| 1105 | CALL read_parts_of_var_list |
---|
| 1106 | CALL close_file( 13 ) |
---|
| 1107 | ENDIF |
---|
| 1108 | #if defined( __parallel ) |
---|
| 1109 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 1110 | #endif |
---|
| 1111 | ENDDO |
---|
[328] | 1112 | |
---|
[767] | 1113 | ENDIF |
---|
| 1114 | |
---|
[151] | 1115 | ! |
---|
[767] | 1116 | !-- Read binary data from restart file |
---|
| 1117 | DO i = 0, io_blocks-1 |
---|
| 1118 | IF ( i == io_group ) THEN |
---|
| 1119 | CALL read_3d_binary |
---|
| 1120 | ENDIF |
---|
| 1121 | #if defined( __parallel ) |
---|
| 1122 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 1123 | #endif |
---|
| 1124 | ENDDO |
---|
| 1125 | |
---|
[328] | 1126 | ! |
---|
[767] | 1127 | !-- Initialization of the turbulence recycling method |
---|
| 1128 | IF ( TRIM( initializing_actions ) == 'cyclic_fill' .AND. & |
---|
| 1129 | turbulent_inflow ) THEN |
---|
| 1130 | ! |
---|
| 1131 | !-- First store the profiles to be used at the inflow. |
---|
| 1132 | !-- These profiles are the (temporally) and horizontally averaged vertical |
---|
| 1133 | !-- profiles from the prerun. Alternatively, prescribed profiles |
---|
| 1134 | !-- for u,v-components can be used. |
---|
| 1135 | ALLOCATE( mean_inflow_profiles(nzb:nzt+1,5) ) |
---|
[151] | 1136 | |
---|
[767] | 1137 | IF ( use_prescribed_profile_data ) THEN |
---|
| 1138 | mean_inflow_profiles(:,1) = u_init ! u |
---|
| 1139 | mean_inflow_profiles(:,2) = v_init ! v |
---|
| 1140 | ELSE |
---|
[328] | 1141 | mean_inflow_profiles(:,1) = hom_sum(:,1,0) ! u |
---|
| 1142 | mean_inflow_profiles(:,2) = hom_sum(:,2,0) ! v |
---|
[767] | 1143 | ENDIF |
---|
| 1144 | mean_inflow_profiles(:,4) = hom_sum(:,4,0) ! pt |
---|
| 1145 | mean_inflow_profiles(:,5) = hom_sum(:,8,0) ! e |
---|
[151] | 1146 | |
---|
| 1147 | ! |
---|
[767] | 1148 | !-- If necessary, adjust the horizontal flow field to the prescribed |
---|
| 1149 | !-- profiles |
---|
| 1150 | IF ( use_prescribed_profile_data ) THEN |
---|
| 1151 | DO i = nxlg, nxrg |
---|
[667] | 1152 | DO j = nysg, nyng |
---|
[328] | 1153 | DO k = nzb, nzt+1 |
---|
[767] | 1154 | u(k,j,i) = u(k,j,i) - hom_sum(k,1,0) + u_init(k) |
---|
| 1155 | v(k,j,i) = v(k,j,i) - hom_sum(k,2,0) + v_init(k) |
---|
[328] | 1156 | ENDDO |
---|
[151] | 1157 | ENDDO |
---|
[767] | 1158 | ENDDO |
---|
| 1159 | ENDIF |
---|
[151] | 1160 | |
---|
| 1161 | ! |
---|
[767] | 1162 | !-- Use these mean profiles at the inflow (provided that Dirichlet |
---|
| 1163 | !-- conditions are used) |
---|
| 1164 | IF ( inflow_l ) THEN |
---|
| 1165 | DO j = nysg, nyng |
---|
| 1166 | DO k = nzb, nzt+1 |
---|
| 1167 | u(k,j,nxlg:-1) = mean_inflow_profiles(k,1) |
---|
| 1168 | v(k,j,nxlg:-1) = mean_inflow_profiles(k,2) |
---|
| 1169 | w(k,j,nxlg:-1) = 0.0 |
---|
| 1170 | pt(k,j,nxlg:-1) = mean_inflow_profiles(k,4) |
---|
| 1171 | e(k,j,nxlg:-1) = mean_inflow_profiles(k,5) |
---|
| 1172 | ENDDO |
---|
| 1173 | ENDDO |
---|
| 1174 | ENDIF |
---|
| 1175 | |
---|
[151] | 1176 | ! |
---|
[767] | 1177 | !-- Calculate the damping factors to be used at the inflow. For a |
---|
| 1178 | !-- turbulent inflow the turbulent fluctuations have to be limited |
---|
| 1179 | !-- vertically because otherwise the turbulent inflow layer will grow |
---|
| 1180 | !-- in time. |
---|
| 1181 | IF ( inflow_damping_height == 9999999.9 ) THEN |
---|
| 1182 | ! |
---|
| 1183 | !-- Default: use the inversion height calculated by the prerun; if |
---|
| 1184 | !-- this is zero, inflow_damping_height must be explicitly |
---|
| 1185 | !-- specified. |
---|
| 1186 | IF ( hom_sum(nzb+6,pr_palm,0) /= 0.0 ) THEN |
---|
| 1187 | inflow_damping_height = hom_sum(nzb+6,pr_palm,0) |
---|
| 1188 | ELSE |
---|
| 1189 | WRITE( message_string, * ) 'inflow_damping_height must be ',& |
---|
| 1190 | 'explicitly specified because&the inversion height ', & |
---|
| 1191 | 'calculated by the prerun is zero.' |
---|
| 1192 | CALL message( 'init_3d_model', 'PA0318', 1, 2, 0, 6, 0 ) |
---|
[292] | 1193 | ENDIF |
---|
[151] | 1194 | |
---|
[767] | 1195 | ENDIF |
---|
| 1196 | |
---|
| 1197 | IF ( inflow_damping_width == 9999999.9 ) THEN |
---|
[151] | 1198 | ! |
---|
[767] | 1199 | !-- Default for the transition range: one tenth of the undamped |
---|
| 1200 | !-- layer |
---|
| 1201 | inflow_damping_width = 0.1 * inflow_damping_height |
---|
[151] | 1202 | |
---|
[767] | 1203 | ENDIF |
---|
[151] | 1204 | |
---|
[767] | 1205 | ALLOCATE( inflow_damping_factor(nzb:nzt+1) ) |
---|
[151] | 1206 | |
---|
[767] | 1207 | DO k = nzb, nzt+1 |
---|
[151] | 1208 | |
---|
[767] | 1209 | IF ( zu(k) <= inflow_damping_height ) THEN |
---|
| 1210 | inflow_damping_factor(k) = 1.0 |
---|
[996] | 1211 | ELSEIF ( zu(k) <= ( inflow_damping_height + inflow_damping_width ) ) THEN |
---|
| 1212 | inflow_damping_factor(k) = 1.0 - & |
---|
| 1213 | ( zu(k) - inflow_damping_height ) / & |
---|
| 1214 | inflow_damping_width |
---|
[767] | 1215 | ELSE |
---|
| 1216 | inflow_damping_factor(k) = 0.0 |
---|
| 1217 | ENDIF |
---|
[151] | 1218 | |
---|
[767] | 1219 | ENDDO |
---|
[151] | 1220 | |
---|
[147] | 1221 | ENDIF |
---|
| 1222 | |
---|
[152] | 1223 | ! |
---|
[359] | 1224 | !-- Inside buildings set velocities and TKE back to zero |
---|
| 1225 | IF ( TRIM( initializing_actions ) == 'cyclic_fill' .AND. & |
---|
| 1226 | topography /= 'flat' ) THEN |
---|
| 1227 | ! |
---|
| 1228 | !-- Inside buildings set velocities and TKE back to zero. |
---|
| 1229 | !-- Other scalars (pt, q, s, km, kh, p, sa, ...) are ignored at present, |
---|
| 1230 | !-- maybe revise later. |
---|
[1001] | 1231 | DO i = nxlg, nxrg |
---|
| 1232 | DO j = nysg, nyng |
---|
| 1233 | u (nzb:nzb_u_inner(j,i),j,i) = 0.0 |
---|
| 1234 | v (nzb:nzb_v_inner(j,i),j,i) = 0.0 |
---|
| 1235 | w (nzb:nzb_w_inner(j,i),j,i) = 0.0 |
---|
| 1236 | e (nzb:nzb_w_inner(j,i),j,i) = 0.0 |
---|
| 1237 | tu_m(nzb:nzb_u_inner(j,i),j,i) = 0.0 |
---|
| 1238 | tv_m(nzb:nzb_v_inner(j,i),j,i) = 0.0 |
---|
| 1239 | tw_m(nzb:nzb_w_inner(j,i),j,i) = 0.0 |
---|
| 1240 | te_m(nzb:nzb_w_inner(j,i),j,i) = 0.0 |
---|
| 1241 | tpt_m(nzb:nzb_w_inner(j,i),j,i) = 0.0 |
---|
[359] | 1242 | ENDDO |
---|
[1001] | 1243 | ENDDO |
---|
[359] | 1244 | |
---|
| 1245 | ENDIF |
---|
| 1246 | |
---|
| 1247 | ! |
---|
[1] | 1248 | !-- Calculate initial temperature field and other constants used in case |
---|
| 1249 | !-- of a sloping surface |
---|
| 1250 | IF ( sloping_surface ) CALL init_slope |
---|
| 1251 | |
---|
| 1252 | ! |
---|
| 1253 | !-- Initialize new time levels (only done in order to set boundary values |
---|
| 1254 | !-- including ghost points) |
---|
| 1255 | e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w |
---|
[1053] | 1256 | IF ( humidity .OR. passive_scalar ) THEN |
---|
| 1257 | q_p = q |
---|
[1115] | 1258 | IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & |
---|
| 1259 | precipitation ) THEN |
---|
[1053] | 1260 | qr_p = qr |
---|
| 1261 | nr_p = nr |
---|
| 1262 | ENDIF |
---|
| 1263 | ENDIF |
---|
[94] | 1264 | IF ( ocean ) sa_p = sa |
---|
[1] | 1265 | |
---|
[181] | 1266 | ! |
---|
| 1267 | !-- Allthough tendency arrays are set in prognostic_equations, they have |
---|
| 1268 | !-- have to be predefined here because they are used (but multiplied with 0) |
---|
| 1269 | !-- there before they are set. |
---|
[1001] | 1270 | te_m = 0.0; tpt_m = 0.0; tu_m = 0.0; tv_m = 0.0; tw_m = 0.0 |
---|
[1053] | 1271 | IF ( humidity .OR. passive_scalar ) THEN |
---|
| 1272 | tq_m = 0.0 |
---|
[1115] | 1273 | IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & |
---|
| 1274 | precipitation ) THEN |
---|
[1053] | 1275 | tqr_m = 0.0 |
---|
| 1276 | tnr_m = 0.0 |
---|
| 1277 | ENDIF |
---|
| 1278 | ENDIF |
---|
[1001] | 1279 | IF ( ocean ) tsa_m = 0.0 |
---|
[181] | 1280 | |
---|
[1] | 1281 | ELSE |
---|
| 1282 | ! |
---|
| 1283 | !-- Actually this part of the programm should not be reached |
---|
[254] | 1284 | message_string = 'unknown initializing problem' |
---|
| 1285 | CALL message( 'init_3d_model', 'PA0193', 1, 2, 0, 6, 0 ) |
---|
[1] | 1286 | ENDIF |
---|
| 1287 | |
---|
[151] | 1288 | |
---|
| 1289 | IF ( TRIM( initializing_actions ) /= 'read_restart_data' ) THEN |
---|
[1] | 1290 | ! |
---|
[151] | 1291 | !-- Initialize old timelevels needed for radiation boundary conditions |
---|
| 1292 | IF ( outflow_l ) THEN |
---|
| 1293 | u_m_l(:,:,:) = u(:,:,1:2) |
---|
| 1294 | v_m_l(:,:,:) = v(:,:,0:1) |
---|
| 1295 | w_m_l(:,:,:) = w(:,:,0:1) |
---|
| 1296 | ENDIF |
---|
| 1297 | IF ( outflow_r ) THEN |
---|
| 1298 | u_m_r(:,:,:) = u(:,:,nx-1:nx) |
---|
| 1299 | v_m_r(:,:,:) = v(:,:,nx-1:nx) |
---|
| 1300 | w_m_r(:,:,:) = w(:,:,nx-1:nx) |
---|
| 1301 | ENDIF |
---|
| 1302 | IF ( outflow_s ) THEN |
---|
| 1303 | u_m_s(:,:,:) = u(:,0:1,:) |
---|
| 1304 | v_m_s(:,:,:) = v(:,1:2,:) |
---|
| 1305 | w_m_s(:,:,:) = w(:,0:1,:) |
---|
| 1306 | ENDIF |
---|
| 1307 | IF ( outflow_n ) THEN |
---|
| 1308 | u_m_n(:,:,:) = u(:,ny-1:ny,:) |
---|
| 1309 | v_m_n(:,:,:) = v(:,ny-1:ny,:) |
---|
| 1310 | w_m_n(:,:,:) = w(:,ny-1:ny,:) |
---|
| 1311 | ENDIF |
---|
[667] | 1312 | |
---|
[151] | 1313 | ENDIF |
---|
[680] | 1314 | |
---|
[667] | 1315 | ! |
---|
| 1316 | !-- Calculate the initial volume flow at the right and north boundary |
---|
[709] | 1317 | IF ( conserve_volume_flow ) THEN |
---|
[151] | 1318 | |
---|
[767] | 1319 | IF ( use_prescribed_profile_data ) THEN |
---|
[667] | 1320 | |
---|
[732] | 1321 | volume_flow_initial_l = 0.0 |
---|
| 1322 | volume_flow_area_l = 0.0 |
---|
| 1323 | |
---|
[667] | 1324 | IF ( nxr == nx ) THEN |
---|
| 1325 | DO j = nys, nyn |
---|
[709] | 1326 | DO k = nzb_2d(j,nx)+1, nzt |
---|
[667] | 1327 | volume_flow_initial_l(1) = volume_flow_initial_l(1) + & |
---|
[767] | 1328 | u_init(k) * dzw(k) |
---|
| 1329 | volume_flow_area_l(1) = volume_flow_area_l(1) + dzw(k) |
---|
| 1330 | ENDDO |
---|
| 1331 | ENDDO |
---|
| 1332 | ENDIF |
---|
| 1333 | |
---|
| 1334 | IF ( nyn == ny ) THEN |
---|
| 1335 | DO i = nxl, nxr |
---|
| 1336 | DO k = nzb_2d(ny,i)+1, nzt |
---|
| 1337 | volume_flow_initial_l(2) = volume_flow_initial_l(2) + & |
---|
| 1338 | v_init(k) * dzw(k) |
---|
| 1339 | volume_flow_area_l(2) = volume_flow_area_l(2) + dzw(k) |
---|
| 1340 | ENDDO |
---|
| 1341 | ENDDO |
---|
| 1342 | ENDIF |
---|
| 1343 | |
---|
| 1344 | #if defined( __parallel ) |
---|
| 1345 | CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),& |
---|
| 1346 | 2, MPI_REAL, MPI_SUM, comm2d, ierr ) |
---|
| 1347 | CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1), & |
---|
| 1348 | 2, MPI_REAL, MPI_SUM, comm2d, ierr ) |
---|
| 1349 | |
---|
| 1350 | #else |
---|
| 1351 | volume_flow_initial = volume_flow_initial_l |
---|
| 1352 | volume_flow_area = volume_flow_area_l |
---|
| 1353 | #endif |
---|
| 1354 | |
---|
| 1355 | ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' ) THEN |
---|
| 1356 | |
---|
| 1357 | volume_flow_initial_l = 0.0 |
---|
| 1358 | volume_flow_area_l = 0.0 |
---|
| 1359 | |
---|
| 1360 | IF ( nxr == nx ) THEN |
---|
| 1361 | DO j = nys, nyn |
---|
| 1362 | DO k = nzb_2d(j,nx)+1, nzt |
---|
| 1363 | volume_flow_initial_l(1) = volume_flow_initial_l(1) + & |
---|
[667] | 1364 | hom_sum(k,1,0) * dzw(k) |
---|
| 1365 | volume_flow_area_l(1) = volume_flow_area_l(1) + dzw(k) |
---|
| 1366 | ENDDO |
---|
| 1367 | ENDDO |
---|
| 1368 | ENDIF |
---|
| 1369 | |
---|
| 1370 | IF ( nyn == ny ) THEN |
---|
| 1371 | DO i = nxl, nxr |
---|
[709] | 1372 | DO k = nzb_2d(ny,i)+1, nzt |
---|
[667] | 1373 | volume_flow_initial_l(2) = volume_flow_initial_l(2) + & |
---|
[709] | 1374 | hom_sum(k,2,0) * dzw(k) |
---|
[667] | 1375 | volume_flow_area_l(2) = volume_flow_area_l(2) + dzw(k) |
---|
| 1376 | ENDDO |
---|
| 1377 | ENDDO |
---|
| 1378 | ENDIF |
---|
| 1379 | |
---|
[732] | 1380 | #if defined( __parallel ) |
---|
| 1381 | CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),& |
---|
| 1382 | 2, MPI_REAL, MPI_SUM, comm2d, ierr ) |
---|
| 1383 | CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1), & |
---|
| 1384 | 2, MPI_REAL, MPI_SUM, comm2d, ierr ) |
---|
| 1385 | |
---|
| 1386 | #else |
---|
| 1387 | volume_flow_initial = volume_flow_initial_l |
---|
| 1388 | volume_flow_area = volume_flow_area_l |
---|
| 1389 | #endif |
---|
| 1390 | |
---|
[667] | 1391 | ELSEIF ( TRIM( initializing_actions ) /= 'read_restart_data' ) THEN |
---|
| 1392 | |
---|
[732] | 1393 | volume_flow_initial_l = 0.0 |
---|
| 1394 | volume_flow_area_l = 0.0 |
---|
| 1395 | |
---|
[667] | 1396 | IF ( nxr == nx ) THEN |
---|
| 1397 | DO j = nys, nyn |
---|
[709] | 1398 | DO k = nzb_2d(j,nx)+1, nzt |
---|
[667] | 1399 | volume_flow_initial_l(1) = volume_flow_initial_l(1) + & |
---|
[709] | 1400 | u(k,j,nx) * dzw(k) |
---|
[667] | 1401 | volume_flow_area_l(1) = volume_flow_area_l(1) + dzw(k) |
---|
| 1402 | ENDDO |
---|
| 1403 | ENDDO |
---|
| 1404 | ENDIF |
---|
| 1405 | |
---|
| 1406 | IF ( nyn == ny ) THEN |
---|
| 1407 | DO i = nxl, nxr |
---|
[709] | 1408 | DO k = nzb_2d(ny,i)+1, nzt |
---|
[667] | 1409 | volume_flow_initial_l(2) = volume_flow_initial_l(2) + & |
---|
| 1410 | v(k,ny,i) * dzw(k) |
---|
| 1411 | volume_flow_area_l(2) = volume_flow_area_l(2) + dzw(k) |
---|
| 1412 | ENDDO |
---|
| 1413 | ENDDO |
---|
| 1414 | ENDIF |
---|
| 1415 | |
---|
| 1416 | #if defined( __parallel ) |
---|
[732] | 1417 | CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),& |
---|
| 1418 | 2, MPI_REAL, MPI_SUM, comm2d, ierr ) |
---|
| 1419 | CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1), & |
---|
| 1420 | 2, MPI_REAL, MPI_SUM, comm2d, ierr ) |
---|
[667] | 1421 | |
---|
| 1422 | #else |
---|
[732] | 1423 | volume_flow_initial = volume_flow_initial_l |
---|
| 1424 | volume_flow_area = volume_flow_area_l |
---|
[667] | 1425 | #endif |
---|
| 1426 | |
---|
[732] | 1427 | ENDIF |
---|
| 1428 | |
---|
[151] | 1429 | ! |
---|
[709] | 1430 | !-- In case of 'bulk_velocity' mode, volume_flow_initial is calculated |
---|
| 1431 | !-- from u|v_bulk instead |
---|
[680] | 1432 | IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' ) THEN |
---|
| 1433 | volume_flow_initial(1) = u_bulk * volume_flow_area(1) |
---|
| 1434 | volume_flow_initial(2) = v_bulk * volume_flow_area(2) |
---|
| 1435 | ENDIF |
---|
[667] | 1436 | |
---|
[680] | 1437 | ENDIF |
---|
| 1438 | |
---|
[787] | 1439 | ! |
---|
| 1440 | !-- Initialize quantities for special advections schemes |
---|
| 1441 | CALL init_advec |
---|
[680] | 1442 | |
---|
[667] | 1443 | ! |
---|
[680] | 1444 | !-- Impose random perturbation on the horizontal velocity field and then |
---|
| 1445 | !-- remove the divergences from the velocity field at the initial stage |
---|
| 1446 | IF ( create_disturbances .AND. & |
---|
| 1447 | TRIM( initializing_actions ) /= 'read_restart_data' .AND. & |
---|
| 1448 | TRIM( initializing_actions ) /= 'cyclic_fill' ) THEN |
---|
| 1449 | |
---|
| 1450 | CALL disturb_field( nzb_u_inner, tend, u ) |
---|
| 1451 | CALL disturb_field( nzb_v_inner, tend, v ) |
---|
| 1452 | n_sor = nsor_ini |
---|
[1221] | 1453 | !$acc data copyin( d, ddzu, ddzw, nzb_s_inner, nzb_u_inner ) & |
---|
| 1454 | !$acc copyin( nzb_v_inner, nzb_w_inner, p, rflags_s_inner, tend ) & |
---|
| 1455 | !$acc copyin( weight_pres, weight_substep ) & |
---|
| 1456 | !$acc copy( tri, tric, u, v, w ) |
---|
[680] | 1457 | CALL pres |
---|
[1111] | 1458 | !$acc end data |
---|
[680] | 1459 | n_sor = nsor |
---|
| 1460 | ENDIF |
---|
| 1461 | |
---|
| 1462 | ! |
---|
[138] | 1463 | !-- Initialization of the leaf area density |
---|
[709] | 1464 | IF ( plant_canopy ) THEN |
---|
[138] | 1465 | |
---|
| 1466 | SELECT CASE ( TRIM( canopy_mode ) ) |
---|
| 1467 | |
---|
| 1468 | CASE( 'block' ) |
---|
| 1469 | |
---|
[667] | 1470 | DO i = nxlg, nxrg |
---|
| 1471 | DO j = nysg, nyng |
---|
[138] | 1472 | lad_s(:,j,i) = lad(:) |
---|
| 1473 | cdc(:,j,i) = drag_coefficient |
---|
[709] | 1474 | IF ( passive_scalar ) THEN |
---|
[153] | 1475 | sls(:,j,i) = leaf_surface_concentration |
---|
| 1476 | sec(:,j,i) = scalar_exchange_coefficient |
---|
| 1477 | ENDIF |
---|
[138] | 1478 | ENDDO |
---|
| 1479 | ENDDO |
---|
| 1480 | |
---|
| 1481 | CASE DEFAULT |
---|
| 1482 | |
---|
| 1483 | ! |
---|
| 1484 | !-- The DEFAULT case is reached either if the parameter |
---|
| 1485 | !-- canopy mode contains a wrong character string or if the |
---|
| 1486 | !-- user has coded a special case in the user interface. |
---|
| 1487 | !-- There, the subroutine user_init_plant_canopy checks |
---|
| 1488 | !-- which of these two conditions applies. |
---|
| 1489 | CALL user_init_plant_canopy |
---|
| 1490 | |
---|
| 1491 | END SELECT |
---|
| 1492 | |
---|
[667] | 1493 | CALL exchange_horiz( lad_s, nbgp ) |
---|
| 1494 | CALL exchange_horiz( cdc, nbgp ) |
---|
[138] | 1495 | |
---|
[709] | 1496 | IF ( passive_scalar ) THEN |
---|
[667] | 1497 | CALL exchange_horiz( sls, nbgp ) |
---|
| 1498 | CALL exchange_horiz( sec, nbgp ) |
---|
[153] | 1499 | ENDIF |
---|
| 1500 | |
---|
| 1501 | ! |
---|
| 1502 | !-- Sharp boundaries of the plant canopy in horizontal directions |
---|
| 1503 | !-- In vertical direction the interpolation is retained, as the leaf |
---|
| 1504 | !-- area density is initialised by prescribing a vertical profile |
---|
| 1505 | !-- consisting of piecewise linear segments. The upper boundary |
---|
| 1506 | !-- of the plant canopy is now defined by lad_w(pch_index,:,:) = 0.0. |
---|
| 1507 | |
---|
[138] | 1508 | DO i = nxl, nxr |
---|
| 1509 | DO j = nys, nyn |
---|
| 1510 | DO k = nzb, nzt+1 |
---|
[709] | 1511 | IF ( lad_s(k,j,i) > 0.0 ) THEN |
---|
[153] | 1512 | lad_u(k,j,i) = lad_s(k,j,i) |
---|
| 1513 | lad_u(k,j,i+1) = lad_s(k,j,i) |
---|
| 1514 | lad_v(k,j,i) = lad_s(k,j,i) |
---|
| 1515 | lad_v(k,j+1,i) = lad_s(k,j,i) |
---|
| 1516 | ENDIF |
---|
[138] | 1517 | ENDDO |
---|
| 1518 | DO k = nzb, nzt |
---|
| 1519 | lad_w(k,j,i) = 0.5 * ( lad_s(k+1,j,i) + lad_s(k,j,i) ) |
---|
| 1520 | ENDDO |
---|
| 1521 | ENDDO |
---|
| 1522 | ENDDO |
---|
| 1523 | |
---|
[153] | 1524 | lad_w(pch_index,:,:) = 0.0 |
---|
| 1525 | lad_w(nzt+1,:,:) = lad_w(nzt,:,:) |
---|
[138] | 1526 | |
---|
[667] | 1527 | CALL exchange_horiz( lad_u, nbgp ) |
---|
| 1528 | CALL exchange_horiz( lad_v, nbgp ) |
---|
| 1529 | CALL exchange_horiz( lad_w, nbgp ) |
---|
[153] | 1530 | |
---|
| 1531 | ! |
---|
| 1532 | !-- Initialisation of the canopy heat source distribution |
---|
[709] | 1533 | IF ( cthf /= 0.0 ) THEN |
---|
[153] | 1534 | ! |
---|
| 1535 | !-- Piecewise evaluation of the leaf area index by |
---|
| 1536 | !-- integration of the leaf area density |
---|
| 1537 | lai(:,:,:) = 0.0 |
---|
[667] | 1538 | DO i = nxlg, nxrg |
---|
| 1539 | DO j = nysg, nyng |
---|
[153] | 1540 | DO k = pch_index-1, 0, -1 |
---|
| 1541 | lai(k,j,i) = lai(k+1,j,i) + & |
---|
| 1542 | ( 0.5 * ( lad_w(k+1,j,i) + & |
---|
| 1543 | lad_s(k+1,j,i) ) * & |
---|
| 1544 | ( zw(k+1) - zu(k+1) ) ) + & |
---|
| 1545 | ( 0.5 * ( lad_w(k,j,i) + & |
---|
| 1546 | lad_s(k+1,j,i) ) * & |
---|
| 1547 | ( zu(k+1) - zw(k) ) ) |
---|
| 1548 | ENDDO |
---|
| 1549 | ENDDO |
---|
| 1550 | ENDDO |
---|
| 1551 | |
---|
| 1552 | ! |
---|
| 1553 | !-- Evaluation of the upward kinematic vertical heat flux within the |
---|
| 1554 | !-- canopy |
---|
[667] | 1555 | DO i = nxlg, nxrg |
---|
| 1556 | DO j = nysg, nyng |
---|
[153] | 1557 | DO k = 0, pch_index |
---|
| 1558 | canopy_heat_flux(k,j,i) = cthf * & |
---|
| 1559 | exp( -0.6 * lai(k,j,i) ) |
---|
| 1560 | ENDDO |
---|
| 1561 | ENDDO |
---|
| 1562 | ENDDO |
---|
| 1563 | |
---|
| 1564 | ! |
---|
| 1565 | !-- The near surface heat flux is derived from the heat flux |
---|
| 1566 | !-- distribution within the canopy |
---|
| 1567 | shf(:,:) = canopy_heat_flux(0,:,:) |
---|
| 1568 | |
---|
| 1569 | ENDIF |
---|
| 1570 | |
---|
[138] | 1571 | ENDIF |
---|
| 1572 | |
---|
| 1573 | ! |
---|
[1] | 1574 | !-- If required, initialize dvrp-software |
---|
| 1575 | IF ( dt_dvrp /= 9999999.9 ) CALL init_dvrp |
---|
| 1576 | |
---|
[96] | 1577 | IF ( ocean ) THEN |
---|
[1] | 1578 | ! |
---|
[96] | 1579 | !-- Initialize quantities needed for the ocean model |
---|
| 1580 | CALL init_ocean |
---|
[388] | 1581 | |
---|
[96] | 1582 | ELSE |
---|
| 1583 | ! |
---|
| 1584 | !-- Initialize quantities for handling cloud physics |
---|
[849] | 1585 | !-- This routine must be called before lpm_init, because |
---|
[96] | 1586 | !-- otherwise, array pt_d_t, needed in data_output_dvrp (called by |
---|
[849] | 1587 | !-- lpm_init) is not defined. |
---|
[96] | 1588 | CALL init_cloud_physics |
---|
| 1589 | ENDIF |
---|
[1] | 1590 | |
---|
| 1591 | ! |
---|
| 1592 | !-- If required, initialize particles |
---|
[849] | 1593 | IF ( particle_advection ) CALL lpm_init |
---|
[1] | 1594 | |
---|
| 1595 | ! |
---|
[673] | 1596 | !-- Initialize the ws-scheme. |
---|
| 1597 | IF ( ws_scheme_sca .OR. ws_scheme_mom ) CALL ws_init |
---|
[1] | 1598 | |
---|
| 1599 | ! |
---|
[709] | 1600 | !-- Setting weighting factors for calculation of perturbation pressure |
---|
| 1601 | !-- and turbulent quantities from the RK substeps |
---|
| 1602 | IF ( TRIM(timestep_scheme) == 'runge-kutta-3' ) THEN ! for RK3-method |
---|
| 1603 | |
---|
[713] | 1604 | weight_substep(1) = 1./6. |
---|
| 1605 | weight_substep(2) = 3./10. |
---|
| 1606 | weight_substep(3) = 8./15. |
---|
[709] | 1607 | |
---|
[713] | 1608 | weight_pres(1) = 1./3. |
---|
| 1609 | weight_pres(2) = 5./12. |
---|
| 1610 | weight_pres(3) = 1./4. |
---|
[709] | 1611 | |
---|
| 1612 | ELSEIF ( TRIM(timestep_scheme) == 'runge-kutta-2' ) THEN ! for RK2-method |
---|
| 1613 | |
---|
[713] | 1614 | weight_substep(1) = 1./2. |
---|
| 1615 | weight_substep(2) = 1./2. |
---|
[673] | 1616 | |
---|
[713] | 1617 | weight_pres(1) = 1./2. |
---|
| 1618 | weight_pres(2) = 1./2. |
---|
[709] | 1619 | |
---|
[1001] | 1620 | ELSE ! for Euler-method |
---|
[709] | 1621 | |
---|
[673] | 1622 | weight_substep(1) = 1.0 |
---|
[709] | 1623 | weight_pres(1) = 1.0 |
---|
| 1624 | |
---|
[673] | 1625 | ENDIF |
---|
| 1626 | |
---|
| 1627 | ! |
---|
[1] | 1628 | !-- Initialize Rayleigh damping factors |
---|
[785] | 1629 | rdf = 0.0 |
---|
| 1630 | rdf_sc = 0.0 |
---|
[1] | 1631 | IF ( rayleigh_damping_factor /= 0.0 ) THEN |
---|
[108] | 1632 | IF ( .NOT. ocean ) THEN |
---|
| 1633 | DO k = nzb+1, nzt |
---|
| 1634 | IF ( zu(k) >= rayleigh_damping_height ) THEN |
---|
| 1635 | rdf(k) = rayleigh_damping_factor * & |
---|
[1] | 1636 | ( SIN( pi * 0.5 * ( zu(k) - rayleigh_damping_height ) & |
---|
| 1637 | / ( zu(nzt) - rayleigh_damping_height ) )& |
---|
| 1638 | )**2 |
---|
[108] | 1639 | ENDIF |
---|
| 1640 | ENDDO |
---|
| 1641 | ELSE |
---|
| 1642 | DO k = nzt, nzb+1, -1 |
---|
| 1643 | IF ( zu(k) <= rayleigh_damping_height ) THEN |
---|
| 1644 | rdf(k) = rayleigh_damping_factor * & |
---|
| 1645 | ( SIN( pi * 0.5 * ( rayleigh_damping_height - zu(k) ) & |
---|
| 1646 | / ( rayleigh_damping_height - zu(nzb+1)))& |
---|
| 1647 | )**2 |
---|
| 1648 | ENDIF |
---|
| 1649 | ENDDO |
---|
| 1650 | ENDIF |
---|
[1] | 1651 | ENDIF |
---|
[785] | 1652 | IF ( scalar_rayleigh_damping ) rdf_sc = rdf |
---|
[1] | 1653 | |
---|
| 1654 | ! |
---|
[240] | 1655 | !-- Initialize the starting level and the vertical smoothing factor used for |
---|
| 1656 | !-- the external pressure gradient |
---|
| 1657 | dp_smooth_factor = 1.0 |
---|
| 1658 | IF ( dp_external ) THEN |
---|
| 1659 | ! |
---|
| 1660 | !-- Set the starting level dp_level_ind_b only if it has not been set before |
---|
| 1661 | !-- (e.g. in init_grid). |
---|
| 1662 | IF ( dp_level_ind_b == 0 ) THEN |
---|
| 1663 | ind_array = MINLOC( ABS( dp_level_b - zu ) ) |
---|
| 1664 | dp_level_ind_b = ind_array(1) - 1 + nzb |
---|
| 1665 | ! MINLOC uses lower array bound 1 |
---|
| 1666 | ENDIF |
---|
| 1667 | IF ( dp_smooth ) THEN |
---|
| 1668 | dp_smooth_factor(:dp_level_ind_b) = 0.0 |
---|
| 1669 | DO k = dp_level_ind_b+1, nzt |
---|
| 1670 | dp_smooth_factor(k) = 0.5 * ( 1.0 + SIN( pi * & |
---|
| 1671 | ( REAL( k - dp_level_ind_b ) / & |
---|
| 1672 | REAL( nzt - dp_level_ind_b ) - 0.5 ) ) ) |
---|
| 1673 | ENDDO |
---|
| 1674 | ENDIF |
---|
| 1675 | ENDIF |
---|
| 1676 | |
---|
| 1677 | ! |
---|
[978] | 1678 | !-- Initialize damping zone for the potential temperature in case of |
---|
| 1679 | !-- non-cyclic lateral boundaries. The damping zone has the maximum value |
---|
| 1680 | !-- at the inflow boundary and decreases to zero at pt_damping_width. |
---|
| 1681 | ptdf_x = 0.0 |
---|
| 1682 | ptdf_y = 0.0 |
---|
[1159] | 1683 | IF ( bc_lr_dirrad ) THEN |
---|
[996] | 1684 | DO i = nxl, nxr |
---|
[978] | 1685 | IF ( ( i * dx ) < pt_damping_width ) THEN |
---|
| 1686 | ptdf_x(i) = pt_damping_factor * ( SIN( pi * 0.5 * & |
---|
| 1687 | REAL( pt_damping_width - i * dx ) / ( & |
---|
| 1688 | REAL( pt_damping_width ) ) ) )**2 |
---|
[73] | 1689 | ENDIF |
---|
| 1690 | ENDDO |
---|
[1159] | 1691 | ELSEIF ( bc_lr_raddir ) THEN |
---|
[996] | 1692 | DO i = nxl, nxr |
---|
[978] | 1693 | IF ( ( i * dx ) > ( nx * dx - pt_damping_width ) ) THEN |
---|
[996] | 1694 | ptdf_x(i) = pt_damping_factor * & |
---|
| 1695 | SIN( pi * 0.5 * ( ( i - nx ) * dx + pt_damping_width ) / & |
---|
| 1696 | REAL( pt_damping_width ) )**2 |
---|
[73] | 1697 | ENDIF |
---|
[978] | 1698 | ENDDO |
---|
[1159] | 1699 | ELSEIF ( bc_ns_dirrad ) THEN |
---|
[996] | 1700 | DO j = nys, nyn |
---|
[978] | 1701 | IF ( ( j * dy ) > ( ny * dy - pt_damping_width ) ) THEN |
---|
[996] | 1702 | ptdf_y(j) = pt_damping_factor * & |
---|
| 1703 | SIN( pi * 0.5 * ( ( j - ny ) * dy + pt_damping_width ) / & |
---|
| 1704 | REAL( pt_damping_width ) )**2 |
---|
[1] | 1705 | ENDIF |
---|
[978] | 1706 | ENDDO |
---|
[1159] | 1707 | ELSEIF ( bc_ns_raddir ) THEN |
---|
[996] | 1708 | DO j = nys, nyn |
---|
[978] | 1709 | IF ( ( j * dy ) < pt_damping_width ) THEN |
---|
[996] | 1710 | ptdf_y(j) = pt_damping_factor * & |
---|
| 1711 | SIN( pi * 0.5 * ( pt_damping_width - j * dy ) / & |
---|
| 1712 | REAL( pt_damping_width ) )**2 |
---|
[1] | 1713 | ENDIF |
---|
[73] | 1714 | ENDDO |
---|
[1] | 1715 | ENDIF |
---|
| 1716 | |
---|
| 1717 | ! |
---|
[709] | 1718 | !-- Initialize local summation arrays for routine flow_statistics. |
---|
| 1719 | !-- This is necessary because they may not yet have been initialized when they |
---|
| 1720 | !-- are called from flow_statistics (or - depending on the chosen model run - |
---|
| 1721 | !-- are never initialized) |
---|
[1] | 1722 | sums_divnew_l = 0.0 |
---|
| 1723 | sums_divold_l = 0.0 |
---|
| 1724 | sums_l_l = 0.0 |
---|
| 1725 | sums_up_fraction_l = 0.0 |
---|
| 1726 | sums_wsts_bc_l = 0.0 |
---|
| 1727 | |
---|
| 1728 | ! |
---|
| 1729 | !-- Pre-set masks for regional statistics. Default is the total model domain. |
---|
[1015] | 1730 | !-- Ghost points are excluded because counting values at the ghost boundaries |
---|
| 1731 | !-- would bias the statistics |
---|
[1] | 1732 | rmask = 1.0 |
---|
[1025] | 1733 | rmask(:,nxlg:nxl-1,:) = 0.0; rmask(:,nxr+1:nxrg,:) = 0.0 |
---|
| 1734 | rmask(nysg:nys-1,:,:) = 0.0; rmask(nyn+1:nyng,:,:) = 0.0 |
---|
[1] | 1735 | |
---|
| 1736 | ! |
---|
[51] | 1737 | !-- User-defined initializing actions. Check afterwards, if maximum number |
---|
[709] | 1738 | !-- of allowed timeseries is exceeded |
---|
[1] | 1739 | CALL user_init |
---|
| 1740 | |
---|
[51] | 1741 | IF ( dots_num > dots_max ) THEN |
---|
[254] | 1742 | WRITE( message_string, * ) 'number of time series quantities exceeds', & |
---|
[274] | 1743 | ' its maximum of dots_max = ', dots_max, & |
---|
[254] | 1744 | ' &Please increase dots_max in modules.f90.' |
---|
| 1745 | CALL message( 'init_3d_model', 'PA0194', 1, 2, 0, 6, 0 ) |
---|
[51] | 1746 | ENDIF |
---|
| 1747 | |
---|
[1] | 1748 | ! |
---|
| 1749 | !-- Input binary data file is not needed anymore. This line must be placed |
---|
| 1750 | !-- after call of user_init! |
---|
| 1751 | CALL close_file( 13 ) |
---|
| 1752 | |
---|
| 1753 | ! |
---|
| 1754 | !-- Compute total sum of active mask grid points |
---|
| 1755 | !-- ngp_2dh: number of grid points of a horizontal cross section through the |
---|
| 1756 | !-- total domain |
---|
| 1757 | !-- ngp_3d: number of grid points of the total domain |
---|
[132] | 1758 | ngp_2dh_outer_l = 0 |
---|
| 1759 | ngp_2dh_outer = 0 |
---|
| 1760 | ngp_2dh_s_inner_l = 0 |
---|
| 1761 | ngp_2dh_s_inner = 0 |
---|
| 1762 | ngp_2dh_l = 0 |
---|
| 1763 | ngp_2dh = 0 |
---|
[485] | 1764 | ngp_3d_inner_l = 0.0 |
---|
[132] | 1765 | ngp_3d_inner = 0 |
---|
| 1766 | ngp_3d = 0 |
---|
| 1767 | ngp_sums = ( nz + 2 ) * ( pr_palm + max_pr_user ) |
---|
[1] | 1768 | |
---|
| 1769 | DO sr = 0, statistic_regions |
---|
| 1770 | DO i = nxl, nxr |
---|
| 1771 | DO j = nys, nyn |
---|
| 1772 | IF ( rmask(j,i,sr) == 1.0 ) THEN |
---|
| 1773 | ! |
---|
| 1774 | !-- All xy-grid points |
---|
| 1775 | ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1 |
---|
| 1776 | ! |
---|
| 1777 | !-- xy-grid points above topography |
---|
| 1778 | DO k = nzb_s_outer(j,i), nz + 1 |
---|
| 1779 | ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr) + 1 |
---|
| 1780 | ENDDO |
---|
[132] | 1781 | DO k = nzb_s_inner(j,i), nz + 1 |
---|
| 1782 | ngp_2dh_s_inner_l(k,sr) = ngp_2dh_s_inner_l(k,sr) + 1 |
---|
| 1783 | ENDDO |
---|
[1] | 1784 | ! |
---|
| 1785 | !-- All grid points of the total domain above topography |
---|
| 1786 | ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr) + & |
---|
| 1787 | ( nz - nzb_s_inner(j,i) + 2 ) |
---|
| 1788 | ENDIF |
---|
| 1789 | ENDDO |
---|
| 1790 | ENDDO |
---|
| 1791 | ENDDO |
---|
| 1792 | |
---|
| 1793 | sr = statistic_regions + 1 |
---|
| 1794 | #if defined( __parallel ) |
---|
[622] | 1795 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
[485] | 1796 | CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM, & |
---|
[1] | 1797 | comm2d, ierr ) |
---|
[622] | 1798 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
[485] | 1799 | CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr, & |
---|
[1] | 1800 | MPI_INTEGER, MPI_SUM, comm2d, ierr ) |
---|
[622] | 1801 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
[485] | 1802 | CALL MPI_ALLREDUCE( ngp_2dh_s_inner_l(0,0), ngp_2dh_s_inner(0,0), & |
---|
[132] | 1803 | (nz+2)*sr, MPI_INTEGER, MPI_SUM, comm2d, ierr ) |
---|
[622] | 1804 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
[485] | 1805 | CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner_tmp(0), sr, MPI_REAL, & |
---|
[1] | 1806 | MPI_SUM, comm2d, ierr ) |
---|
[485] | 1807 | ngp_3d_inner = INT( ngp_3d_inner_tmp, KIND = SELECTED_INT_KIND( 18 ) ) |
---|
[1] | 1808 | #else |
---|
[132] | 1809 | ngp_2dh = ngp_2dh_l |
---|
| 1810 | ngp_2dh_outer = ngp_2dh_outer_l |
---|
| 1811 | ngp_2dh_s_inner = ngp_2dh_s_inner_l |
---|
[485] | 1812 | ngp_3d_inner = INT( ngp_3d_inner_l, KIND = SELECTED_INT_KIND( 18 ) ) |
---|
[1] | 1813 | #endif |
---|
| 1814 | |
---|
[560] | 1815 | ngp_3d = INT ( ngp_2dh, KIND = SELECTED_INT_KIND( 18 ) ) * & |
---|
| 1816 | INT ( (nz + 2 ), KIND = SELECTED_INT_KIND( 18 ) ) |
---|
[1] | 1817 | |
---|
| 1818 | ! |
---|
| 1819 | !-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics, |
---|
| 1820 | !-- buoyancy, etc. A zero value will occur for cases where all grid points of |
---|
| 1821 | !-- the respective subdomain lie below the surface topography |
---|
[667] | 1822 | ngp_2dh_outer = MAX( 1, ngp_2dh_outer(:,:) ) |
---|
[631] | 1823 | ngp_3d_inner = MAX( INT(1, KIND = SELECTED_INT_KIND( 18 )), & |
---|
| 1824 | ngp_3d_inner(:) ) |
---|
[667] | 1825 | ngp_2dh_s_inner = MAX( 1, ngp_2dh_s_inner(:,:) ) |
---|
[1] | 1826 | |
---|
[485] | 1827 | DEALLOCATE( ngp_2dh_l, ngp_2dh_outer_l, ngp_3d_inner_l, ngp_3d_inner_tmp ) |
---|
[1] | 1828 | |
---|
| 1829 | |
---|
| 1830 | END SUBROUTINE init_3d_model |
---|