1 | !> @file prognostic_equations.f90 |
---|
2 | !--------------------------------------------------------------------------------! |
---|
3 | ! This file is part of PALM. |
---|
4 | ! |
---|
5 | ! PALM is free software: you can redistribute it and/or modify it under the terms |
---|
6 | ! of the GNU General Public License as published by the Free Software Foundation, |
---|
7 | ! either version 3 of the License, or (at your option) any later version. |
---|
8 | ! |
---|
9 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
10 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
11 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
12 | ! |
---|
13 | ! You should have received a copy of the GNU General Public License along with |
---|
14 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
15 | ! |
---|
16 | ! Copyright 1997-2016 Leibniz Universitaet Hannover |
---|
17 | !--------------------------------------------------------------------------------! |
---|
18 | ! |
---|
19 | ! Current revisions: |
---|
20 | ! ------------------ |
---|
21 | ! |
---|
22 | ! |
---|
23 | ! Former revisions: |
---|
24 | ! ----------------- |
---|
25 | ! $Id: prognostic_equations.f90 1914 2016-05-26 14:44:07Z suehring $ |
---|
26 | ! |
---|
27 | ! 1873 2016-04-18 14:50:06Z maronga |
---|
28 | ! Module renamed (removed _mod) |
---|
29 | ! |
---|
30 | ! 1850 2016-04-08 13:29:27Z maronga |
---|
31 | ! Module renamed |
---|
32 | ! |
---|
33 | ! 1826 2016-04-07 12:01:39Z maronga |
---|
34 | ! Renamed canopy model calls. |
---|
35 | ! |
---|
36 | ! 1822 2016-04-07 07:49:42Z hoffmann |
---|
37 | ! Kessler microphysics scheme moved to microphysics. |
---|
38 | |
---|
39 | ! 1819 2016-04-06 21:05:15Z witha |
---|
40 | ! Added calls for wind turbine model |
---|
41 | ! |
---|
42 | ! 1757 2016-02-22 15:49:32Z maronga |
---|
43 | ! |
---|
44 | ! |
---|
45 | ! 1691 2015-10-26 16:17:44Z maronga |
---|
46 | ! Added optional model spin-up without radiation / land surface model calls. |
---|
47 | ! Formatting corrections. |
---|
48 | ! |
---|
49 | ! 1682 2015-10-07 23:56:08Z knoop |
---|
50 | ! Code annotations made doxygen readable |
---|
51 | ! |
---|
52 | ! 1585 2015-04-30 07:05:52Z maronga |
---|
53 | ! Added call for temperature tendency calculation due to radiative flux divergence |
---|
54 | ! |
---|
55 | ! 1517 2015-01-07 19:12:25Z hoffmann |
---|
56 | ! advec_s_bc_mod addded, since advec_s_bc is now a module |
---|
57 | ! |
---|
58 | ! 1496 2014-12-02 17:25:50Z maronga |
---|
59 | ! Renamed "radiation" -> "cloud_top_radiation" |
---|
60 | ! |
---|
61 | ! 1484 2014-10-21 10:53:05Z kanani |
---|
62 | ! Changes due to new module structure of the plant canopy model: |
---|
63 | ! parameters cthf and plant_canopy moved to module plant_canopy_model_mod. |
---|
64 | ! Removed double-listing of use_upstream_for_tke in ONLY-list of module |
---|
65 | ! control_parameters |
---|
66 | ! |
---|
67 | ! 1409 2014-05-23 12:11:32Z suehring |
---|
68 | ! Bugfix: i_omp_start changed for advec_u_ws at left inflow and outflow boundary. |
---|
69 | ! This ensures that left-hand side fluxes are also calculated for nxl in that |
---|
70 | ! case, even though the solution at nxl is overwritten in boundary_conds() |
---|
71 | ! |
---|
72 | ! 1398 2014-05-07 11:15:00Z heinze |
---|
73 | ! Rayleigh-damping for horizontal velocity components changed: instead of damping |
---|
74 | ! against ug and vg, damping against u_init and v_init is used to allow for a |
---|
75 | ! homogenized treatment in case of nudging |
---|
76 | ! |
---|
77 | ! 1380 2014-04-28 12:40:45Z heinze |
---|
78 | ! Change order of calls for scalar prognostic quantities: |
---|
79 | ! ls_advec -> nudging -> subsidence since initial profiles |
---|
80 | ! |
---|
81 | ! 1374 2014-04-25 12:55:07Z raasch |
---|
82 | ! missing variables added to ONLY lists |
---|
83 | ! |
---|
84 | ! 1365 2014-04-22 15:03:56Z boeske |
---|
85 | ! Calls of ls_advec for large scale advection added, |
---|
86 | ! subroutine subsidence is only called if use_subsidence_tendencies = .F., |
---|
87 | ! new argument ls_index added to the calls of subsidence |
---|
88 | ! +ls_index |
---|
89 | ! |
---|
90 | ! 1361 2014-04-16 15:17:48Z hoffmann |
---|
91 | ! Two-moment microphysics moved to the start of prognostic equations. This makes |
---|
92 | ! the 3d arrays for tend_q, tend_qr, tend_pt and tend_pt redundant. |
---|
93 | ! Additionally, it is allowed to call the microphysics just once during the time |
---|
94 | ! step (not at each sub-time step). |
---|
95 | ! |
---|
96 | ! Two-moment cloud physics added for vector and accelerator optimization. |
---|
97 | ! |
---|
98 | ! 1353 2014-04-08 15:21:23Z heinze |
---|
99 | ! REAL constants provided with KIND-attribute |
---|
100 | ! |
---|
101 | ! 1337 2014-03-25 15:11:48Z heinze |
---|
102 | ! Bugfix: REAL constants provided with KIND-attribute |
---|
103 | ! |
---|
104 | ! 1332 2014-03-25 11:59:43Z suehring |
---|
105 | ! Bugfix: call advec_ws or advec_pw for TKE only if NOT use_upstream_for_tke |
---|
106 | ! |
---|
107 | ! 1330 2014-03-24 17:29:32Z suehring |
---|
108 | ! In case of SGS-particle velocity advection of TKE is also allowed with |
---|
109 | ! dissipative 5th-order scheme. |
---|
110 | ! |
---|
111 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
112 | ! ONLY-attribute added to USE-statements, |
---|
113 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
114 | ! kinds are defined in new module kinds, |
---|
115 | ! old module precision_kind is removed, |
---|
116 | ! revision history before 2012 removed, |
---|
117 | ! comment fields (!:) to be used for variable explanations added to |
---|
118 | ! all variable declaration statements |
---|
119 | ! |
---|
120 | ! 1318 2014-03-17 13:35:16Z raasch |
---|
121 | ! module interfaces removed |
---|
122 | ! |
---|
123 | ! 1257 2013-11-08 15:18:40Z raasch |
---|
124 | ! openacc loop vector clauses removed, independent clauses added |
---|
125 | ! |
---|
126 | ! 1246 2013-11-01 08:59:45Z heinze |
---|
127 | ! enable nudging also for accelerator version |
---|
128 | ! |
---|
129 | ! 1241 2013-10-30 11:36:58Z heinze |
---|
130 | ! usage of nudging enabled (so far not implemented for accelerator version) |
---|
131 | ! |
---|
132 | ! 1179 2013-06-14 05:57:58Z raasch |
---|
133 | ! two arguments removed from routine buoyancy, ref_state updated on device |
---|
134 | ! |
---|
135 | ! 1128 2013-04-12 06:19:32Z raasch |
---|
136 | ! those parts requiring global communication moved to time_integration, |
---|
137 | ! loop index bounds in accelerator version replaced by i_left, i_right, j_south, |
---|
138 | ! j_north |
---|
139 | ! |
---|
140 | ! 1115 2013-03-26 18:16:16Z hoffmann |
---|
141 | ! optimized cloud physics: calculation of microphysical tendencies transfered |
---|
142 | ! to microphysics.f90; qr and nr are only calculated if precipitation is required |
---|
143 | ! |
---|
144 | ! 1111 2013-03-08 23:54:10Z raasch |
---|
145 | ! update directives for prognostic quantities removed |
---|
146 | ! |
---|
147 | ! 1106 2013-03-04 05:31:38Z raasch |
---|
148 | ! small changes in code formatting |
---|
149 | ! |
---|
150 | ! 1092 2013-02-02 11:24:22Z raasch |
---|
151 | ! unused variables removed |
---|
152 | ! |
---|
153 | ! 1053 2012-11-13 17:11:03Z hoffmann |
---|
154 | ! implementation of two new prognostic equations for rain drop concentration (nr) |
---|
155 | ! and rain water content (qr) |
---|
156 | ! |
---|
157 | ! currently, only available for cache loop optimization |
---|
158 | ! |
---|
159 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
160 | ! code put under GPL (PALM 3.9) |
---|
161 | ! |
---|
162 | ! 1019 2012-09-28 06:46:45Z raasch |
---|
163 | ! non-optimized version of prognostic_equations removed |
---|
164 | ! |
---|
165 | ! 1015 2012-09-27 09:23:24Z raasch |
---|
166 | ! new branch prognostic_equations_acc |
---|
167 | ! OpenACC statements added + code changes required for GPU optimization |
---|
168 | ! |
---|
169 | ! 1001 2012-09-13 14:08:46Z raasch |
---|
170 | ! all actions concerning leapfrog- and upstream-spline-scheme removed |
---|
171 | ! |
---|
172 | ! 978 2012-08-09 08:28:32Z fricke |
---|
173 | ! km_damp_x and km_damp_y removed in calls of diffusion_u and diffusion_v |
---|
174 | ! add ptdf_x, ptdf_y for damping the potential temperature at the inflow |
---|
175 | ! boundary in case of non-cyclic lateral boundaries |
---|
176 | ! Bugfix: first thread index changes for WS-scheme at the inflow |
---|
177 | ! |
---|
178 | ! 940 2012-07-09 14:31:00Z raasch |
---|
179 | ! temperature equation can be switched off |
---|
180 | ! |
---|
181 | ! Revision 1.1 2000/04/13 14:56:27 schroeter |
---|
182 | ! Initial revision |
---|
183 | ! |
---|
184 | ! |
---|
185 | ! Description: |
---|
186 | ! ------------ |
---|
187 | !> Solving the prognostic equations. |
---|
188 | !------------------------------------------------------------------------------! |
---|
189 | MODULE prognostic_equations_mod |
---|
190 | |
---|
191 | |
---|
192 | |
---|
193 | USE arrays_3d, & |
---|
194 | ONLY: diss_l_e, diss_l_nr, diss_l_pt, diss_l_q, diss_l_qr, & |
---|
195 | diss_l_sa, diss_s_e, diss_s_nr, diss_s_pt, diss_s_q, & |
---|
196 | diss_s_qr, diss_s_sa, e, e_p, flux_s_e, flux_s_nr, flux_s_pt, & |
---|
197 | flux_s_q, flux_s_qr, flux_s_sa, flux_l_e, flux_l_nr, & |
---|
198 | flux_l_pt, flux_l_q, flux_l_qr, flux_l_sa, nr, nr_p, nrsws, & |
---|
199 | nrswst, pt, ptdf_x, ptdf_y, pt_init, pt_p, prho, q, q_init, & |
---|
200 | q_p, qsws, qswst, qr, qr_p, qrsws, qrswst, rdf, rdf_sc, & |
---|
201 | ref_state, rho, sa, sa_init, sa_p, saswsb, saswst, shf, tend, & |
---|
202 | te_m, tnr_m, tpt_m, tq_m, tqr_m, tsa_m, tswst, tu_m, tv_m, & |
---|
203 | tw_m, u, ug, u_init, u_p, v, vg, vpt, v_init, v_p, w, w_p |
---|
204 | |
---|
205 | USE control_parameters, & |
---|
206 | ONLY: call_microphysics_at_all_substeps, cloud_physics, & |
---|
207 | cloud_top_radiation, constant_diffusion, dp_external, & |
---|
208 | dp_level_ind_b, dp_smooth_factor, dpdxy, dt_3d, humidity, & |
---|
209 | inflow_l, intermediate_timestep_count, & |
---|
210 | intermediate_timestep_count_max, large_scale_forcing, & |
---|
211 | large_scale_subsidence, microphysics_seifert, & |
---|
212 | microphysics_sat_adjust, neutral, nudging, ocean, outflow_l, & |
---|
213 | outflow_s, passive_scalar, prho_reference, prho_reference, & |
---|
214 | prho_reference, pt_reference, pt_reference, pt_reference, & |
---|
215 | scalar_advec, scalar_advec, simulated_time, sloping_surface, & |
---|
216 | timestep_scheme, tsc, use_subsidence_tendencies, & |
---|
217 | use_upstream_for_tke, wall_heatflux, & |
---|
218 | wall_nrflux, wall_qflux, wall_qflux, wall_qflux, wall_qrflux, & |
---|
219 | wall_salinityflux, ws_scheme_mom, ws_scheme_sca |
---|
220 | |
---|
221 | USE cpulog, & |
---|
222 | ONLY: cpu_log, log_point |
---|
223 | |
---|
224 | USE eqn_state_seawater_mod, & |
---|
225 | ONLY: eqn_state_seawater |
---|
226 | |
---|
227 | USE indices, & |
---|
228 | ONLY: i_left, i_right, j_north, j_south, nxl, nxlu, nxr, nyn, nys, & |
---|
229 | nysv, nzb_s_inner, nzb_u_inner, nzb_v_inner, nzb_w_inner, nzt |
---|
230 | |
---|
231 | USE advec_ws, & |
---|
232 | ONLY: advec_s_ws, advec_s_ws_acc, advec_u_ws, advec_u_ws_acc, & |
---|
233 | advec_v_ws, advec_v_ws_acc, advec_w_ws, advec_w_ws_acc |
---|
234 | |
---|
235 | USE advec_s_bc_mod, & |
---|
236 | ONLY: advec_s_bc |
---|
237 | |
---|
238 | USE advec_s_pw_mod, & |
---|
239 | ONLY: advec_s_pw |
---|
240 | |
---|
241 | USE advec_s_up_mod, & |
---|
242 | ONLY: advec_s_up |
---|
243 | |
---|
244 | USE advec_u_pw_mod, & |
---|
245 | ONLY: advec_u_pw |
---|
246 | |
---|
247 | USE advec_u_up_mod, & |
---|
248 | ONLY: advec_u_up |
---|
249 | |
---|
250 | USE advec_v_pw_mod, & |
---|
251 | ONLY: advec_v_pw |
---|
252 | |
---|
253 | USE advec_v_up_mod, & |
---|
254 | ONLY: advec_v_up |
---|
255 | |
---|
256 | USE advec_w_pw_mod, & |
---|
257 | ONLY: advec_w_pw |
---|
258 | |
---|
259 | USE advec_w_up_mod, & |
---|
260 | ONLY: advec_w_up |
---|
261 | |
---|
262 | USE buoyancy_mod, & |
---|
263 | ONLY: buoyancy, buoyancy_acc |
---|
264 | |
---|
265 | USE calc_radiation_mod, & |
---|
266 | ONLY: calc_radiation |
---|
267 | |
---|
268 | USE coriolis_mod, & |
---|
269 | ONLY: coriolis, coriolis_acc |
---|
270 | |
---|
271 | USE diffusion_e_mod, & |
---|
272 | ONLY: diffusion_e, diffusion_e_acc |
---|
273 | |
---|
274 | USE diffusion_s_mod, & |
---|
275 | ONLY: diffusion_s, diffusion_s_acc |
---|
276 | |
---|
277 | USE diffusion_u_mod, & |
---|
278 | ONLY: diffusion_u, diffusion_u_acc |
---|
279 | |
---|
280 | USE diffusion_v_mod, & |
---|
281 | ONLY: diffusion_v, diffusion_v_acc |
---|
282 | |
---|
283 | USE diffusion_w_mod, & |
---|
284 | ONLY: diffusion_w, diffusion_w_acc |
---|
285 | |
---|
286 | USE kinds |
---|
287 | |
---|
288 | USE ls_forcing_mod, & |
---|
289 | ONLY: ls_advec |
---|
290 | |
---|
291 | USE microphysics_mod, & |
---|
292 | ONLY: microphysics_control |
---|
293 | |
---|
294 | USE nudge_mod, & |
---|
295 | ONLY: nudge |
---|
296 | |
---|
297 | USE plant_canopy_model_mod, & |
---|
298 | ONLY: cthf, plant_canopy, pcm_tendency |
---|
299 | |
---|
300 | USE production_e_mod, & |
---|
301 | ONLY: production_e, production_e_acc |
---|
302 | |
---|
303 | USE radiation_model_mod, & |
---|
304 | ONLY: radiation, radiation_scheme, radiation_tendency, & |
---|
305 | skip_time_do_radiation |
---|
306 | |
---|
307 | USE statistics, & |
---|
308 | ONLY: hom |
---|
309 | |
---|
310 | USE subsidence_mod, & |
---|
311 | ONLY: subsidence |
---|
312 | |
---|
313 | USE user_actions_mod, & |
---|
314 | ONLY: user_actions |
---|
315 | |
---|
316 | USE wind_turbine_model_mod, & |
---|
317 | ONLY: wind_turbine, wtm_tendencies |
---|
318 | |
---|
319 | |
---|
320 | PRIVATE |
---|
321 | PUBLIC prognostic_equations_cache, prognostic_equations_vector, & |
---|
322 | prognostic_equations_acc |
---|
323 | |
---|
324 | INTERFACE prognostic_equations_cache |
---|
325 | MODULE PROCEDURE prognostic_equations_cache |
---|
326 | END INTERFACE prognostic_equations_cache |
---|
327 | |
---|
328 | INTERFACE prognostic_equations_vector |
---|
329 | MODULE PROCEDURE prognostic_equations_vector |
---|
330 | END INTERFACE prognostic_equations_vector |
---|
331 | |
---|
332 | INTERFACE prognostic_equations_acc |
---|
333 | MODULE PROCEDURE prognostic_equations_acc |
---|
334 | END INTERFACE prognostic_equations_acc |
---|
335 | |
---|
336 | |
---|
337 | CONTAINS |
---|
338 | |
---|
339 | |
---|
340 | !------------------------------------------------------------------------------! |
---|
341 | ! Description: |
---|
342 | ! ------------ |
---|
343 | !> Version with one optimized loop over all equations. It is only allowed to |
---|
344 | !> be called for the Wicker and Skamarock or Piascek-Williams advection scheme. |
---|
345 | !> |
---|
346 | !> Here the calls of most subroutines are embedded in two DO loops over i and j, |
---|
347 | !> so communication between CPUs is not allowed (does not make sense) within |
---|
348 | !> these loops. |
---|
349 | !> |
---|
350 | !> (Optimized to avoid cache missings, i.e. for Power4/5-architectures.) |
---|
351 | !------------------------------------------------------------------------------! |
---|
352 | |
---|
353 | SUBROUTINE prognostic_equations_cache |
---|
354 | |
---|
355 | |
---|
356 | IMPLICIT NONE |
---|
357 | |
---|
358 | INTEGER(iwp) :: i !< |
---|
359 | INTEGER(iwp) :: i_omp_start !< |
---|
360 | INTEGER(iwp) :: j !< |
---|
361 | INTEGER(iwp) :: k !< |
---|
362 | INTEGER(iwp) :: omp_get_thread_num !< |
---|
363 | INTEGER(iwp) :: tn = 0 !< |
---|
364 | |
---|
365 | LOGICAL :: loop_start !< |
---|
366 | |
---|
367 | |
---|
368 | ! |
---|
369 | !-- Time measurement can only be performed for the whole set of equations |
---|
370 | CALL cpu_log( log_point(32), 'all progn.equations', 'start' ) |
---|
371 | |
---|
372 | ! |
---|
373 | !-- Loop over all prognostic equations |
---|
374 | !$OMP PARALLEL private (i,i_omp_start,j,k,loop_start,tn) |
---|
375 | |
---|
376 | !$ tn = omp_get_thread_num() |
---|
377 | loop_start = .TRUE. |
---|
378 | !$OMP DO |
---|
379 | DO i = nxl, nxr |
---|
380 | |
---|
381 | ! |
---|
382 | !-- Store the first loop index. It differs for each thread and is required |
---|
383 | !-- later in advec_ws |
---|
384 | IF ( loop_start ) THEN |
---|
385 | loop_start = .FALSE. |
---|
386 | i_omp_start = i |
---|
387 | ENDIF |
---|
388 | |
---|
389 | DO j = nys, nyn |
---|
390 | ! |
---|
391 | !-- If required, calculate cloud microphysics |
---|
392 | IF ( cloud_physics .AND. .NOT. microphysics_sat_adjust .AND. & |
---|
393 | ( intermediate_timestep_count == 1 .OR. & |
---|
394 | call_microphysics_at_all_substeps ) & |
---|
395 | ) THEN |
---|
396 | CALL microphysics_control( i, j ) |
---|
397 | ENDIF |
---|
398 | ! |
---|
399 | !-- Tendency terms for u-velocity component |
---|
400 | IF ( .NOT. outflow_l .OR. i > nxl ) THEN |
---|
401 | |
---|
402 | tend(:,j,i) = 0.0_wp |
---|
403 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
404 | IF ( ws_scheme_mom ) THEN |
---|
405 | CALL advec_u_ws( i, j, i_omp_start, tn ) |
---|
406 | ELSE |
---|
407 | CALL advec_u_pw( i, j ) |
---|
408 | ENDIF |
---|
409 | ELSE |
---|
410 | CALL advec_u_up( i, j ) |
---|
411 | ENDIF |
---|
412 | CALL diffusion_u( i, j ) |
---|
413 | CALL coriolis( i, j, 1 ) |
---|
414 | IF ( sloping_surface .AND. .NOT. neutral ) THEN |
---|
415 | CALL buoyancy( i, j, pt, 1 ) |
---|
416 | ENDIF |
---|
417 | |
---|
418 | ! |
---|
419 | !-- Drag by plant canopy |
---|
420 | IF ( plant_canopy ) CALL pcm_tendency( i, j, 1 ) |
---|
421 | |
---|
422 | ! |
---|
423 | !-- External pressure gradient |
---|
424 | IF ( dp_external ) THEN |
---|
425 | DO k = dp_level_ind_b+1, nzt |
---|
426 | tend(k,j,i) = tend(k,j,i) - dpdxy(1) * dp_smooth_factor(k) |
---|
427 | ENDDO |
---|
428 | ENDIF |
---|
429 | |
---|
430 | ! |
---|
431 | !-- Nudging |
---|
432 | IF ( nudging ) CALL nudge( i, j, simulated_time, 'u' ) |
---|
433 | |
---|
434 | ! |
---|
435 | !-- Forces by wind turbines |
---|
436 | IF ( wind_turbine ) CALL wtm_tendencies( i, j, 1 ) |
---|
437 | |
---|
438 | CALL user_actions( i, j, 'u-tendency' ) |
---|
439 | ! |
---|
440 | !-- Prognostic equation for u-velocity component |
---|
441 | DO k = nzb_u_inner(j,i)+1, nzt |
---|
442 | u_p(k,j,i) = u(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
443 | tsc(3) * tu_m(k,j,i) ) & |
---|
444 | - tsc(5) * rdf(k) * ( u(k,j,i) - u_init(k) ) |
---|
445 | ENDDO |
---|
446 | |
---|
447 | ! |
---|
448 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
449 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
450 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
451 | DO k = nzb_u_inner(j,i)+1, nzt |
---|
452 | tu_m(k,j,i) = tend(k,j,i) |
---|
453 | ENDDO |
---|
454 | ELSEIF ( intermediate_timestep_count < & |
---|
455 | intermediate_timestep_count_max ) THEN |
---|
456 | DO k = nzb_u_inner(j,i)+1, nzt |
---|
457 | tu_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tu_m(k,j,i) |
---|
458 | ENDDO |
---|
459 | ENDIF |
---|
460 | ENDIF |
---|
461 | |
---|
462 | ENDIF |
---|
463 | |
---|
464 | ! |
---|
465 | !-- Tendency terms for v-velocity component |
---|
466 | IF ( .NOT. outflow_s .OR. j > nys ) THEN |
---|
467 | |
---|
468 | tend(:,j,i) = 0.0_wp |
---|
469 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
470 | IF ( ws_scheme_mom ) THEN |
---|
471 | CALL advec_v_ws( i, j, i_omp_start, tn ) |
---|
472 | ELSE |
---|
473 | CALL advec_v_pw( i, j ) |
---|
474 | ENDIF |
---|
475 | ELSE |
---|
476 | CALL advec_v_up( i, j ) |
---|
477 | ENDIF |
---|
478 | CALL diffusion_v( i, j ) |
---|
479 | CALL coriolis( i, j, 2 ) |
---|
480 | |
---|
481 | ! |
---|
482 | !-- Drag by plant canopy |
---|
483 | IF ( plant_canopy ) CALL pcm_tendency( i, j, 2 ) |
---|
484 | |
---|
485 | ! |
---|
486 | !-- External pressure gradient |
---|
487 | IF ( dp_external ) THEN |
---|
488 | DO k = dp_level_ind_b+1, nzt |
---|
489 | tend(k,j,i) = tend(k,j,i) - dpdxy(2) * dp_smooth_factor(k) |
---|
490 | ENDDO |
---|
491 | ENDIF |
---|
492 | |
---|
493 | ! |
---|
494 | !-- Nudging |
---|
495 | IF ( nudging ) CALL nudge( i, j, simulated_time, 'v' ) |
---|
496 | |
---|
497 | ! |
---|
498 | !-- Forces by wind turbines |
---|
499 | IF ( wind_turbine ) CALL wtm_tendencies( i, j, 2 ) |
---|
500 | |
---|
501 | CALL user_actions( i, j, 'v-tendency' ) |
---|
502 | ! |
---|
503 | !-- Prognostic equation for v-velocity component |
---|
504 | DO k = nzb_v_inner(j,i)+1, nzt |
---|
505 | v_p(k,j,i) = v(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
506 | tsc(3) * tv_m(k,j,i) ) & |
---|
507 | - tsc(5) * rdf(k) * ( v(k,j,i) - v_init(k) ) |
---|
508 | ENDDO |
---|
509 | |
---|
510 | ! |
---|
511 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
512 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
513 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
514 | DO k = nzb_v_inner(j,i)+1, nzt |
---|
515 | tv_m(k,j,i) = tend(k,j,i) |
---|
516 | ENDDO |
---|
517 | ELSEIF ( intermediate_timestep_count < & |
---|
518 | intermediate_timestep_count_max ) THEN |
---|
519 | DO k = nzb_v_inner(j,i)+1, nzt |
---|
520 | tv_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tv_m(k,j,i) |
---|
521 | ENDDO |
---|
522 | ENDIF |
---|
523 | ENDIF |
---|
524 | |
---|
525 | ENDIF |
---|
526 | |
---|
527 | ! |
---|
528 | !-- Tendency terms for w-velocity component |
---|
529 | tend(:,j,i) = 0.0_wp |
---|
530 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
531 | IF ( ws_scheme_mom ) THEN |
---|
532 | CALL advec_w_ws( i, j, i_omp_start, tn ) |
---|
533 | ELSE |
---|
534 | CALL advec_w_pw( i, j ) |
---|
535 | END IF |
---|
536 | ELSE |
---|
537 | CALL advec_w_up( i, j ) |
---|
538 | ENDIF |
---|
539 | CALL diffusion_w( i, j ) |
---|
540 | CALL coriolis( i, j, 3 ) |
---|
541 | |
---|
542 | IF ( .NOT. neutral ) THEN |
---|
543 | IF ( ocean ) THEN |
---|
544 | CALL buoyancy( i, j, rho, 3 ) |
---|
545 | ELSE |
---|
546 | IF ( .NOT. humidity ) THEN |
---|
547 | CALL buoyancy( i, j, pt, 3 ) |
---|
548 | ELSE |
---|
549 | CALL buoyancy( i, j, vpt, 3 ) |
---|
550 | ENDIF |
---|
551 | ENDIF |
---|
552 | ENDIF |
---|
553 | |
---|
554 | ! |
---|
555 | !-- Drag by plant canopy |
---|
556 | IF ( plant_canopy ) CALL pcm_tendency( i, j, 3 ) |
---|
557 | |
---|
558 | ! |
---|
559 | !-- Forces by wind turbines |
---|
560 | IF ( wind_turbine ) CALL wtm_tendencies( i, j, 3 ) |
---|
561 | |
---|
562 | CALL user_actions( i, j, 'w-tendency' ) |
---|
563 | |
---|
564 | ! |
---|
565 | !-- Prognostic equation for w-velocity component |
---|
566 | DO k = nzb_w_inner(j,i)+1, nzt-1 |
---|
567 | w_p(k,j,i) = w(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
568 | tsc(3) * tw_m(k,j,i) ) & |
---|
569 | - tsc(5) * rdf(k) * w(k,j,i) |
---|
570 | ENDDO |
---|
571 | |
---|
572 | ! |
---|
573 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
574 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
575 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
576 | DO k = nzb_w_inner(j,i)+1, nzt-1 |
---|
577 | tw_m(k,j,i) = tend(k,j,i) |
---|
578 | ENDDO |
---|
579 | ELSEIF ( intermediate_timestep_count < & |
---|
580 | intermediate_timestep_count_max ) THEN |
---|
581 | DO k = nzb_w_inner(j,i)+1, nzt-1 |
---|
582 | tw_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tw_m(k,j,i) |
---|
583 | ENDDO |
---|
584 | ENDIF |
---|
585 | ENDIF |
---|
586 | |
---|
587 | ! |
---|
588 | !-- If required, compute prognostic equation for potential temperature |
---|
589 | IF ( .NOT. neutral ) THEN |
---|
590 | ! |
---|
591 | !-- Tendency terms for potential temperature |
---|
592 | tend(:,j,i) = 0.0_wp |
---|
593 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
594 | IF ( ws_scheme_sca ) THEN |
---|
595 | CALL advec_s_ws( i, j, pt, 'pt', flux_s_pt, diss_s_pt, & |
---|
596 | flux_l_pt, diss_l_pt, i_omp_start, tn ) |
---|
597 | ELSE |
---|
598 | CALL advec_s_pw( i, j, pt ) |
---|
599 | ENDIF |
---|
600 | ELSE |
---|
601 | CALL advec_s_up( i, j, pt ) |
---|
602 | ENDIF |
---|
603 | CALL diffusion_s( i, j, pt, shf, tswst, wall_heatflux ) |
---|
604 | |
---|
605 | ! |
---|
606 | !-- If required compute heating/cooling due to long wave radiation |
---|
607 | !-- processes |
---|
608 | IF ( cloud_top_radiation ) THEN |
---|
609 | CALL calc_radiation( i, j ) |
---|
610 | ENDIF |
---|
611 | |
---|
612 | ! |
---|
613 | !-- Consideration of heat sources within the plant canopy |
---|
614 | IF ( plant_canopy .AND. cthf /= 0.0_wp ) THEN |
---|
615 | CALL pcm_tendency( i, j, 4 ) |
---|
616 | ENDIF |
---|
617 | |
---|
618 | ! |
---|
619 | !-- Large scale advection |
---|
620 | IF ( large_scale_forcing ) THEN |
---|
621 | CALL ls_advec( i, j, simulated_time, 'pt' ) |
---|
622 | ENDIF |
---|
623 | |
---|
624 | ! |
---|
625 | !-- Nudging |
---|
626 | IF ( nudging ) CALL nudge( i, j, simulated_time, 'pt' ) |
---|
627 | |
---|
628 | ! |
---|
629 | !-- If required, compute effect of large-scale subsidence/ascent |
---|
630 | IF ( large_scale_subsidence .AND. & |
---|
631 | .NOT. use_subsidence_tendencies ) THEN |
---|
632 | CALL subsidence( i, j, tend, pt, pt_init, 2 ) |
---|
633 | ENDIF |
---|
634 | |
---|
635 | ! |
---|
636 | !-- If required, add tendency due to radiative heating/cooling |
---|
637 | IF ( radiation .AND. radiation_scheme == 'rrtmg' .AND. & |
---|
638 | simulated_time > skip_time_do_radiation ) THEN |
---|
639 | CALL radiation_tendency ( i, j, tend ) |
---|
640 | ENDIF |
---|
641 | |
---|
642 | |
---|
643 | CALL user_actions( i, j, 'pt-tendency' ) |
---|
644 | |
---|
645 | ! |
---|
646 | !-- Prognostic equation for potential temperature |
---|
647 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
648 | pt_p(k,j,i) = pt(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
649 | tsc(3) * tpt_m(k,j,i) ) & |
---|
650 | - tsc(5) * ( pt(k,j,i) - pt_init(k) ) *& |
---|
651 | ( rdf_sc(k) + ptdf_x(i) + ptdf_y(j) ) |
---|
652 | ENDDO |
---|
653 | |
---|
654 | ! |
---|
655 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
656 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
657 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
658 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
659 | tpt_m(k,j,i) = tend(k,j,i) |
---|
660 | ENDDO |
---|
661 | ELSEIF ( intermediate_timestep_count < & |
---|
662 | intermediate_timestep_count_max ) THEN |
---|
663 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
664 | tpt_m(k,j,i) = -9.5625_wp * tend(k,j,i) + & |
---|
665 | 5.3125_wp * tpt_m(k,j,i) |
---|
666 | ENDDO |
---|
667 | ENDIF |
---|
668 | ENDIF |
---|
669 | |
---|
670 | ENDIF |
---|
671 | |
---|
672 | ! |
---|
673 | !-- If required, compute prognostic equation for salinity |
---|
674 | IF ( ocean ) THEN |
---|
675 | |
---|
676 | ! |
---|
677 | !-- Tendency-terms for salinity |
---|
678 | tend(:,j,i) = 0.0_wp |
---|
679 | IF ( timestep_scheme(1:5) == 'runge' ) & |
---|
680 | THEN |
---|
681 | IF ( ws_scheme_sca ) THEN |
---|
682 | CALL advec_s_ws( i, j, sa, 'sa', flux_s_sa, & |
---|
683 | diss_s_sa, flux_l_sa, diss_l_sa, i_omp_start, tn ) |
---|
684 | ELSE |
---|
685 | CALL advec_s_pw( i, j, sa ) |
---|
686 | ENDIF |
---|
687 | ELSE |
---|
688 | CALL advec_s_up( i, j, sa ) |
---|
689 | ENDIF |
---|
690 | CALL diffusion_s( i, j, sa, saswsb, saswst, wall_salinityflux ) |
---|
691 | |
---|
692 | CALL user_actions( i, j, 'sa-tendency' ) |
---|
693 | |
---|
694 | ! |
---|
695 | !-- Prognostic equation for salinity |
---|
696 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
697 | sa_p(k,j,i) = sa(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
698 | tsc(3) * tsa_m(k,j,i) ) & |
---|
699 | - tsc(5) * rdf_sc(k) * & |
---|
700 | ( sa(k,j,i) - sa_init(k) ) |
---|
701 | IF ( sa_p(k,j,i) < 0.0_wp ) sa_p(k,j,i) = 0.1_wp * sa(k,j,i) |
---|
702 | ENDDO |
---|
703 | |
---|
704 | ! |
---|
705 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
706 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
707 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
708 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
709 | tsa_m(k,j,i) = tend(k,j,i) |
---|
710 | ENDDO |
---|
711 | ELSEIF ( intermediate_timestep_count < & |
---|
712 | intermediate_timestep_count_max ) THEN |
---|
713 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
714 | tsa_m(k,j,i) = -9.5625_wp * tend(k,j,i) + & |
---|
715 | 5.3125_wp * tsa_m(k,j,i) |
---|
716 | ENDDO |
---|
717 | ENDIF |
---|
718 | ENDIF |
---|
719 | |
---|
720 | ! |
---|
721 | !-- Calculate density by the equation of state for seawater |
---|
722 | CALL eqn_state_seawater( i, j ) |
---|
723 | |
---|
724 | ENDIF |
---|
725 | |
---|
726 | ! |
---|
727 | !-- If required, compute prognostic equation for total water content / |
---|
728 | !-- scalar |
---|
729 | IF ( humidity .OR. passive_scalar ) THEN |
---|
730 | |
---|
731 | ! |
---|
732 | !-- Tendency-terms for total water content / scalar |
---|
733 | tend(:,j,i) = 0.0_wp |
---|
734 | IF ( timestep_scheme(1:5) == 'runge' ) & |
---|
735 | THEN |
---|
736 | IF ( ws_scheme_sca ) THEN |
---|
737 | CALL advec_s_ws( i, j, q, 'q', flux_s_q, & |
---|
738 | diss_s_q, flux_l_q, diss_l_q, i_omp_start, tn ) |
---|
739 | ELSE |
---|
740 | CALL advec_s_pw( i, j, q ) |
---|
741 | ENDIF |
---|
742 | ELSE |
---|
743 | CALL advec_s_up( i, j, q ) |
---|
744 | ENDIF |
---|
745 | CALL diffusion_s( i, j, q, qsws, qswst, wall_qflux ) |
---|
746 | |
---|
747 | ! |
---|
748 | !-- Sink or source of scalar concentration due to canopy elements |
---|
749 | IF ( plant_canopy ) CALL pcm_tendency( i, j, 5 ) |
---|
750 | |
---|
751 | ! |
---|
752 | !-- Large scale advection |
---|
753 | IF ( large_scale_forcing ) THEN |
---|
754 | CALL ls_advec( i, j, simulated_time, 'q' ) |
---|
755 | ENDIF |
---|
756 | |
---|
757 | ! |
---|
758 | !-- Nudging |
---|
759 | IF ( nudging ) CALL nudge( i, j, simulated_time, 'q' ) |
---|
760 | |
---|
761 | ! |
---|
762 | !-- If required compute influence of large-scale subsidence/ascent |
---|
763 | IF ( large_scale_subsidence .AND. & |
---|
764 | .NOT. use_subsidence_tendencies ) THEN |
---|
765 | CALL subsidence( i, j, tend, q, q_init, 3 ) |
---|
766 | ENDIF |
---|
767 | |
---|
768 | CALL user_actions( i, j, 'q-tendency' ) |
---|
769 | |
---|
770 | ! |
---|
771 | !-- Prognostic equation for total water content / scalar |
---|
772 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
773 | q_p(k,j,i) = q(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
774 | tsc(3) * tq_m(k,j,i) ) & |
---|
775 | - tsc(5) * rdf_sc(k) * & |
---|
776 | ( q(k,j,i) - q_init(k) ) |
---|
777 | IF ( q_p(k,j,i) < 0.0_wp ) q_p(k,j,i) = 0.1_wp * q(k,j,i) |
---|
778 | ENDDO |
---|
779 | |
---|
780 | ! |
---|
781 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
782 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
783 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
784 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
785 | tq_m(k,j,i) = tend(k,j,i) |
---|
786 | ENDDO |
---|
787 | ELSEIF ( intermediate_timestep_count < & |
---|
788 | intermediate_timestep_count_max ) THEN |
---|
789 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
790 | tq_m(k,j,i) = -9.5625_wp * tend(k,j,i) + & |
---|
791 | 5.3125_wp * tq_m(k,j,i) |
---|
792 | ENDDO |
---|
793 | ENDIF |
---|
794 | ENDIF |
---|
795 | |
---|
796 | ! |
---|
797 | !-- If required, calculate prognostic equations for rain water content |
---|
798 | !-- and rain drop concentration |
---|
799 | IF ( cloud_physics .AND. microphysics_seifert ) THEN |
---|
800 | ! |
---|
801 | !-- Calculate prognostic equation for rain water content |
---|
802 | tend(:,j,i) = 0.0_wp |
---|
803 | IF ( timestep_scheme(1:5) == 'runge' ) & |
---|
804 | THEN |
---|
805 | IF ( ws_scheme_sca ) THEN |
---|
806 | CALL advec_s_ws( i, j, qr, 'qr', flux_s_qr, & |
---|
807 | diss_s_qr, flux_l_qr, diss_l_qr, & |
---|
808 | i_omp_start, tn ) |
---|
809 | ELSE |
---|
810 | CALL advec_s_pw( i, j, qr ) |
---|
811 | ENDIF |
---|
812 | ELSE |
---|
813 | CALL advec_s_up( i, j, qr ) |
---|
814 | ENDIF |
---|
815 | CALL diffusion_s( i, j, qr, qrsws, qrswst, wall_qrflux ) |
---|
816 | |
---|
817 | ! |
---|
818 | !-- Prognostic equation for rain water content |
---|
819 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
820 | qr_p(k,j,i) = qr(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
821 | tsc(3) * tqr_m(k,j,i) ) & |
---|
822 | - tsc(5) * rdf_sc(k) * qr(k,j,i) |
---|
823 | IF ( qr_p(k,j,i) < 0.0_wp ) qr_p(k,j,i) = 0.0_wp |
---|
824 | ENDDO |
---|
825 | ! |
---|
826 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
827 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
828 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
829 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
830 | tqr_m(k,j,i) = tend(k,j,i) |
---|
831 | ENDDO |
---|
832 | ELSEIF ( intermediate_timestep_count < & |
---|
833 | intermediate_timestep_count_max ) THEN |
---|
834 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
835 | tqr_m(k,j,i) = -9.5625_wp * tend(k,j,i) + & |
---|
836 | 5.3125_wp * tqr_m(k,j,i) |
---|
837 | ENDDO |
---|
838 | ENDIF |
---|
839 | ENDIF |
---|
840 | |
---|
841 | ! |
---|
842 | !-- Calculate prognostic equation for rain drop concentration. |
---|
843 | tend(:,j,i) = 0.0_wp |
---|
844 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
845 | IF ( ws_scheme_sca ) THEN |
---|
846 | CALL advec_s_ws( i, j, nr, 'nr', flux_s_nr, & |
---|
847 | diss_s_nr, flux_l_nr, diss_l_nr, & |
---|
848 | i_omp_start, tn ) |
---|
849 | ELSE |
---|
850 | CALL advec_s_pw( i, j, nr ) |
---|
851 | ENDIF |
---|
852 | ELSE |
---|
853 | CALL advec_s_up( i, j, nr ) |
---|
854 | ENDIF |
---|
855 | CALL diffusion_s( i, j, nr, nrsws, nrswst, wall_nrflux ) |
---|
856 | |
---|
857 | ! |
---|
858 | !-- Prognostic equation for rain drop concentration |
---|
859 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
860 | nr_p(k,j,i) = nr(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
861 | tsc(3) * tnr_m(k,j,i) ) & |
---|
862 | - tsc(5) * rdf_sc(k) * nr(k,j,i) |
---|
863 | IF ( nr_p(k,j,i) < 0.0_wp ) nr_p(k,j,i) = 0.0_wp |
---|
864 | ENDDO |
---|
865 | ! |
---|
866 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
867 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
868 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
869 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
870 | tnr_m(k,j,i) = tend(k,j,i) |
---|
871 | ENDDO |
---|
872 | ELSEIF ( intermediate_timestep_count < & |
---|
873 | intermediate_timestep_count_max ) THEN |
---|
874 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
875 | tnr_m(k,j,i) = -9.5625_wp * tend(k,j,i) + & |
---|
876 | 5.3125_wp * tnr_m(k,j,i) |
---|
877 | ENDDO |
---|
878 | ENDIF |
---|
879 | ENDIF |
---|
880 | |
---|
881 | ENDIF |
---|
882 | |
---|
883 | ENDIF |
---|
884 | |
---|
885 | ! |
---|
886 | !-- If required, compute prognostic equation for turbulent kinetic |
---|
887 | !-- energy (TKE) |
---|
888 | IF ( .NOT. constant_diffusion ) THEN |
---|
889 | |
---|
890 | ! |
---|
891 | !-- Tendency-terms for TKE |
---|
892 | tend(:,j,i) = 0.0_wp |
---|
893 | IF ( timestep_scheme(1:5) == 'runge' & |
---|
894 | .AND. .NOT. use_upstream_for_tke ) THEN |
---|
895 | IF ( ws_scheme_sca ) THEN |
---|
896 | CALL advec_s_ws( i, j, e, 'e', flux_s_e, diss_s_e, & |
---|
897 | flux_l_e, diss_l_e , i_omp_start, tn ) |
---|
898 | ELSE |
---|
899 | CALL advec_s_pw( i, j, e ) |
---|
900 | ENDIF |
---|
901 | ELSE |
---|
902 | CALL advec_s_up( i, j, e ) |
---|
903 | ENDIF |
---|
904 | IF ( .NOT. humidity ) THEN |
---|
905 | IF ( ocean ) THEN |
---|
906 | CALL diffusion_e( i, j, prho, prho_reference ) |
---|
907 | ELSE |
---|
908 | CALL diffusion_e( i, j, pt, pt_reference ) |
---|
909 | ENDIF |
---|
910 | ELSE |
---|
911 | CALL diffusion_e( i, j, vpt, pt_reference ) |
---|
912 | ENDIF |
---|
913 | CALL production_e( i, j ) |
---|
914 | |
---|
915 | ! |
---|
916 | !-- Additional sink term for flows through plant canopies |
---|
917 | IF ( plant_canopy ) CALL pcm_tendency( i, j, 6 ) |
---|
918 | |
---|
919 | CALL user_actions( i, j, 'e-tendency' ) |
---|
920 | |
---|
921 | ! |
---|
922 | !-- Prognostic equation for TKE. |
---|
923 | !-- Eliminate negative TKE values, which can occur due to numerical |
---|
924 | !-- reasons in the course of the integration. In such cases the old |
---|
925 | !-- TKE value is reduced by 90%. |
---|
926 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
927 | e_p(k,j,i) = e(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
928 | tsc(3) * te_m(k,j,i) ) |
---|
929 | IF ( e_p(k,j,i) < 0.0_wp ) e_p(k,j,i) = 0.1_wp * e(k,j,i) |
---|
930 | ENDDO |
---|
931 | |
---|
932 | ! |
---|
933 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
934 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
935 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
936 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
937 | te_m(k,j,i) = tend(k,j,i) |
---|
938 | ENDDO |
---|
939 | ELSEIF ( intermediate_timestep_count < & |
---|
940 | intermediate_timestep_count_max ) THEN |
---|
941 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
942 | te_m(k,j,i) = -9.5625_wp * tend(k,j,i) + & |
---|
943 | 5.3125_wp * te_m(k,j,i) |
---|
944 | ENDDO |
---|
945 | ENDIF |
---|
946 | ENDIF |
---|
947 | |
---|
948 | ENDIF ! TKE equation |
---|
949 | |
---|
950 | ENDDO |
---|
951 | ENDDO |
---|
952 | !$OMP END PARALLEL |
---|
953 | |
---|
954 | CALL cpu_log( log_point(32), 'all progn.equations', 'stop' ) |
---|
955 | |
---|
956 | |
---|
957 | END SUBROUTINE prognostic_equations_cache |
---|
958 | |
---|
959 | |
---|
960 | !------------------------------------------------------------------------------! |
---|
961 | ! Description: |
---|
962 | ! ------------ |
---|
963 | !> Version for vector machines |
---|
964 | !------------------------------------------------------------------------------! |
---|
965 | |
---|
966 | SUBROUTINE prognostic_equations_vector |
---|
967 | |
---|
968 | |
---|
969 | IMPLICIT NONE |
---|
970 | |
---|
971 | INTEGER(iwp) :: i !< |
---|
972 | INTEGER(iwp) :: j !< |
---|
973 | INTEGER(iwp) :: k !< |
---|
974 | |
---|
975 | REAL(wp) :: sbt !< |
---|
976 | |
---|
977 | |
---|
978 | ! |
---|
979 | !-- If required, calculate cloud microphysical impacts |
---|
980 | IF ( cloud_physics .AND. .NOT. microphysics_sat_adjust .AND. & |
---|
981 | ( intermediate_timestep_count == 1 .OR. & |
---|
982 | call_microphysics_at_all_substeps ) & |
---|
983 | ) THEN |
---|
984 | CALL cpu_log( log_point(51), 'microphysics', 'start' ) |
---|
985 | CALL microphysics_control |
---|
986 | CALL cpu_log( log_point(51), 'microphysics', 'stop' ) |
---|
987 | ENDIF |
---|
988 | |
---|
989 | ! |
---|
990 | !-- u-velocity component |
---|
991 | CALL cpu_log( log_point(5), 'u-equation', 'start' ) |
---|
992 | |
---|
993 | tend = 0.0_wp |
---|
994 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
995 | IF ( ws_scheme_mom ) THEN |
---|
996 | CALL advec_u_ws |
---|
997 | ELSE |
---|
998 | CALL advec_u_pw |
---|
999 | ENDIF |
---|
1000 | ELSE |
---|
1001 | CALL advec_u_up |
---|
1002 | ENDIF |
---|
1003 | CALL diffusion_u |
---|
1004 | CALL coriolis( 1 ) |
---|
1005 | IF ( sloping_surface .AND. .NOT. neutral ) THEN |
---|
1006 | CALL buoyancy( pt, 1 ) |
---|
1007 | ENDIF |
---|
1008 | |
---|
1009 | ! |
---|
1010 | !-- Drag by plant canopy |
---|
1011 | IF ( plant_canopy ) CALL pcm_tendency( 1 ) |
---|
1012 | |
---|
1013 | ! |
---|
1014 | !-- External pressure gradient |
---|
1015 | IF ( dp_external ) THEN |
---|
1016 | DO i = nxlu, nxr |
---|
1017 | DO j = nys, nyn |
---|
1018 | DO k = dp_level_ind_b+1, nzt |
---|
1019 | tend(k,j,i) = tend(k,j,i) - dpdxy(1) * dp_smooth_factor(k) |
---|
1020 | ENDDO |
---|
1021 | ENDDO |
---|
1022 | ENDDO |
---|
1023 | ENDIF |
---|
1024 | |
---|
1025 | ! |
---|
1026 | !-- Nudging |
---|
1027 | IF ( nudging ) CALL nudge( simulated_time, 'u' ) |
---|
1028 | |
---|
1029 | ! |
---|
1030 | !-- Forces by wind turbines |
---|
1031 | IF ( wind_turbine ) CALL wtm_tendencies( 1 ) |
---|
1032 | |
---|
1033 | CALL user_actions( 'u-tendency' ) |
---|
1034 | |
---|
1035 | ! |
---|
1036 | !-- Prognostic equation for u-velocity component |
---|
1037 | DO i = nxlu, nxr |
---|
1038 | DO j = nys, nyn |
---|
1039 | DO k = nzb_u_inner(j,i)+1, nzt |
---|
1040 | u_p(k,j,i) = u(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
1041 | tsc(3) * tu_m(k,j,i) ) & |
---|
1042 | - tsc(5) * rdf(k) * ( u(k,j,i) - u_init(k) ) |
---|
1043 | ENDDO |
---|
1044 | ENDDO |
---|
1045 | ENDDO |
---|
1046 | |
---|
1047 | ! |
---|
1048 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
1049 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1050 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
1051 | DO i = nxlu, nxr |
---|
1052 | DO j = nys, nyn |
---|
1053 | DO k = nzb_u_inner(j,i)+1, nzt |
---|
1054 | tu_m(k,j,i) = tend(k,j,i) |
---|
1055 | ENDDO |
---|
1056 | ENDDO |
---|
1057 | ENDDO |
---|
1058 | ELSEIF ( intermediate_timestep_count < & |
---|
1059 | intermediate_timestep_count_max ) THEN |
---|
1060 | DO i = nxlu, nxr |
---|
1061 | DO j = nys, nyn |
---|
1062 | DO k = nzb_u_inner(j,i)+1, nzt |
---|
1063 | tu_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tu_m(k,j,i) |
---|
1064 | ENDDO |
---|
1065 | ENDDO |
---|
1066 | ENDDO |
---|
1067 | ENDIF |
---|
1068 | ENDIF |
---|
1069 | |
---|
1070 | CALL cpu_log( log_point(5), 'u-equation', 'stop' ) |
---|
1071 | |
---|
1072 | ! |
---|
1073 | !-- v-velocity component |
---|
1074 | CALL cpu_log( log_point(6), 'v-equation', 'start' ) |
---|
1075 | |
---|
1076 | tend = 0.0_wp |
---|
1077 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1078 | IF ( ws_scheme_mom ) THEN |
---|
1079 | CALL advec_v_ws |
---|
1080 | ELSE |
---|
1081 | CALL advec_v_pw |
---|
1082 | END IF |
---|
1083 | ELSE |
---|
1084 | CALL advec_v_up |
---|
1085 | ENDIF |
---|
1086 | CALL diffusion_v |
---|
1087 | CALL coriolis( 2 ) |
---|
1088 | |
---|
1089 | ! |
---|
1090 | !-- Drag by plant canopy |
---|
1091 | IF ( plant_canopy ) CALL pcm_tendency( 2 ) |
---|
1092 | |
---|
1093 | ! |
---|
1094 | !-- External pressure gradient |
---|
1095 | IF ( dp_external ) THEN |
---|
1096 | DO i = nxl, nxr |
---|
1097 | DO j = nysv, nyn |
---|
1098 | DO k = dp_level_ind_b+1, nzt |
---|
1099 | tend(k,j,i) = tend(k,j,i) - dpdxy(2) * dp_smooth_factor(k) |
---|
1100 | ENDDO |
---|
1101 | ENDDO |
---|
1102 | ENDDO |
---|
1103 | ENDIF |
---|
1104 | |
---|
1105 | ! |
---|
1106 | !-- Nudging |
---|
1107 | IF ( nudging ) CALL nudge( simulated_time, 'v' ) |
---|
1108 | |
---|
1109 | ! |
---|
1110 | !-- Forces by wind turbines |
---|
1111 | IF ( wind_turbine ) CALL wtm_tendencies( 2 ) |
---|
1112 | |
---|
1113 | CALL user_actions( 'v-tendency' ) |
---|
1114 | |
---|
1115 | ! |
---|
1116 | !-- Prognostic equation for v-velocity component |
---|
1117 | DO i = nxl, nxr |
---|
1118 | DO j = nysv, nyn |
---|
1119 | DO k = nzb_v_inner(j,i)+1, nzt |
---|
1120 | v_p(k,j,i) = v(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
1121 | tsc(3) * tv_m(k,j,i) ) & |
---|
1122 | - tsc(5) * rdf(k) * ( v(k,j,i) - v_init(k) ) |
---|
1123 | ENDDO |
---|
1124 | ENDDO |
---|
1125 | ENDDO |
---|
1126 | |
---|
1127 | ! |
---|
1128 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
1129 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1130 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
1131 | DO i = nxl, nxr |
---|
1132 | DO j = nysv, nyn |
---|
1133 | DO k = nzb_v_inner(j,i)+1, nzt |
---|
1134 | tv_m(k,j,i) = tend(k,j,i) |
---|
1135 | ENDDO |
---|
1136 | ENDDO |
---|
1137 | ENDDO |
---|
1138 | ELSEIF ( intermediate_timestep_count < & |
---|
1139 | intermediate_timestep_count_max ) THEN |
---|
1140 | DO i = nxl, nxr |
---|
1141 | DO j = nysv, nyn |
---|
1142 | DO k = nzb_v_inner(j,i)+1, nzt |
---|
1143 | tv_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tv_m(k,j,i) |
---|
1144 | ENDDO |
---|
1145 | ENDDO |
---|
1146 | ENDDO |
---|
1147 | ENDIF |
---|
1148 | ENDIF |
---|
1149 | |
---|
1150 | CALL cpu_log( log_point(6), 'v-equation', 'stop' ) |
---|
1151 | |
---|
1152 | ! |
---|
1153 | !-- w-velocity component |
---|
1154 | CALL cpu_log( log_point(7), 'w-equation', 'start' ) |
---|
1155 | |
---|
1156 | tend = 0.0_wp |
---|
1157 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1158 | IF ( ws_scheme_mom ) THEN |
---|
1159 | CALL advec_w_ws |
---|
1160 | ELSE |
---|
1161 | CALL advec_w_pw |
---|
1162 | ENDIF |
---|
1163 | ELSE |
---|
1164 | CALL advec_w_up |
---|
1165 | ENDIF |
---|
1166 | CALL diffusion_w |
---|
1167 | CALL coriolis( 3 ) |
---|
1168 | |
---|
1169 | IF ( .NOT. neutral ) THEN |
---|
1170 | IF ( ocean ) THEN |
---|
1171 | CALL buoyancy( rho, 3 ) |
---|
1172 | ELSE |
---|
1173 | IF ( .NOT. humidity ) THEN |
---|
1174 | CALL buoyancy( pt, 3 ) |
---|
1175 | ELSE |
---|
1176 | CALL buoyancy( vpt, 3 ) |
---|
1177 | ENDIF |
---|
1178 | ENDIF |
---|
1179 | ENDIF |
---|
1180 | |
---|
1181 | ! |
---|
1182 | !-- Drag by plant canopy |
---|
1183 | IF ( plant_canopy ) CALL pcm_tendency( 3 ) |
---|
1184 | |
---|
1185 | ! |
---|
1186 | !-- Forces by wind turbines |
---|
1187 | IF ( wind_turbine ) CALL wtm_tendencies( 3 ) |
---|
1188 | |
---|
1189 | CALL user_actions( 'w-tendency' ) |
---|
1190 | |
---|
1191 | ! |
---|
1192 | !-- Prognostic equation for w-velocity component |
---|
1193 | DO i = nxl, nxr |
---|
1194 | DO j = nys, nyn |
---|
1195 | DO k = nzb_w_inner(j,i)+1, nzt-1 |
---|
1196 | w_p(k,j,i) = w(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
1197 | tsc(3) * tw_m(k,j,i) ) & |
---|
1198 | - tsc(5) * rdf(k) * w(k,j,i) |
---|
1199 | ENDDO |
---|
1200 | ENDDO |
---|
1201 | ENDDO |
---|
1202 | |
---|
1203 | ! |
---|
1204 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
1205 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1206 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
1207 | DO i = nxl, nxr |
---|
1208 | DO j = nys, nyn |
---|
1209 | DO k = nzb_w_inner(j,i)+1, nzt-1 |
---|
1210 | tw_m(k,j,i) = tend(k,j,i) |
---|
1211 | ENDDO |
---|
1212 | ENDDO |
---|
1213 | ENDDO |
---|
1214 | ELSEIF ( intermediate_timestep_count < & |
---|
1215 | intermediate_timestep_count_max ) THEN |
---|
1216 | DO i = nxl, nxr |
---|
1217 | DO j = nys, nyn |
---|
1218 | DO k = nzb_w_inner(j,i)+1, nzt-1 |
---|
1219 | tw_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tw_m(k,j,i) |
---|
1220 | ENDDO |
---|
1221 | ENDDO |
---|
1222 | ENDDO |
---|
1223 | ENDIF |
---|
1224 | ENDIF |
---|
1225 | |
---|
1226 | CALL cpu_log( log_point(7), 'w-equation', 'stop' ) |
---|
1227 | |
---|
1228 | |
---|
1229 | ! |
---|
1230 | !-- If required, compute prognostic equation for potential temperature |
---|
1231 | IF ( .NOT. neutral ) THEN |
---|
1232 | |
---|
1233 | CALL cpu_log( log_point(13), 'pt-equation', 'start' ) |
---|
1234 | |
---|
1235 | ! |
---|
1236 | !-- pt-tendency terms with communication |
---|
1237 | sbt = tsc(2) |
---|
1238 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
1239 | |
---|
1240 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
1241 | ! |
---|
1242 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
1243 | sbt = 1.0_wp |
---|
1244 | ENDIF |
---|
1245 | tend = 0.0_wp |
---|
1246 | CALL advec_s_bc( pt, 'pt' ) |
---|
1247 | |
---|
1248 | ENDIF |
---|
1249 | |
---|
1250 | ! |
---|
1251 | !-- pt-tendency terms with no communication |
---|
1252 | IF ( scalar_advec /= 'bc-scheme' ) THEN |
---|
1253 | tend = 0.0_wp |
---|
1254 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1255 | IF ( ws_scheme_sca ) THEN |
---|
1256 | CALL advec_s_ws( pt, 'pt' ) |
---|
1257 | ELSE |
---|
1258 | CALL advec_s_pw( pt ) |
---|
1259 | ENDIF |
---|
1260 | ELSE |
---|
1261 | CALL advec_s_up( pt ) |
---|
1262 | ENDIF |
---|
1263 | ENDIF |
---|
1264 | |
---|
1265 | CALL diffusion_s( pt, shf, tswst, wall_heatflux ) |
---|
1266 | |
---|
1267 | ! |
---|
1268 | !-- If required compute heating/cooling due to long wave radiation processes |
---|
1269 | IF ( cloud_top_radiation ) THEN |
---|
1270 | CALL calc_radiation |
---|
1271 | ENDIF |
---|
1272 | |
---|
1273 | ! |
---|
1274 | !-- Consideration of heat sources within the plant canopy |
---|
1275 | IF ( plant_canopy .AND. ( cthf /= 0.0_wp ) ) THEN |
---|
1276 | CALL pcm_tendency( 4 ) |
---|
1277 | ENDIF |
---|
1278 | |
---|
1279 | ! |
---|
1280 | !-- Large scale advection |
---|
1281 | IF ( large_scale_forcing ) THEN |
---|
1282 | CALL ls_advec( simulated_time, 'pt' ) |
---|
1283 | ENDIF |
---|
1284 | |
---|
1285 | ! |
---|
1286 | !-- Nudging |
---|
1287 | IF ( nudging ) CALL nudge( simulated_time, 'pt' ) |
---|
1288 | |
---|
1289 | ! |
---|
1290 | !-- If required compute influence of large-scale subsidence/ascent |
---|
1291 | IF ( large_scale_subsidence .AND. & |
---|
1292 | .NOT. use_subsidence_tendencies ) THEN |
---|
1293 | CALL subsidence( tend, pt, pt_init, 2 ) |
---|
1294 | ENDIF |
---|
1295 | |
---|
1296 | ! |
---|
1297 | !-- If required, add tendency due to radiative heating/cooling |
---|
1298 | IF ( radiation .AND. radiation_scheme == 'rrtmg' .AND. & |
---|
1299 | simulated_time > skip_time_do_radiation ) THEN |
---|
1300 | CALL radiation_tendency ( tend ) |
---|
1301 | ENDIF |
---|
1302 | |
---|
1303 | CALL user_actions( 'pt-tendency' ) |
---|
1304 | |
---|
1305 | ! |
---|
1306 | !-- Prognostic equation for potential temperature |
---|
1307 | DO i = nxl, nxr |
---|
1308 | DO j = nys, nyn |
---|
1309 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1310 | pt_p(k,j,i) = pt(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
1311 | tsc(3) * tpt_m(k,j,i) ) & |
---|
1312 | - tsc(5) * ( pt(k,j,i) - pt_init(k) ) *& |
---|
1313 | ( rdf_sc(k) + ptdf_x(i) + ptdf_y(j) ) |
---|
1314 | ENDDO |
---|
1315 | ENDDO |
---|
1316 | ENDDO |
---|
1317 | |
---|
1318 | ! |
---|
1319 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
1320 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1321 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
1322 | DO i = nxl, nxr |
---|
1323 | DO j = nys, nyn |
---|
1324 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1325 | tpt_m(k,j,i) = tend(k,j,i) |
---|
1326 | ENDDO |
---|
1327 | ENDDO |
---|
1328 | ENDDO |
---|
1329 | ELSEIF ( intermediate_timestep_count < & |
---|
1330 | intermediate_timestep_count_max ) THEN |
---|
1331 | DO i = nxl, nxr |
---|
1332 | DO j = nys, nyn |
---|
1333 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1334 | tpt_m(k,j,i) = -9.5625_wp * tend(k,j,i) + & |
---|
1335 | 5.3125_wp * tpt_m(k,j,i) |
---|
1336 | ENDDO |
---|
1337 | ENDDO |
---|
1338 | ENDDO |
---|
1339 | ENDIF |
---|
1340 | ENDIF |
---|
1341 | |
---|
1342 | CALL cpu_log( log_point(13), 'pt-equation', 'stop' ) |
---|
1343 | |
---|
1344 | ENDIF |
---|
1345 | |
---|
1346 | ! |
---|
1347 | !-- If required, compute prognostic equation for salinity |
---|
1348 | IF ( ocean ) THEN |
---|
1349 | |
---|
1350 | CALL cpu_log( log_point(37), 'sa-equation', 'start' ) |
---|
1351 | |
---|
1352 | ! |
---|
1353 | !-- sa-tendency terms with communication |
---|
1354 | sbt = tsc(2) |
---|
1355 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
1356 | |
---|
1357 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
1358 | ! |
---|
1359 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
1360 | sbt = 1.0_wp |
---|
1361 | ENDIF |
---|
1362 | tend = 0.0_wp |
---|
1363 | CALL advec_s_bc( sa, 'sa' ) |
---|
1364 | |
---|
1365 | ENDIF |
---|
1366 | |
---|
1367 | ! |
---|
1368 | !-- sa-tendency terms with no communication |
---|
1369 | IF ( scalar_advec /= 'bc-scheme' ) THEN |
---|
1370 | tend = 0.0_wp |
---|
1371 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1372 | IF ( ws_scheme_sca ) THEN |
---|
1373 | CALL advec_s_ws( sa, 'sa' ) |
---|
1374 | ELSE |
---|
1375 | CALL advec_s_pw( sa ) |
---|
1376 | ENDIF |
---|
1377 | ELSE |
---|
1378 | CALL advec_s_up( sa ) |
---|
1379 | ENDIF |
---|
1380 | ENDIF |
---|
1381 | |
---|
1382 | CALL diffusion_s( sa, saswsb, saswst, wall_salinityflux ) |
---|
1383 | |
---|
1384 | CALL user_actions( 'sa-tendency' ) |
---|
1385 | |
---|
1386 | ! |
---|
1387 | !-- Prognostic equation for salinity |
---|
1388 | DO i = nxl, nxr |
---|
1389 | DO j = nys, nyn |
---|
1390 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1391 | sa_p(k,j,i) = sa(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
1392 | tsc(3) * tsa_m(k,j,i) ) & |
---|
1393 | - tsc(5) * rdf_sc(k) * & |
---|
1394 | ( sa(k,j,i) - sa_init(k) ) |
---|
1395 | IF ( sa_p(k,j,i) < 0.0_wp ) sa_p(k,j,i) = 0.1_wp * sa(k,j,i) |
---|
1396 | ENDDO |
---|
1397 | ENDDO |
---|
1398 | ENDDO |
---|
1399 | |
---|
1400 | ! |
---|
1401 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
1402 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1403 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
1404 | DO i = nxl, nxr |
---|
1405 | DO j = nys, nyn |
---|
1406 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1407 | tsa_m(k,j,i) = tend(k,j,i) |
---|
1408 | ENDDO |
---|
1409 | ENDDO |
---|
1410 | ENDDO |
---|
1411 | ELSEIF ( intermediate_timestep_count < & |
---|
1412 | intermediate_timestep_count_max ) THEN |
---|
1413 | DO i = nxl, nxr |
---|
1414 | DO j = nys, nyn |
---|
1415 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1416 | tsa_m(k,j,i) = -9.5625_wp * tend(k,j,i) + & |
---|
1417 | 5.3125_wp * tsa_m(k,j,i) |
---|
1418 | ENDDO |
---|
1419 | ENDDO |
---|
1420 | ENDDO |
---|
1421 | ENDIF |
---|
1422 | ENDIF |
---|
1423 | |
---|
1424 | CALL cpu_log( log_point(37), 'sa-equation', 'stop' ) |
---|
1425 | |
---|
1426 | ! |
---|
1427 | !-- Calculate density by the equation of state for seawater |
---|
1428 | CALL cpu_log( log_point(38), 'eqns-seawater', 'start' ) |
---|
1429 | CALL eqn_state_seawater |
---|
1430 | CALL cpu_log( log_point(38), 'eqns-seawater', 'stop' ) |
---|
1431 | |
---|
1432 | ENDIF |
---|
1433 | |
---|
1434 | ! |
---|
1435 | !-- If required, compute prognostic equation for total water content / scalar |
---|
1436 | IF ( humidity .OR. passive_scalar ) THEN |
---|
1437 | |
---|
1438 | CALL cpu_log( log_point(29), 'q/s-equation', 'start' ) |
---|
1439 | |
---|
1440 | ! |
---|
1441 | !-- Scalar/q-tendency terms with communication |
---|
1442 | sbt = tsc(2) |
---|
1443 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
1444 | |
---|
1445 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
1446 | ! |
---|
1447 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
1448 | sbt = 1.0_wp |
---|
1449 | ENDIF |
---|
1450 | tend = 0.0_wp |
---|
1451 | CALL advec_s_bc( q, 'q' ) |
---|
1452 | |
---|
1453 | ENDIF |
---|
1454 | |
---|
1455 | ! |
---|
1456 | !-- Scalar/q-tendency terms with no communication |
---|
1457 | IF ( scalar_advec /= 'bc-scheme' ) THEN |
---|
1458 | tend = 0.0_wp |
---|
1459 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1460 | IF ( ws_scheme_sca ) THEN |
---|
1461 | CALL advec_s_ws( q, 'q' ) |
---|
1462 | ELSE |
---|
1463 | CALL advec_s_pw( q ) |
---|
1464 | ENDIF |
---|
1465 | ELSE |
---|
1466 | CALL advec_s_up( q ) |
---|
1467 | ENDIF |
---|
1468 | ENDIF |
---|
1469 | |
---|
1470 | CALL diffusion_s( q, qsws, qswst, wall_qflux ) |
---|
1471 | |
---|
1472 | ! |
---|
1473 | !-- Sink or source of scalar concentration due to canopy elements |
---|
1474 | IF ( plant_canopy ) CALL pcm_tendency( 5 ) |
---|
1475 | |
---|
1476 | ! |
---|
1477 | !-- Large scale advection |
---|
1478 | IF ( large_scale_forcing ) THEN |
---|
1479 | CALL ls_advec( simulated_time, 'q' ) |
---|
1480 | ENDIF |
---|
1481 | |
---|
1482 | ! |
---|
1483 | !-- Nudging |
---|
1484 | IF ( nudging ) CALL nudge( simulated_time, 'q' ) |
---|
1485 | |
---|
1486 | ! |
---|
1487 | !-- If required compute influence of large-scale subsidence/ascent |
---|
1488 | IF ( large_scale_subsidence .AND. & |
---|
1489 | .NOT. use_subsidence_tendencies ) THEN |
---|
1490 | CALL subsidence( tend, q, q_init, 3 ) |
---|
1491 | ENDIF |
---|
1492 | |
---|
1493 | CALL user_actions( 'q-tendency' ) |
---|
1494 | |
---|
1495 | ! |
---|
1496 | !-- Prognostic equation for total water content / scalar |
---|
1497 | DO i = nxl, nxr |
---|
1498 | DO j = nys, nyn |
---|
1499 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1500 | q_p(k,j,i) = q(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
1501 | tsc(3) * tq_m(k,j,i) ) & |
---|
1502 | - tsc(5) * rdf_sc(k) * & |
---|
1503 | ( q(k,j,i) - q_init(k) ) |
---|
1504 | IF ( q_p(k,j,i) < 0.0_wp ) q_p(k,j,i) = 0.1_wp * q(k,j,i) |
---|
1505 | ENDDO |
---|
1506 | ENDDO |
---|
1507 | ENDDO |
---|
1508 | |
---|
1509 | ! |
---|
1510 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
1511 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1512 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
1513 | DO i = nxl, nxr |
---|
1514 | DO j = nys, nyn |
---|
1515 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1516 | tq_m(k,j,i) = tend(k,j,i) |
---|
1517 | ENDDO |
---|
1518 | ENDDO |
---|
1519 | ENDDO |
---|
1520 | ELSEIF ( intermediate_timestep_count < & |
---|
1521 | intermediate_timestep_count_max ) THEN |
---|
1522 | DO i = nxl, nxr |
---|
1523 | DO j = nys, nyn |
---|
1524 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1525 | tq_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tq_m(k,j,i) |
---|
1526 | ENDDO |
---|
1527 | ENDDO |
---|
1528 | ENDDO |
---|
1529 | ENDIF |
---|
1530 | ENDIF |
---|
1531 | |
---|
1532 | CALL cpu_log( log_point(29), 'q/s-equation', 'stop' ) |
---|
1533 | |
---|
1534 | ! |
---|
1535 | !-- If required, calculate prognostic equations for rain water content |
---|
1536 | !-- and rain drop concentration |
---|
1537 | IF ( cloud_physics .AND. microphysics_seifert ) THEN |
---|
1538 | |
---|
1539 | CALL cpu_log( log_point(52), 'qr-equation', 'start' ) |
---|
1540 | |
---|
1541 | ! |
---|
1542 | !-- Calculate prognostic equation for rain water content |
---|
1543 | sbt = tsc(2) |
---|
1544 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
1545 | |
---|
1546 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
1547 | ! |
---|
1548 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
1549 | sbt = 1.0_wp |
---|
1550 | ENDIF |
---|
1551 | tend = 0.0_wp |
---|
1552 | CALL advec_s_bc( qr, 'qr' ) |
---|
1553 | |
---|
1554 | ENDIF |
---|
1555 | |
---|
1556 | ! |
---|
1557 | !-- qr-tendency terms with no communication |
---|
1558 | IF ( scalar_advec /= 'bc-scheme' ) THEN |
---|
1559 | tend = 0.0_wp |
---|
1560 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1561 | IF ( ws_scheme_sca ) THEN |
---|
1562 | CALL advec_s_ws( qr, 'qr' ) |
---|
1563 | ELSE |
---|
1564 | CALL advec_s_pw( qr ) |
---|
1565 | ENDIF |
---|
1566 | ELSE |
---|
1567 | CALL advec_s_up( qr ) |
---|
1568 | ENDIF |
---|
1569 | ENDIF |
---|
1570 | |
---|
1571 | CALL diffusion_s( qr, qrsws, qrswst, wall_qrflux ) |
---|
1572 | |
---|
1573 | CALL user_actions( 'qr-tendency' ) |
---|
1574 | |
---|
1575 | ! |
---|
1576 | !-- Prognostic equation for rain water content |
---|
1577 | DO i = nxl, nxr |
---|
1578 | DO j = nys, nyn |
---|
1579 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1580 | qr_p(k,j,i) = qr(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
1581 | tsc(3) * tqr_m(k,j,i) ) & |
---|
1582 | - tsc(5) * rdf_sc(k) * qr(k,j,i) |
---|
1583 | IF ( qr_p(k,j,i) < 0.0_wp ) qr_p(k,j,i) = 0.0_wp |
---|
1584 | ENDDO |
---|
1585 | ENDDO |
---|
1586 | ENDDO |
---|
1587 | |
---|
1588 | ! |
---|
1589 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
1590 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1591 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
1592 | DO i = nxl, nxr |
---|
1593 | DO j = nys, nyn |
---|
1594 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1595 | tqr_m(k,j,i) = tend(k,j,i) |
---|
1596 | ENDDO |
---|
1597 | ENDDO |
---|
1598 | ENDDO |
---|
1599 | ELSEIF ( intermediate_timestep_count < & |
---|
1600 | intermediate_timestep_count_max ) THEN |
---|
1601 | DO i = nxl, nxr |
---|
1602 | DO j = nys, nyn |
---|
1603 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1604 | tqr_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * & |
---|
1605 | tqr_m(k,j,i) |
---|
1606 | ENDDO |
---|
1607 | ENDDO |
---|
1608 | ENDDO |
---|
1609 | ENDIF |
---|
1610 | ENDIF |
---|
1611 | |
---|
1612 | CALL cpu_log( log_point(52), 'qr-equation', 'stop' ) |
---|
1613 | CALL cpu_log( log_point(53), 'nr-equation', 'start' ) |
---|
1614 | |
---|
1615 | ! |
---|
1616 | !-- Calculate prognostic equation for rain drop concentration |
---|
1617 | sbt = tsc(2) |
---|
1618 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
1619 | |
---|
1620 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
1621 | ! |
---|
1622 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
1623 | sbt = 1.0_wp |
---|
1624 | ENDIF |
---|
1625 | tend = 0.0_wp |
---|
1626 | CALL advec_s_bc( nr, 'nr' ) |
---|
1627 | |
---|
1628 | ENDIF |
---|
1629 | |
---|
1630 | ! |
---|
1631 | !-- nr-tendency terms with no communication |
---|
1632 | IF ( scalar_advec /= 'bc-scheme' ) THEN |
---|
1633 | tend = 0.0_wp |
---|
1634 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1635 | IF ( ws_scheme_sca ) THEN |
---|
1636 | CALL advec_s_ws( nr, 'nr' ) |
---|
1637 | ELSE |
---|
1638 | CALL advec_s_pw( nr ) |
---|
1639 | ENDIF |
---|
1640 | ELSE |
---|
1641 | CALL advec_s_up( nr ) |
---|
1642 | ENDIF |
---|
1643 | ENDIF |
---|
1644 | |
---|
1645 | CALL diffusion_s( nr, nrsws, nrswst, wall_nrflux ) |
---|
1646 | |
---|
1647 | ! |
---|
1648 | !-- Prognostic equation for rain drop concentration |
---|
1649 | DO i = nxl, nxr |
---|
1650 | DO j = nys, nyn |
---|
1651 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1652 | nr_p(k,j,i) = nr(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
1653 | tsc(3) * tnr_m(k,j,i) ) & |
---|
1654 | - tsc(5) * rdf_sc(k) * nr(k,j,i) |
---|
1655 | IF ( nr_p(k,j,i) < 0.0_wp ) nr_p(k,j,i) = 0.0_wp |
---|
1656 | ENDDO |
---|
1657 | ENDDO |
---|
1658 | ENDDO |
---|
1659 | |
---|
1660 | ! |
---|
1661 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
1662 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1663 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
1664 | DO i = nxl, nxr |
---|
1665 | DO j = nys, nyn |
---|
1666 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1667 | tnr_m(k,j,i) = tend(k,j,i) |
---|
1668 | ENDDO |
---|
1669 | ENDDO |
---|
1670 | ENDDO |
---|
1671 | ELSEIF ( intermediate_timestep_count < & |
---|
1672 | intermediate_timestep_count_max ) THEN |
---|
1673 | DO i = nxl, nxr |
---|
1674 | DO j = nys, nyn |
---|
1675 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1676 | tnr_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * & |
---|
1677 | tnr_m(k,j,i) |
---|
1678 | ENDDO |
---|
1679 | ENDDO |
---|
1680 | ENDDO |
---|
1681 | ENDIF |
---|
1682 | ENDIF |
---|
1683 | |
---|
1684 | CALL cpu_log( log_point(53), 'nr-equation', 'stop' ) |
---|
1685 | |
---|
1686 | ENDIF |
---|
1687 | |
---|
1688 | ENDIF |
---|
1689 | |
---|
1690 | ! |
---|
1691 | !-- If required, compute prognostic equation for turbulent kinetic |
---|
1692 | !-- energy (TKE) |
---|
1693 | IF ( .NOT. constant_diffusion ) THEN |
---|
1694 | |
---|
1695 | CALL cpu_log( log_point(16), 'tke-equation', 'start' ) |
---|
1696 | |
---|
1697 | sbt = tsc(2) |
---|
1698 | IF ( .NOT. use_upstream_for_tke ) THEN |
---|
1699 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
1700 | |
---|
1701 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
1702 | ! |
---|
1703 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
1704 | sbt = 1.0_wp |
---|
1705 | ENDIF |
---|
1706 | tend = 0.0_wp |
---|
1707 | CALL advec_s_bc( e, 'e' ) |
---|
1708 | |
---|
1709 | ENDIF |
---|
1710 | ENDIF |
---|
1711 | |
---|
1712 | ! |
---|
1713 | !-- TKE-tendency terms with no communication |
---|
1714 | IF ( scalar_advec /= 'bc-scheme' .OR. use_upstream_for_tke ) THEN |
---|
1715 | IF ( use_upstream_for_tke ) THEN |
---|
1716 | tend = 0.0_wp |
---|
1717 | CALL advec_s_up( e ) |
---|
1718 | ELSE |
---|
1719 | tend = 0.0_wp |
---|
1720 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1721 | IF ( ws_scheme_sca ) THEN |
---|
1722 | CALL advec_s_ws( e, 'e' ) |
---|
1723 | ELSE |
---|
1724 | CALL advec_s_pw( e ) |
---|
1725 | ENDIF |
---|
1726 | ELSE |
---|
1727 | CALL advec_s_up( e ) |
---|
1728 | ENDIF |
---|
1729 | ENDIF |
---|
1730 | ENDIF |
---|
1731 | |
---|
1732 | IF ( .NOT. humidity ) THEN |
---|
1733 | IF ( ocean ) THEN |
---|
1734 | CALL diffusion_e( prho, prho_reference ) |
---|
1735 | ELSE |
---|
1736 | CALL diffusion_e( pt, pt_reference ) |
---|
1737 | ENDIF |
---|
1738 | ELSE |
---|
1739 | CALL diffusion_e( vpt, pt_reference ) |
---|
1740 | ENDIF |
---|
1741 | |
---|
1742 | CALL production_e |
---|
1743 | |
---|
1744 | ! |
---|
1745 | !-- Additional sink term for flows through plant canopies |
---|
1746 | IF ( plant_canopy ) CALL pcm_tendency( 6 ) |
---|
1747 | CALL user_actions( 'e-tendency' ) |
---|
1748 | |
---|
1749 | ! |
---|
1750 | !-- Prognostic equation for TKE. |
---|
1751 | !-- Eliminate negative TKE values, which can occur due to numerical |
---|
1752 | !-- reasons in the course of the integration. In such cases the old TKE |
---|
1753 | !-- value is reduced by 90%. |
---|
1754 | DO i = nxl, nxr |
---|
1755 | DO j = nys, nyn |
---|
1756 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1757 | e_p(k,j,i) = e(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
1758 | tsc(3) * te_m(k,j,i) ) |
---|
1759 | IF ( e_p(k,j,i) < 0.0_wp ) e_p(k,j,i) = 0.1_wp * e(k,j,i) |
---|
1760 | ENDDO |
---|
1761 | ENDDO |
---|
1762 | ENDDO |
---|
1763 | |
---|
1764 | ! |
---|
1765 | !-- Calculate tendencies for the next Runge-Kutta step |
---|
1766 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1767 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
1768 | DO i = nxl, nxr |
---|
1769 | DO j = nys, nyn |
---|
1770 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1771 | te_m(k,j,i) = tend(k,j,i) |
---|
1772 | ENDDO |
---|
1773 | ENDDO |
---|
1774 | ENDDO |
---|
1775 | ELSEIF ( intermediate_timestep_count < & |
---|
1776 | intermediate_timestep_count_max ) THEN |
---|
1777 | DO i = nxl, nxr |
---|
1778 | DO j = nys, nyn |
---|
1779 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
1780 | te_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * te_m(k,j,i) |
---|
1781 | ENDDO |
---|
1782 | ENDDO |
---|
1783 | ENDDO |
---|
1784 | ENDIF |
---|
1785 | ENDIF |
---|
1786 | |
---|
1787 | CALL cpu_log( log_point(16), 'tke-equation', 'stop' ) |
---|
1788 | |
---|
1789 | ENDIF |
---|
1790 | |
---|
1791 | END SUBROUTINE prognostic_equations_vector |
---|
1792 | |
---|
1793 | |
---|
1794 | !------------------------------------------------------------------------------! |
---|
1795 | ! Description: |
---|
1796 | ! ------------ |
---|
1797 | !> Version for accelerator boards |
---|
1798 | !------------------------------------------------------------------------------! |
---|
1799 | |
---|
1800 | SUBROUTINE prognostic_equations_acc |
---|
1801 | |
---|
1802 | |
---|
1803 | IMPLICIT NONE |
---|
1804 | |
---|
1805 | INTEGER(iwp) :: i !< |
---|
1806 | INTEGER(iwp) :: j !< |
---|
1807 | INTEGER(iwp) :: k !< |
---|
1808 | INTEGER(iwp) :: runge_step !< |
---|
1809 | |
---|
1810 | REAL(wp) :: sbt !< |
---|
1811 | |
---|
1812 | ! |
---|
1813 | !-- Set switch for intermediate Runge-Kutta step |
---|
1814 | runge_step = 0 |
---|
1815 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1816 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
1817 | runge_step = 1 |
---|
1818 | ELSEIF ( intermediate_timestep_count < & |
---|
1819 | intermediate_timestep_count_max ) THEN |
---|
1820 | runge_step = 2 |
---|
1821 | ENDIF |
---|
1822 | ENDIF |
---|
1823 | |
---|
1824 | ! |
---|
1825 | !-- If required, calculate cloud microphysical impacts (two-moment scheme) |
---|
1826 | IF ( cloud_physics .AND. .NOT. microphysics_sat_adjust .AND. & |
---|
1827 | ( intermediate_timestep_count == 1 .OR. & |
---|
1828 | call_microphysics_at_all_substeps ) & |
---|
1829 | ) THEN |
---|
1830 | CALL cpu_log( log_point(51), 'microphysics', 'start' ) |
---|
1831 | CALL microphysics_control |
---|
1832 | CALL cpu_log( log_point(51), 'microphysics', 'stop' ) |
---|
1833 | ENDIF |
---|
1834 | |
---|
1835 | ! |
---|
1836 | !-- u-velocity component |
---|
1837 | !++ Statistics still not completely ported to accelerators |
---|
1838 | !$acc update device( hom, ref_state ) |
---|
1839 | CALL cpu_log( log_point(5), 'u-equation', 'start' ) |
---|
1840 | |
---|
1841 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1842 | IF ( ws_scheme_mom ) THEN |
---|
1843 | CALL advec_u_ws_acc |
---|
1844 | ELSE |
---|
1845 | tend = 0.0_wp ! to be removed later?? |
---|
1846 | CALL advec_u_pw |
---|
1847 | ENDIF |
---|
1848 | ELSE |
---|
1849 | CALL advec_u_up |
---|
1850 | ENDIF |
---|
1851 | CALL diffusion_u_acc |
---|
1852 | CALL coriolis_acc( 1 ) |
---|
1853 | IF ( sloping_surface .AND. .NOT. neutral ) THEN |
---|
1854 | CALL buoyancy( pt, 1 ) |
---|
1855 | ENDIF |
---|
1856 | |
---|
1857 | ! |
---|
1858 | !-- Drag by plant canopy |
---|
1859 | IF ( plant_canopy ) CALL pcm_tendency( 1 ) |
---|
1860 | |
---|
1861 | ! |
---|
1862 | !-- External pressure gradient |
---|
1863 | IF ( dp_external ) THEN |
---|
1864 | DO i = i_left, i_right |
---|
1865 | DO j = j_south, j_north |
---|
1866 | DO k = dp_level_ind_b+1, nzt |
---|
1867 | tend(k,j,i) = tend(k,j,i) - dpdxy(1) * dp_smooth_factor(k) |
---|
1868 | ENDDO |
---|
1869 | ENDDO |
---|
1870 | ENDDO |
---|
1871 | ENDIF |
---|
1872 | |
---|
1873 | ! |
---|
1874 | !-- Nudging |
---|
1875 | IF ( nudging ) CALL nudge( simulated_time, 'u' ) |
---|
1876 | |
---|
1877 | ! |
---|
1878 | !-- Forces by wind turbines |
---|
1879 | IF ( wind_turbine ) CALL wtm_tendencies( 1 ) |
---|
1880 | |
---|
1881 | CALL user_actions( 'u-tendency' ) |
---|
1882 | |
---|
1883 | ! |
---|
1884 | !-- Prognostic equation for u-velocity component |
---|
1885 | !$acc kernels present( nzb_u_inner, rdf, tend, tu_m, u, u_init, u_p ) |
---|
1886 | !$acc loop independent |
---|
1887 | DO i = i_left, i_right |
---|
1888 | !$acc loop independent |
---|
1889 | DO j = j_south, j_north |
---|
1890 | !$acc loop independent |
---|
1891 | DO k = 1, nzt |
---|
1892 | IF ( k > nzb_u_inner(j,i) ) THEN |
---|
1893 | u_p(k,j,i) = u(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
1894 | tsc(3) * tu_m(k,j,i) ) & |
---|
1895 | - tsc(5) * rdf(k) * ( u(k,j,i) - u_init(k) ) |
---|
1896 | ! |
---|
1897 | !-- Tendencies for the next Runge-Kutta step |
---|
1898 | IF ( runge_step == 1 ) THEN |
---|
1899 | tu_m(k,j,i) = tend(k,j,i) |
---|
1900 | ELSEIF ( runge_step == 2 ) THEN |
---|
1901 | tu_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tu_m(k,j,i) |
---|
1902 | ENDIF |
---|
1903 | ENDIF |
---|
1904 | ENDDO |
---|
1905 | ENDDO |
---|
1906 | ENDDO |
---|
1907 | !$acc end kernels |
---|
1908 | |
---|
1909 | CALL cpu_log( log_point(5), 'u-equation', 'stop' ) |
---|
1910 | |
---|
1911 | ! |
---|
1912 | !-- v-velocity component |
---|
1913 | CALL cpu_log( log_point(6), 'v-equation', 'start' ) |
---|
1914 | |
---|
1915 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1916 | IF ( ws_scheme_mom ) THEN |
---|
1917 | CALL advec_v_ws_acc |
---|
1918 | ELSE |
---|
1919 | tend = 0.0_wp ! to be removed later?? |
---|
1920 | CALL advec_v_pw |
---|
1921 | END IF |
---|
1922 | ELSE |
---|
1923 | CALL advec_v_up |
---|
1924 | ENDIF |
---|
1925 | CALL diffusion_v_acc |
---|
1926 | CALL coriolis_acc( 2 ) |
---|
1927 | |
---|
1928 | ! |
---|
1929 | !-- Drag by plant canopy |
---|
1930 | IF ( plant_canopy ) CALL pcm_tendency( 2 ) |
---|
1931 | |
---|
1932 | ! |
---|
1933 | !-- External pressure gradient |
---|
1934 | IF ( dp_external ) THEN |
---|
1935 | DO i = i_left, i_right |
---|
1936 | DO j = j_south, j_north |
---|
1937 | DO k = dp_level_ind_b+1, nzt |
---|
1938 | tend(k,j,i) = tend(k,j,i) - dpdxy(2) * dp_smooth_factor(k) |
---|
1939 | ENDDO |
---|
1940 | ENDDO |
---|
1941 | ENDDO |
---|
1942 | ENDIF |
---|
1943 | |
---|
1944 | ! |
---|
1945 | !-- Nudging |
---|
1946 | IF ( nudging ) CALL nudge( simulated_time, 'v' ) |
---|
1947 | |
---|
1948 | ! |
---|
1949 | !-- Forces by wind turbines |
---|
1950 | IF ( wind_turbine ) CALL wtm_tendencies( 2 ) |
---|
1951 | |
---|
1952 | CALL user_actions( 'v-tendency' ) |
---|
1953 | |
---|
1954 | ! |
---|
1955 | !-- Prognostic equation for v-velocity component |
---|
1956 | !$acc kernels present( nzb_v_inner, rdf, tend, tv_m, v, v_init, v_p ) |
---|
1957 | !$acc loop independent |
---|
1958 | DO i = i_left, i_right |
---|
1959 | !$acc loop independent |
---|
1960 | DO j = j_south, j_north |
---|
1961 | !$acc loop independent |
---|
1962 | DO k = 1, nzt |
---|
1963 | IF ( k > nzb_v_inner(j,i) ) THEN |
---|
1964 | v_p(k,j,i) = v(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
1965 | tsc(3) * tv_m(k,j,i) ) & |
---|
1966 | - tsc(5) * rdf(k) * ( v(k,j,i) - v_init(k) ) |
---|
1967 | ! |
---|
1968 | !-- Tendencies for the next Runge-Kutta step |
---|
1969 | IF ( runge_step == 1 ) THEN |
---|
1970 | tv_m(k,j,i) = tend(k,j,i) |
---|
1971 | ELSEIF ( runge_step == 2 ) THEN |
---|
1972 | tv_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tv_m(k,j,i) |
---|
1973 | ENDIF |
---|
1974 | ENDIF |
---|
1975 | ENDDO |
---|
1976 | ENDDO |
---|
1977 | ENDDO |
---|
1978 | !$acc end kernels |
---|
1979 | |
---|
1980 | CALL cpu_log( log_point(6), 'v-equation', 'stop' ) |
---|
1981 | |
---|
1982 | ! |
---|
1983 | !-- w-velocity component |
---|
1984 | CALL cpu_log( log_point(7), 'w-equation', 'start' ) |
---|
1985 | |
---|
1986 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
1987 | IF ( ws_scheme_mom ) THEN |
---|
1988 | CALL advec_w_ws_acc |
---|
1989 | ELSE |
---|
1990 | tend = 0.0_wp ! to be removed later?? |
---|
1991 | CALL advec_w_pw |
---|
1992 | ENDIF |
---|
1993 | ELSE |
---|
1994 | CALL advec_w_up |
---|
1995 | ENDIF |
---|
1996 | CALL diffusion_w_acc |
---|
1997 | CALL coriolis_acc( 3 ) |
---|
1998 | |
---|
1999 | IF ( .NOT. neutral ) THEN |
---|
2000 | IF ( ocean ) THEN |
---|
2001 | CALL buoyancy( rho, 3 ) |
---|
2002 | ELSE |
---|
2003 | IF ( .NOT. humidity ) THEN |
---|
2004 | CALL buoyancy_acc( pt, 3 ) |
---|
2005 | ELSE |
---|
2006 | CALL buoyancy( vpt, 3 ) |
---|
2007 | ENDIF |
---|
2008 | ENDIF |
---|
2009 | ENDIF |
---|
2010 | |
---|
2011 | ! |
---|
2012 | !-- Drag by plant canopy |
---|
2013 | IF ( plant_canopy ) CALL pcm_tendency( 3 ) |
---|
2014 | |
---|
2015 | ! |
---|
2016 | !-- Forces by wind turbines |
---|
2017 | IF ( wind_turbine ) CALL wtm_tendencies( 3 ) |
---|
2018 | |
---|
2019 | CALL user_actions( 'w-tendency' ) |
---|
2020 | |
---|
2021 | ! |
---|
2022 | !-- Prognostic equation for w-velocity component |
---|
2023 | !$acc kernels present( nzb_w_inner, rdf, tend, tw_m, w, w_p ) |
---|
2024 | !$acc loop independent |
---|
2025 | DO i = i_left, i_right |
---|
2026 | !$acc loop independent |
---|
2027 | DO j = j_south, j_north |
---|
2028 | !$acc loop independent |
---|
2029 | DO k = 1, nzt-1 |
---|
2030 | IF ( k > nzb_w_inner(j,i) ) THEN |
---|
2031 | w_p(k,j,i) = w(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & |
---|
2032 | tsc(3) * tw_m(k,j,i) ) & |
---|
2033 | - tsc(5) * rdf(k) * w(k,j,i) |
---|
2034 | ! |
---|
2035 | !-- Tendencies for the next Runge-Kutta step |
---|
2036 | IF ( runge_step == 1 ) THEN |
---|
2037 | tw_m(k,j,i) = tend(k,j,i) |
---|
2038 | ELSEIF ( runge_step == 2 ) THEN |
---|
2039 | tw_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tw_m(k,j,i) |
---|
2040 | ENDIF |
---|
2041 | ENDIF |
---|
2042 | ENDDO |
---|
2043 | ENDDO |
---|
2044 | ENDDO |
---|
2045 | !$acc end kernels |
---|
2046 | |
---|
2047 | CALL cpu_log( log_point(7), 'w-equation', 'stop' ) |
---|
2048 | |
---|
2049 | |
---|
2050 | ! |
---|
2051 | !-- If required, compute prognostic equation for potential temperature |
---|
2052 | IF ( .NOT. neutral ) THEN |
---|
2053 | |
---|
2054 | CALL cpu_log( log_point(13), 'pt-equation', 'start' ) |
---|
2055 | |
---|
2056 | ! |
---|
2057 | !-- pt-tendency terms with communication |
---|
2058 | sbt = tsc(2) |
---|
2059 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
2060 | |
---|
2061 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
2062 | ! |
---|
2063 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
2064 | sbt = 1.0_wp |
---|
2065 | ENDIF |
---|
2066 | tend = 0.0_wp |
---|
2067 | CALL advec_s_bc( pt, 'pt' ) |
---|
2068 | |
---|
2069 | ENDIF |
---|
2070 | |
---|
2071 | ! |
---|
2072 | !-- pt-tendency terms with no communication |
---|
2073 | IF ( scalar_advec /= 'bc-scheme' ) THEN |
---|
2074 | tend = 0.0_wp |
---|
2075 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
2076 | IF ( ws_scheme_sca ) THEN |
---|
2077 | CALL advec_s_ws_acc( pt, 'pt' ) |
---|
2078 | ELSE |
---|
2079 | tend = 0.0_wp ! to be removed later?? |
---|
2080 | CALL advec_s_pw( pt ) |
---|
2081 | ENDIF |
---|
2082 | ELSE |
---|
2083 | CALL advec_s_up( pt ) |
---|
2084 | ENDIF |
---|
2085 | ENDIF |
---|
2086 | |
---|
2087 | CALL diffusion_s_acc( pt, shf, tswst, wall_heatflux ) |
---|
2088 | |
---|
2089 | ! |
---|
2090 | !-- If required compute heating/cooling due to long wave radiation processes |
---|
2091 | IF ( cloud_top_radiation ) THEN |
---|
2092 | CALL calc_radiation |
---|
2093 | ENDIF |
---|
2094 | |
---|
2095 | ! |
---|
2096 | !-- Consideration of heat sources within the plant canopy |
---|
2097 | IF ( plant_canopy .AND. ( cthf /= 0.0_wp ) ) THEN |
---|
2098 | CALL pcm_tendency( 4 ) |
---|
2099 | ENDIF |
---|
2100 | |
---|
2101 | ! |
---|
2102 | !-- Large scale advection |
---|
2103 | IF ( large_scale_forcing ) THEN |
---|
2104 | CALL ls_advec( simulated_time, 'pt' ) |
---|
2105 | ENDIF |
---|
2106 | |
---|
2107 | ! |
---|
2108 | !-- Nudging |
---|
2109 | IF ( nudging ) CALL nudge( simulated_time, 'pt' ) |
---|
2110 | |
---|
2111 | ! |
---|
2112 | !-- If required compute influence of large-scale subsidence/ascent |
---|
2113 | IF ( large_scale_subsidence .AND. & |
---|
2114 | .NOT. use_subsidence_tendencies ) THEN |
---|
2115 | CALL subsidence( tend, pt, pt_init, 2 ) |
---|
2116 | ENDIF |
---|
2117 | |
---|
2118 | IF ( radiation .AND. radiation_scheme == 'rrtmg' .AND. & |
---|
2119 | simulated_time > skip_time_do_radiation ) THEN |
---|
2120 | CALL radiation_tendency ( tend ) |
---|
2121 | ENDIF |
---|
2122 | |
---|
2123 | CALL user_actions( 'pt-tendency' ) |
---|
2124 | |
---|
2125 | ! |
---|
2126 | !-- Prognostic equation for potential temperature |
---|
2127 | !$acc kernels present( nzb_s_inner, rdf_sc, ptdf_x, ptdf_y, pt_init ) & |
---|
2128 | !$acc present( tend, tpt_m, pt, pt_p ) |
---|
2129 | !$acc loop independent |
---|
2130 | DO i = i_left, i_right |
---|
2131 | !$acc loop independent |
---|
2132 | DO j = j_south, j_north |
---|
2133 | !$acc loop independent |
---|
2134 | DO k = 1, nzt |
---|
2135 | IF ( k > nzb_s_inner(j,i) ) THEN |
---|
2136 | pt_p(k,j,i) = pt(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
2137 | tsc(3) * tpt_m(k,j,i) ) & |
---|
2138 | - tsc(5) * ( pt(k,j,i) - pt_init(k) ) *& |
---|
2139 | ( rdf_sc(k) + ptdf_x(i) + ptdf_y(j) ) |
---|
2140 | ! |
---|
2141 | !-- Tendencies for the next Runge-Kutta step |
---|
2142 | IF ( runge_step == 1 ) THEN |
---|
2143 | tpt_m(k,j,i) = tend(k,j,i) |
---|
2144 | ELSEIF ( runge_step == 2 ) THEN |
---|
2145 | tpt_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tpt_m(k,j,i) |
---|
2146 | ENDIF |
---|
2147 | ENDIF |
---|
2148 | ENDDO |
---|
2149 | ENDDO |
---|
2150 | ENDDO |
---|
2151 | !$acc end kernels |
---|
2152 | |
---|
2153 | CALL cpu_log( log_point(13), 'pt-equation', 'stop' ) |
---|
2154 | |
---|
2155 | ENDIF |
---|
2156 | |
---|
2157 | ! |
---|
2158 | !-- If required, compute prognostic equation for salinity |
---|
2159 | IF ( ocean ) THEN |
---|
2160 | |
---|
2161 | CALL cpu_log( log_point(37), 'sa-equation', 'start' ) |
---|
2162 | |
---|
2163 | ! |
---|
2164 | !-- sa-tendency terms with communication |
---|
2165 | sbt = tsc(2) |
---|
2166 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
2167 | |
---|
2168 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
2169 | ! |
---|
2170 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
2171 | sbt = 1.0_wp |
---|
2172 | ENDIF |
---|
2173 | tend = 0.0_wp |
---|
2174 | CALL advec_s_bc( sa, 'sa' ) |
---|
2175 | |
---|
2176 | ENDIF |
---|
2177 | |
---|
2178 | ! |
---|
2179 | !-- sa-tendency terms with no communication |
---|
2180 | IF ( scalar_advec /= 'bc-scheme' ) THEN |
---|
2181 | tend = 0.0_wp |
---|
2182 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
2183 | IF ( ws_scheme_sca ) THEN |
---|
2184 | CALL advec_s_ws( sa, 'sa' ) |
---|
2185 | ELSE |
---|
2186 | CALL advec_s_pw( sa ) |
---|
2187 | ENDIF |
---|
2188 | ELSE |
---|
2189 | CALL advec_s_up( sa ) |
---|
2190 | ENDIF |
---|
2191 | ENDIF |
---|
2192 | |
---|
2193 | CALL diffusion_s( sa, saswsb, saswst, wall_salinityflux ) |
---|
2194 | |
---|
2195 | CALL user_actions( 'sa-tendency' ) |
---|
2196 | |
---|
2197 | ! |
---|
2198 | !-- Prognostic equation for salinity |
---|
2199 | DO i = i_left, i_right |
---|
2200 | DO j = j_south, j_north |
---|
2201 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
2202 | sa_p(k,j,i) = sa(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
2203 | tsc(3) * tsa_m(k,j,i) ) & |
---|
2204 | - tsc(5) * rdf_sc(k) * & |
---|
2205 | ( sa(k,j,i) - sa_init(k) ) |
---|
2206 | IF ( sa_p(k,j,i) < 0.0_wp ) sa_p(k,j,i) = 0.1_wp * sa(k,j,i) |
---|
2207 | ! |
---|
2208 | !-- Tendencies for the next Runge-Kutta step |
---|
2209 | IF ( runge_step == 1 ) THEN |
---|
2210 | tsa_m(k,j,i) = tend(k,j,i) |
---|
2211 | ELSEIF ( runge_step == 2 ) THEN |
---|
2212 | tsa_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tsa_m(k,j,i) |
---|
2213 | ENDIF |
---|
2214 | ENDDO |
---|
2215 | ENDDO |
---|
2216 | ENDDO |
---|
2217 | |
---|
2218 | CALL cpu_log( log_point(37), 'sa-equation', 'stop' ) |
---|
2219 | |
---|
2220 | ! |
---|
2221 | !-- Calculate density by the equation of state for seawater |
---|
2222 | CALL cpu_log( log_point(38), 'eqns-seawater', 'start' ) |
---|
2223 | CALL eqn_state_seawater |
---|
2224 | CALL cpu_log( log_point(38), 'eqns-seawater', 'stop' ) |
---|
2225 | |
---|
2226 | ENDIF |
---|
2227 | |
---|
2228 | ! |
---|
2229 | !-- If required, compute prognostic equation for total water content / scalar |
---|
2230 | IF ( humidity .OR. passive_scalar ) THEN |
---|
2231 | |
---|
2232 | CALL cpu_log( log_point(29), 'q/s-equation', 'start' ) |
---|
2233 | |
---|
2234 | ! |
---|
2235 | !-- Scalar/q-tendency terms with communication |
---|
2236 | sbt = tsc(2) |
---|
2237 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
2238 | |
---|
2239 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
2240 | ! |
---|
2241 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
2242 | sbt = 1.0_wp |
---|
2243 | ENDIF |
---|
2244 | tend = 0.0_wp |
---|
2245 | CALL advec_s_bc( q, 'q' ) |
---|
2246 | |
---|
2247 | ENDIF |
---|
2248 | |
---|
2249 | ! |
---|
2250 | !-- Scalar/q-tendency terms with no communication |
---|
2251 | IF ( scalar_advec /= 'bc-scheme' ) THEN |
---|
2252 | tend = 0.0_wp |
---|
2253 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
2254 | IF ( ws_scheme_sca ) THEN |
---|
2255 | CALL advec_s_ws( q, 'q' ) |
---|
2256 | ELSE |
---|
2257 | CALL advec_s_pw( q ) |
---|
2258 | ENDIF |
---|
2259 | ELSE |
---|
2260 | CALL advec_s_up( q ) |
---|
2261 | ENDIF |
---|
2262 | ENDIF |
---|
2263 | |
---|
2264 | CALL diffusion_s( q, qsws, qswst, wall_qflux ) |
---|
2265 | |
---|
2266 | ! |
---|
2267 | !-- Sink or source of scalar concentration due to canopy elements |
---|
2268 | IF ( plant_canopy ) CALL pcm_tendency( 5 ) |
---|
2269 | |
---|
2270 | ! |
---|
2271 | !-- Large scale advection |
---|
2272 | IF ( large_scale_forcing ) THEN |
---|
2273 | CALL ls_advec( simulated_time, 'q' ) |
---|
2274 | ENDIF |
---|
2275 | |
---|
2276 | ! |
---|
2277 | !-- Nudging |
---|
2278 | IF ( nudging ) CALL nudge( simulated_time, 'q' ) |
---|
2279 | |
---|
2280 | ! |
---|
2281 | !-- If required compute influence of large-scale subsidence/ascent |
---|
2282 | IF ( large_scale_subsidence .AND. & |
---|
2283 | .NOT. use_subsidence_tendencies ) THEN |
---|
2284 | CALL subsidence( tend, q, q_init, 3 ) |
---|
2285 | ENDIF |
---|
2286 | |
---|
2287 | CALL user_actions( 'q-tendency' ) |
---|
2288 | |
---|
2289 | ! |
---|
2290 | !-- Prognostic equation for total water content / scalar |
---|
2291 | DO i = i_left, i_right |
---|
2292 | DO j = j_south, j_north |
---|
2293 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
2294 | q_p(k,j,i) = q(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
2295 | tsc(3) * tq_m(k,j,i) ) & |
---|
2296 | - tsc(5) * rdf_sc(k) * & |
---|
2297 | ( q(k,j,i) - q_init(k) ) |
---|
2298 | IF ( q_p(k,j,i) < 0.0_wp ) q_p(k,j,i) = 0.1_wp * q(k,j,i) |
---|
2299 | ! |
---|
2300 | !-- Tendencies for the next Runge-Kutta step |
---|
2301 | IF ( runge_step == 1 ) THEN |
---|
2302 | tq_m(k,j,i) = tend(k,j,i) |
---|
2303 | ELSEIF ( runge_step == 2 ) THEN |
---|
2304 | tq_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * tq_m(k,j,i) |
---|
2305 | ENDIF |
---|
2306 | ENDDO |
---|
2307 | ENDDO |
---|
2308 | ENDDO |
---|
2309 | |
---|
2310 | CALL cpu_log( log_point(29), 'q/s-equation', 'stop' ) |
---|
2311 | |
---|
2312 | ! |
---|
2313 | !-- If required, calculate prognostic equations for rain water content |
---|
2314 | !-- and rain drop concentration |
---|
2315 | IF ( cloud_physics .AND. microphysics_seifert ) THEN |
---|
2316 | |
---|
2317 | CALL cpu_log( log_point(52), 'qr-equation', 'start' ) |
---|
2318 | ! |
---|
2319 | !-- qr-tendency terms with communication |
---|
2320 | sbt = tsc(2) |
---|
2321 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
2322 | |
---|
2323 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
2324 | ! |
---|
2325 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
2326 | sbt = 1.0_wp |
---|
2327 | ENDIF |
---|
2328 | tend = 0.0_wp |
---|
2329 | CALL advec_s_bc( qr, 'qr' ) |
---|
2330 | |
---|
2331 | ENDIF |
---|
2332 | |
---|
2333 | ! |
---|
2334 | !-- qr-tendency terms with no communication |
---|
2335 | IF ( scalar_advec /= 'bc-scheme' ) THEN |
---|
2336 | tend = 0.0_wp |
---|
2337 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
2338 | IF ( ws_scheme_sca ) THEN |
---|
2339 | CALL advec_s_ws( qr, 'qr' ) |
---|
2340 | ELSE |
---|
2341 | CALL advec_s_pw( qr ) |
---|
2342 | ENDIF |
---|
2343 | ELSE |
---|
2344 | CALL advec_s_up( qr ) |
---|
2345 | ENDIF |
---|
2346 | ENDIF |
---|
2347 | |
---|
2348 | CALL diffusion_s( qr, qrsws, qrswst, wall_qrflux ) |
---|
2349 | |
---|
2350 | ! |
---|
2351 | !-- Prognostic equation for rain water content |
---|
2352 | DO i = i_left, i_right |
---|
2353 | DO j = j_south, j_north |
---|
2354 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
2355 | qr_p(k,j,i) = qr(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
2356 | tsc(3) * tqr_m(k,j,i) ) & |
---|
2357 | - tsc(5) * rdf_sc(k) * qr(k,j,i) |
---|
2358 | IF ( qr_p(k,j,i) < 0.0_wp ) qr_p(k,j,i) = 0.0_wp |
---|
2359 | ! |
---|
2360 | !-- Tendencies for the next Runge-Kutta step |
---|
2361 | IF ( runge_step == 1 ) THEN |
---|
2362 | tqr_m(k,j,i) = tend(k,j,i) |
---|
2363 | ELSEIF ( runge_step == 2 ) THEN |
---|
2364 | tqr_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * & |
---|
2365 | tqr_m(k,j,i) |
---|
2366 | ENDIF |
---|
2367 | ENDDO |
---|
2368 | ENDDO |
---|
2369 | ENDDO |
---|
2370 | |
---|
2371 | CALL cpu_log( log_point(52), 'qr-equation', 'stop' ) |
---|
2372 | CALL cpu_log( log_point(53), 'nr-equation', 'start' ) |
---|
2373 | |
---|
2374 | ! |
---|
2375 | !-- nr-tendency terms with communication |
---|
2376 | sbt = tsc(2) |
---|
2377 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
2378 | |
---|
2379 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
2380 | ! |
---|
2381 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
2382 | sbt = 1.0_wp |
---|
2383 | ENDIF |
---|
2384 | tend = 0.0_wp |
---|
2385 | CALL advec_s_bc( nr, 'nr' ) |
---|
2386 | |
---|
2387 | ENDIF |
---|
2388 | |
---|
2389 | ! |
---|
2390 | !-- nr-tendency terms with no communication |
---|
2391 | IF ( scalar_advec /= 'bc-scheme' ) THEN |
---|
2392 | tend = 0.0_wp |
---|
2393 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
2394 | IF ( ws_scheme_sca ) THEN |
---|
2395 | CALL advec_s_ws( nr, 'nr' ) |
---|
2396 | ELSE |
---|
2397 | CALL advec_s_pw( nr ) |
---|
2398 | ENDIF |
---|
2399 | ELSE |
---|
2400 | CALL advec_s_up( nr ) |
---|
2401 | ENDIF |
---|
2402 | ENDIF |
---|
2403 | |
---|
2404 | CALL diffusion_s( nr, nrsws, nrswst, wall_nrflux ) |
---|
2405 | |
---|
2406 | ! |
---|
2407 | !-- Prognostic equation for rain drop concentration |
---|
2408 | DO i = i_left, i_right |
---|
2409 | DO j = j_south, j_north |
---|
2410 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
2411 | nr_p(k,j,i) = nr(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
2412 | tsc(3) * tnr_m(k,j,i) ) & |
---|
2413 | - tsc(5) * rdf_sc(k) * nr(k,j,i) |
---|
2414 | IF ( nr_p(k,j,i) < 0.0_wp ) nr_p(k,j,i) = 0.0_wp |
---|
2415 | ! |
---|
2416 | !-- Tendencies for the next Runge-Kutta step |
---|
2417 | IF ( runge_step == 1 ) THEN |
---|
2418 | tnr_m(k,j,i) = tend(k,j,i) |
---|
2419 | ELSEIF ( runge_step == 2 ) THEN |
---|
2420 | tnr_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * & |
---|
2421 | tnr_m(k,j,i) |
---|
2422 | ENDIF |
---|
2423 | ENDDO |
---|
2424 | ENDDO |
---|
2425 | ENDDO |
---|
2426 | |
---|
2427 | CALL cpu_log( log_point(53), 'nr-equation', 'stop' ) |
---|
2428 | |
---|
2429 | ENDIF |
---|
2430 | |
---|
2431 | ENDIF |
---|
2432 | |
---|
2433 | ! |
---|
2434 | !-- If required, compute prognostic equation for turbulent kinetic |
---|
2435 | !-- energy (TKE) |
---|
2436 | IF ( .NOT. constant_diffusion ) THEN |
---|
2437 | |
---|
2438 | CALL cpu_log( log_point(16), 'tke-equation', 'start' ) |
---|
2439 | |
---|
2440 | sbt = tsc(2) |
---|
2441 | IF ( .NOT. use_upstream_for_tke ) THEN |
---|
2442 | IF ( scalar_advec == 'bc-scheme' ) THEN |
---|
2443 | |
---|
2444 | IF ( timestep_scheme(1:5) /= 'runge' ) THEN |
---|
2445 | ! |
---|
2446 | !-- Bott-Chlond scheme always uses Euler time step. Thus: |
---|
2447 | sbt = 1.0_wp |
---|
2448 | ENDIF |
---|
2449 | tend = 0.0_wp |
---|
2450 | CALL advec_s_bc( e, 'e' ) |
---|
2451 | |
---|
2452 | ENDIF |
---|
2453 | ENDIF |
---|
2454 | |
---|
2455 | ! |
---|
2456 | !-- TKE-tendency terms with no communication |
---|
2457 | IF ( scalar_advec /= 'bc-scheme' .OR. use_upstream_for_tke ) THEN |
---|
2458 | IF ( use_upstream_for_tke ) THEN |
---|
2459 | tend = 0.0_wp |
---|
2460 | CALL advec_s_up( e ) |
---|
2461 | ELSE |
---|
2462 | IF ( timestep_scheme(1:5) == 'runge' ) THEN |
---|
2463 | IF ( ws_scheme_sca ) THEN |
---|
2464 | CALL advec_s_ws_acc( e, 'e' ) |
---|
2465 | ELSE |
---|
2466 | tend = 0.0_wp ! to be removed later?? |
---|
2467 | CALL advec_s_pw( e ) |
---|
2468 | ENDIF |
---|
2469 | ELSE |
---|
2470 | tend = 0.0_wp ! to be removed later?? |
---|
2471 | CALL advec_s_up( e ) |
---|
2472 | ENDIF |
---|
2473 | ENDIF |
---|
2474 | ENDIF |
---|
2475 | |
---|
2476 | IF ( .NOT. humidity ) THEN |
---|
2477 | IF ( ocean ) THEN |
---|
2478 | CALL diffusion_e( prho, prho_reference ) |
---|
2479 | ELSE |
---|
2480 | CALL diffusion_e_acc( pt, pt_reference ) |
---|
2481 | ENDIF |
---|
2482 | ELSE |
---|
2483 | CALL diffusion_e( vpt, pt_reference ) |
---|
2484 | ENDIF |
---|
2485 | |
---|
2486 | CALL production_e_acc |
---|
2487 | |
---|
2488 | ! |
---|
2489 | !-- Additional sink term for flows through plant canopies |
---|
2490 | IF ( plant_canopy ) CALL pcm_tendency( 6 ) |
---|
2491 | CALL user_actions( 'e-tendency' ) |
---|
2492 | |
---|
2493 | ! |
---|
2494 | !-- Prognostic equation for TKE. |
---|
2495 | !-- Eliminate negative TKE values, which can occur due to numerical |
---|
2496 | !-- reasons in the course of the integration. In such cases the old TKE |
---|
2497 | !-- value is reduced by 90%. |
---|
2498 | !$acc kernels present( e, e_p, nzb_s_inner, tend, te_m ) |
---|
2499 | !$acc loop independent |
---|
2500 | DO i = i_left, i_right |
---|
2501 | !$acc loop independent |
---|
2502 | DO j = j_south, j_north |
---|
2503 | !$acc loop independent |
---|
2504 | DO k = 1, nzt |
---|
2505 | IF ( k > nzb_s_inner(j,i) ) THEN |
---|
2506 | e_p(k,j,i) = e(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & |
---|
2507 | tsc(3) * te_m(k,j,i) ) |
---|
2508 | IF ( e_p(k,j,i) < 0.0_wp ) e_p(k,j,i) = 0.1_wp * e(k,j,i) |
---|
2509 | ! |
---|
2510 | !-- Tendencies for the next Runge-Kutta step |
---|
2511 | IF ( runge_step == 1 ) THEN |
---|
2512 | te_m(k,j,i) = tend(k,j,i) |
---|
2513 | ELSEIF ( runge_step == 2 ) THEN |
---|
2514 | te_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * te_m(k,j,i) |
---|
2515 | ENDIF |
---|
2516 | ENDIF |
---|
2517 | ENDDO |
---|
2518 | ENDDO |
---|
2519 | ENDDO |
---|
2520 | !$acc end kernels |
---|
2521 | |
---|
2522 | CALL cpu_log( log_point(16), 'tke-equation', 'stop' ) |
---|
2523 | |
---|
2524 | ENDIF |
---|
2525 | |
---|
2526 | END SUBROUTINE prognostic_equations_acc |
---|
2527 | |
---|
2528 | |
---|
2529 | END MODULE prognostic_equations_mod |
---|