source: palm/trunk/SOURCE/time_integration.f90 @ 2270

Last change on this file since 2270 was 2259, checked in by gronemeier, 7 years ago

Implemented synthetic turbulence generator

  • Property svn:keywords set to Id
File size: 42.1 KB
RevLine 
[1682]1!> @file time_integration.f90
[2000]2!------------------------------------------------------------------------------!
[1036]3! This file is part of PALM.
4!
[2000]5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
[1036]9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
[2101]17! Copyright 1997-2017 Leibniz Universitaet Hannover
[2000]18!------------------------------------------------------------------------------!
[1036]19!
[484]20! Current revisions:
[1092]21! ------------------
[1928]22!
[2233]23!
[1366]24! Former revisions:
25! -----------------
26! $Id: time_integration.f90 2259 2017-06-08 09:09:11Z maronga $
[2259]27! Implemented synthetic turbulence generator
[1366]28!
[2259]29! 2233 2017-05-30 18:08:54Z suehring
30!
[2233]31! 2232 2017-05-30 17:47:52Z suehring
32! Adjustments to new topography and surface concept
33! Modify passed parameters for disturb_field
34!
[2179]35! 2178 2017-03-17 11:07:39Z hellstea
36! Setting perturbations at all times near inflow boundary is removed
37! in case of nested boundaries
38!
[2175]39! 2174 2017-03-13 08:18:57Z maronga
40! Added support for nesting with cloud microphysics
41!
[2119]42! 2118 2017-01-17 16:38:49Z raasch
43! OpenACC directives and related code removed
44!
[2051]45! 2050 2016-11-08 15:00:55Z gronemeier
46! Implement turbulent outflow condition
47!
[2032]48! 2031 2016-10-21 15:11:58Z knoop
49! renamed variable rho to rho_ocean
50!
[2012]51! 2011 2016-09-19 17:29:57Z kanani
52! Flag urban_surface is now defined in module control_parameters,
53! removed commented CALLs of global_min_max.
54!
[2008]55! 2007 2016-08-24 15:47:17Z kanani
56! Added CALLs for new urban surface model
57!
[2001]58! 2000 2016-08-20 18:09:15Z knoop
59! Forced header and separation lines into 80 columns
60!
[1977]61! 1976 2016-07-27 13:28:04Z maronga
62! Simplified calls to radiation model
63!
[1961]64! 1960 2016-07-12 16:34:24Z suehring
65! Separate humidity and passive scalar
66!
[1958]67! 1957 2016-07-07 10:43:48Z suehring
68! flight module added
69!
[1933]70! 1919 2016-05-27 14:51:23Z raasch
71! Initial version of purely vertical nesting introduced.
[1928]72!
[1919]73! 1918 2016-05-27 14:35:57Z raasch
74! determination of time step moved to the end of the time step loop,
75! the first time step is now always calculated before the time step loop (i.e.
76! also in case of restart runs)
77!
[1917]78! 1914 2016-05-26 14:44:07Z witha
79! Added call for wind turbine model
80!
[1879]81! 1878 2016-04-19 12:30:36Z hellstea
82! Synchronization for nested runs rewritten
83!
[1854]84! 1853 2016-04-11 09:00:35Z maronga
85! Adjusted for use with radiation_scheme = constant
86!
[1851]87! 1849 2016-04-08 11:33:18Z hoffmann
88! Adapted for modularization of microphysics
89!
[1834]90! 1833 2016-04-07 14:23:03Z raasch
91! spectrum renamed spectra_mod, spectra related variables moved to spectra_mod
92!
[1832]93! 1831 2016-04-07 13:15:51Z hoffmann
94! turbulence renamed collision_turbulence
95!
[1823]96! 1822 2016-04-07 07:49:42Z hoffmann
97! icloud_scheme replaced by microphysics_*
98!
[1809]99! 1808 2016-04-05 19:44:00Z raasch
100! output message in case unscheduled radiation calls removed
101!
[1798]102! 1797 2016-03-21 16:50:28Z raasch
103! introduction of different datatransfer modes
104!
[1792]105! 1791 2016-03-11 10:41:25Z raasch
106! call of pmci_update_new removed
107!
[1787]108! 1786 2016-03-08 05:49:27Z raasch
109! +module spectrum
110!
[1784]111! 1783 2016-03-06 18:36:17Z raasch
112! switch back of netcdf data format for mask output moved to the mask output
113! routine
114!
[1782]115! 1781 2016-03-03 15:12:23Z raasch
116! some pmc calls removed at the beginning (before timeloop),
117! pmc initialization moved to the main program
118!
[1765]119! 1764 2016-02-28 12:45:19Z raasch
120! PMC_ACTIVE flags removed,
121! bugfix: nest synchronization after first call of timestep
122!
[1763]123! 1762 2016-02-25 12:31:13Z hellstea
124! Introduction of nested domain feature
125!
[1737]126! 1736 2015-12-04 08:56:33Z raasch
127! no perturbations added to total domain if energy limit has been set zero
128!
[1692]129! 1691 2015-10-26 16:17:44Z maronga
130! Added option for spin-ups without land surface and radiation models. Moved calls
131! for radiation and lan surface schemes.
132!
[1683]133! 1682 2015-10-07 23:56:08Z knoop
134! Code annotations made doxygen readable
135!
[1672]136! 1671 2015-09-25 03:29:37Z raasch
137! bugfix: ghostpoint exchange for array diss in case that sgs velocities are used
138! for particles
139!
[1586]140! 1585 2015-04-30 07:05:52Z maronga
141! Moved call of radiation scheme. Added support for RRTM
142!
[1552]143! 1551 2015-03-03 14:18:16Z maronga
144! Added interface for different radiation schemes.
145!
[1497]146! 1496 2014-12-02 17:25:50Z maronga
147! Added calls for the land surface model and radiation scheme
148!
[1403]149! 1402 2014-05-09 14:25:13Z raasch
150! location messages modified
151!
[1385]152! 1384 2014-05-02 14:31:06Z raasch
153! location messages added
154!
[1381]155! 1380 2014-04-28 12:40:45Z heinze
156! CALL of nudge_ref added
157! bc_pt_t_val and bc_q_t_val are updated in case nudging is used
158!
[1366]159! 1365 2014-04-22 15:03:56Z boeske
[1365]160! Reset sums_ls_l to zero at each timestep
161! +sums_ls_l
162! Calculation of reference state (previously in subroutine calc_mean_profile)
163
[1343]164! 1342 2014-03-26 17:04:47Z kanani
165! REAL constants defined as wp-kind
166!
[1321]167! 1320 2014-03-20 08:40:49Z raasch
[1320]168! ONLY-attribute added to USE-statements,
169! kind-parameters added to all INTEGER and REAL declaration statements,
170! kinds are defined in new module kinds,
171! old module precision_kind is removed,
172! revision history before 2012 removed,
173! comment fields (!:) to be used for variable explanations added to
174! all variable declaration statements
[1319]175! 1318 2014-03-17 13:35:16Z raasch
176! module interfaces removed
177!
[1309]178! 1308 2014-03-13 14:58:42Z fricke
179! +netcdf_data_format_save
180! For masked data, parallel netcdf output is not tested so far, hence
181! netcdf_data_format is switched back to non-paralell output.
182!
[1277]183! 1276 2014-01-15 13:40:41Z heinze
184! Use LSF_DATA also in case of Dirichlet bottom boundary condition for scalars
185!
[1258]186! 1257 2013-11-08 15:18:40Z raasch
187! acc-update-host directive for timestep removed
188!
[1242]189! 1241 2013-10-30 11:36:58Z heinze
190! Generalize calc_mean_profile for wider use
191! Determine shf and qsws in dependence on data from LSF_DATA
192! Determine ug and vg in dependence on data from LSF_DATA
[1222]193! 1221 2013-09-10 08:59:13Z raasch
194! host update of arrays before timestep is called
195!
[1182]196! 1179 2013-06-14 05:57:58Z raasch
197! mean profiles for reference state are only calculated if required,
198! small bugfix for background communication
199!
[1172]200! 1171 2013-05-30 11:27:45Z raasch
201! split of prognostic_equations deactivated (comment lines), for the time being
202!
[1132]203! 1128 2013-04-12 06:19:32Z raasch
[1128]204! asynchronous transfer of ghost point data realized for acc-optimized version:
205! prognostic_equations are first called two times for those points required for
206! the left-right and north-south exchange, respectively, and then for the
207! remaining points,
208! those parts requiring global communication moved from prognostic_equations to
209! here
[392]210!
[1116]211! 1115 2013-03-26 18:16:16Z hoffmann
212! calculation of qr and nr is restricted to precipitation
213!
[1114]214! 1113 2013-03-10 02:48:14Z raasch
215! GPU-porting of boundary conditions,
216! openACC directives updated
217! formal parameter removed from routine boundary_conds
218!
[1112]219! 1111 2013-03-08 23:54:10Z raasch
220! +internal timestep counter for cpu statistics added,
221! openACC directives updated
222!
[1093]223! 1092 2013-02-02 11:24:22Z raasch
224! unused variables removed
225!
[1066]226! 1065 2012-11-22 17:42:36Z hoffmann
227! exchange of diss (dissipation rate) in case of turbulence = .TRUE. added
228!
[1054]229! 1053 2012-11-13 17:11:03Z hoffmann
230! exchange of ghost points for nr, qr added
231!
[1037]232! 1036 2012-10-22 13:43:42Z raasch
233! code put under GPL (PALM 3.9)
234!
[1020]235! 1019 2012-09-28 06:46:45Z raasch
236! non-optimized version of prognostic_equations removed
237!
[1017]238! 1015 2012-09-27 09:23:24Z raasch
239! +call of prognostic_equations_acc
240!
[1002]241! 1001 2012-09-13 14:08:46Z raasch
242! all actions concerning leapfrog- and upstream-spline-scheme removed
243!
[850]244! 849 2012-03-15 10:35:09Z raasch
245! advec_particles renamed lpm, first_call_advec_particles renamed first_call_lpm
246!
[826]247! 825 2012-02-19 03:03:44Z raasch
248! wang_collision_kernel renamed wang_kernel
249!
[1]250! Revision 1.1  1997/08/11 06:19:04  raasch
251! Initial revision
252!
253!
254! Description:
255! ------------
[1682]256!> Integration in time of the model equations, statistical analysis and graphic
257!> output
[1]258!------------------------------------------------------------------------------!
[1682]259 SUBROUTINE time_integration
260 
[1]261
[1320]262    USE advec_ws,                                                              &
263        ONLY:  ws_statistics
264
265    USE arrays_3d,                                                             &
[2174]266        ONLY:  diss, dzu, e, e_p, nr, nr_p, prho, pt, pt_p, pt_init, q_init, q,&
267               ql, ql_c, ql_v, ql_vp, qr, qr_p, q_p, ref_state, rho_ocean,     &
268               s, s_p, sa_p, tend, u, u_p, v, vpt, v_p, w, w_p
[1320]269
[1365]270    USE calc_mean_profile_mod,                                                 &
[1320]271        ONLY:  calc_mean_profile
272
273    USE control_parameters,                                                    &
274        ONLY:  advected_distance_x, advected_distance_y, average_count_3d,     &
[1833]275               averaging_interval, averaging_interval_pr,                      &
276               bc_lr_cyc, bc_ns_cyc, bc_pt_t_val,                              &
[1380]277               bc_q_t_val, call_psolver_at_all_substeps, cloud_droplets,       &
[1691]278               cloud_physics, constant_flux_layer, constant_heatflux,          &
279               create_disturbances, dopr_n, constant_diffusion, coupling_mode, &
280               coupling_start_time, current_timestep_number,                   &
281               disturbance_created, disturbance_energy_limit, dist_range,      &
282               do_sum, dt_3d, dt_averaging_input, dt_averaging_input_pr,       &
283               dt_coupling, dt_data_output_av, dt_disturb, dt_do2d_xy,         &
284               dt_do2d_xz, dt_do2d_yz, dt_do3d, dt_domask,dt_dopts, dt_dopr,   &
[1833]285               dt_dopr_listing, dt_dots, dt_dvrp, dt_run_control,              &
[1320]286               end_time, first_call_lpm, galilei_transformation, humidity,     &
[2232]287               intermediate_timestep_count, intermediate_timestep_count_max,   &
288               land_surface, large_scale_forcing,                              &
[1822]289               loop_optimization, lsf_surf, lsf_vert, masks,                   &
290               microphysics_seifert, mid, nest_domain,                         &
[1783]291               neutral, nr_timesteps_this_run, nudging,                        &
[2118]292               ocean, passive_scalar,                                          &
[1320]293               prho_reference, pt_reference, pt_slope_offset, random_heatflux, &
294               run_coupled, simulated_time, simulated_time_chr,                &
295               skip_time_do2d_xy, skip_time_do2d_xz, skip_time_do2d_yz,        &
296               skip_time_do3d, skip_time_domask, skip_time_dopr,               &
[1833]297               skip_time_data_output_av, sloping_surface,                      &
[1320]298               stop_dt, terminate_coupled, terminate_run, timestep_scheme,     &
299               time_coupling, time_do2d_xy, time_do2d_xz, time_do2d_yz,        &
300               time_do3d, time_domask, time_dopr, time_dopr_av,                &
301               time_dopr_listing, time_dopts, time_dosp, time_dosp_av,         &
302               time_dots, time_do_av, time_do_sla, time_disturb, time_dvrp,    &
[1496]303               time_run_control, time_since_reference_point,                   &
[2050]304               turbulent_inflow, turbulent_outflow, urban_surface,             &
[2011]305               use_initial_profile_as_reference,                               &
[1957]306               use_single_reference_value, u_gtrans, v_gtrans, virtual_flight, &
307               ws_scheme_mom, ws_scheme_sca
[1320]308
309    USE cpulog,                                                                &
310        ONLY:  cpu_log, log_point, log_point_s
311
[1957]312    USE flight_mod,                                                            &
313        ONLY:  flight_measurement
314
315
[1320]316    USE indices,                                                               &
[2232]317        ONLY:  nbgp, nx, nxl, nxlg, nxr, nxrg, nyn, nyng, nys, nysg, nzb, nzt
[1320]318
319    USE interaction_droplets_ptq_mod,                                          &
320        ONLY:  interaction_droplets_ptq
321
[1918]322    USE interfaces
323
[1320]324    USE kinds
325
[1496]326    USE land_surface_model_mod,                                                &
[2232]327        ONLY:  lsm_energy_balance, lsm_soil_model,                             &
[1691]328               skip_time_do_lsm
[1496]329
[1320]330    USE ls_forcing_mod,                                                        &
331        ONLY:  ls_forcing_surf, ls_forcing_vert
332
[1849]333    USE microphysics_mod,                                                      &
334        ONLY: collision_turbulence
335
[1365]336    USE nudge_mod,                                                             &
[1380]337        ONLY:  calc_tnudge, nudge_ref
[1365]338
[1320]339    USE particle_attributes,                                                   &
[1671]340        ONLY:  particle_advection, particle_advection_start,                   &
341               use_sgs_for_particles, wang_kernel
[1320]342
[1]343    USE pegrid
344
[1762]345    USE pmc_interface,                                                         &
[1933]346        ONLY:  nested_run, nesting_mode, pmci_datatrans,                       &
347               pmci_ensure_nest_mass_conservation, pmci_synchronize
[1762]348
[1320]349    USE production_e_mod,                                                      &
350        ONLY:  production_e_init
351
[1402]352    USE progress_bar,                                                          &
353        ONLY:  finish_progress_bar, output_progress_bar
354
[1320]355    USE prognostic_equations_mod,                                              &
[2118]356        ONLY:  prognostic_equations_cache, prognostic_equations_vector
[1320]357
[1496]358    USE radiation_model_mod,                                                   &
[1976]359        ONLY: dt_radiation, force_radiation_call, radiation, radiation_control,&
360              skip_time_do_radiation, time_radiation
[1496]361
[1833]362    USE spectra_mod,                                                           &
363        ONLY: average_count_sp, averaging_interval_sp, calc_spectra, dt_dosp,  &
364              skip_time_dosp
[1786]365
[1320]366    USE statistics,                                                            &
[1918]367        ONLY:  flow_statistics_called, hom, pr_palm, sums_ls_l, u_max,         &
368               u_max_ijk, v_max, v_max_ijk, w_max, w_max_ijk
[1320]369
[1691]370    USE surface_layer_fluxes_mod,                                              &
371        ONLY:  surface_layer_fluxes
372
[2232]373    USE surface_mod,                                                           &
374        ONLY:  surf_def_h, surf_lsm_h, surf_usm_h
375
[2007]376    USE urban_surface_mod,                                                     &
[2232]377        ONLY:  usm_material_heat_model, usm_material_model,                    &
[2259]378               usm_radiation, usm_surface_energy_balance
[2007]379
[2259]380    USE synthetic_turbulence_generator_mod,                                    &
381        ONLY:  stg_main, use_synthetic_turbulence_generator
382
[1320]383    USE user_actions_mod,                                                      &
384        ONLY:  user_actions
385
[1914]386    USE wind_turbine_model_mod,                                                &
387        ONLY:  wind_turbine, wtm_forces
388
[1]389    IMPLICIT NONE
390
[1682]391    CHARACTER (LEN=9) ::  time_to_string          !<
[1]392
[1918]393    REAL(wp) ::  dt_3d_old  !< temporary storage of timestep to be used for
394                            !< steering of run control output interval
395
[1]396!
[1918]397!-- At beginning determine the first time step
398    CALL timestep
[667]399
[1764]400!
401!-- Synchronize the timestep in case of nested run.
402    IF ( nested_run )  THEN
[1878]403!
404!--    Synchronization by unifying the time step.
405!--    Global minimum of all time-steps is used for all.
406       CALL pmci_synchronize
[1764]407    ENDIF
408
[1918]409!
410!-- Determine and print out the run control quantities before the first time
411!-- step of this run. For the initial run, some statistics (e.g. divergence)
412!-- need to be determined first.
413    IF ( simulated_time == 0.0_wp )  CALL flow_statistics
[1]414    CALL run_control
415
[108]416!
417!-- Data exchange between coupled models in case that a call has been omitted
418!-- at the end of the previous run of a job chain.
[291]419    IF ( coupling_mode /= 'uncoupled'  .AND.  run_coupled )  THEN
[108]420!
421!--    In case of model termination initiated by the local model the coupler
422!--    must not be called because this would again cause an MPI hang.
[1918]423       DO WHILE ( time_coupling >= dt_coupling  .AND.  terminate_coupled == 0 )
[108]424          CALL surface_coupler
425          time_coupling = time_coupling - dt_coupling
426       ENDDO
[1918]427       IF (time_coupling == 0.0_wp  .AND.                                      &
428           time_since_reference_point < dt_coupling )                          &
[348]429       THEN
430          time_coupling = time_since_reference_point
431       ENDIF
[108]432    ENDIF
433
[1]434#if defined( __dvrp_graphics )
435!
436!-- Time measurement with dvrp software 
437    CALL DVRP_LOG_EVENT( 2, current_timestep_number )
438#endif
439
[1402]440    CALL location_message( 'start with time-stepping', .TRUE. )
[1]441!
442!-- Start of the time loop
443    DO  WHILE ( simulated_time < end_time  .AND.  .NOT. stop_dt  .AND. &
444                .NOT. terminate_run )
445
446       CALL cpu_log( log_point_s(10), 'timesteps', 'start' )
447!
[1241]448!--    Determine ug, vg and w_subs in dependence on data from external file
449!--    LSF_DATA
[1365]450       IF ( large_scale_forcing .AND. lsf_vert )  THEN
[1241]451           CALL ls_forcing_vert ( simulated_time )
[1365]452           sums_ls_l = 0.0_wp
[1241]453       ENDIF
454
455!
[1380]456!--    Set pt_init and q_init to the current profiles taken from
457!--    NUDGING_DATA
458       IF ( nudging )  THEN
459           CALL nudge_ref ( simulated_time )
460!
461!--        Store temperature gradient at the top boundary for possible Neumann
462!--        boundary condition
463           bc_pt_t_val = ( pt_init(nzt+1) - pt_init(nzt) ) / dzu(nzt+1)
464           bc_q_t_val  = ( q_init(nzt+1) - q_init(nzt) ) / dzu(nzt+1)
465       ENDIF
466
467!
[1]468!--    Execute the user-defined actions
469       CALL user_actions( 'before_timestep' )
470
471!
[1914]472!--    Calculate forces by wind turbines
473       IF ( wind_turbine )  THEN
474
475          CALL cpu_log( log_point(55), 'wind_turbine', 'start' )
476
477          CALL wtm_forces
478
479          CALL cpu_log( log_point(55), 'wind_turbine', 'stop' )
480
481       ENDIF       
482       
483!
[1]484!--    Start of intermediate step loop
485       intermediate_timestep_count = 0
486       DO  WHILE ( intermediate_timestep_count < &
487                   intermediate_timestep_count_max )
488
489          intermediate_timestep_count = intermediate_timestep_count + 1
490
491!
492!--       Set the steering factors for the prognostic equations which depend
493!--       on the timestep scheme
494          CALL timestep_scheme_steering
495
496!
[1128]497!--       Calculate those variables needed in the tendency terms which need
498!--       global communication
[1179]499          IF ( .NOT. use_single_reference_value  .AND. &
500               .NOT. use_initial_profile_as_reference )  THEN
501!
502!--          Horizontally averaged profiles to be used as reference state in
503!--          buoyancy terms (WARNING: only the respective last call of
504!--          calc_mean_profile defines the reference state!)
[1365]505             IF ( .NOT. neutral )  THEN
506                CALL calc_mean_profile( pt, 4 )
507                ref_state(:)  = hom(:,1,4,0) ! this is used in the buoyancy term
508             ENDIF
509             IF ( ocean )  THEN
[2031]510                CALL calc_mean_profile( rho_ocean, 64 )
[1365]511                ref_state(:)  = hom(:,1,64,0)
512             ENDIF
513             IF ( humidity )  THEN
514                CALL calc_mean_profile( vpt, 44 )
515                ref_state(:)  = hom(:,1,44,0)
516             ENDIF
517
[1179]518          ENDIF
519
[1128]520          IF ( .NOT. constant_diffusion )  CALL production_e_init
521          IF ( ( ws_scheme_mom .OR. ws_scheme_sca )  .AND.  &
522               intermediate_timestep_count == 1 )  CALL ws_statistics
[1365]523!
524!--       In case of nudging calculate current nudging time scale and horizontal
[1380]525!--       means of u, v, pt and q
[1365]526          IF ( nudging )  THEN
527             CALL calc_tnudge( simulated_time )
528             CALL calc_mean_profile( u, 1 )
529             CALL calc_mean_profile( v, 2 )
530             CALL calc_mean_profile( pt, 4 )
531             CALL calc_mean_profile( q, 41 )
532          ENDIF
[1128]533
534!
[1]535!--       Solve the prognostic equations. A fast cache optimized version with
536!--       only one single loop is used in case of Piascek-Williams advection
537!--       scheme. NEC vector machines use a different version, because
538!--       in the other versions a good vectorization is prohibited due to
539!--       inlining problems.
[1019]540          IF ( loop_optimization == 'cache' )  THEN
541             CALL prognostic_equations_cache
542          ELSEIF ( loop_optimization == 'vector' )  THEN
[63]543             CALL prognostic_equations_vector
[1]544          ENDIF
545
546!
[849]547!--       Particle transport/physics with the Lagrangian particle model
548!--       (only once during intermediate steps, because it uses an Euler-step)
[1128]549!--       ### particle model should be moved before prognostic_equations, in order
550!--       to regard droplet interactions directly
[63]551          IF ( particle_advection  .AND.                         &
552               simulated_time >= particle_advection_start  .AND. &
[1]553               intermediate_timestep_count == 1 )  THEN
[849]554             CALL lpm
555             first_call_lpm = .FALSE.
[1]556          ENDIF
557
558!
559!--       Interaction of droplets with temperature and specific humidity.
560!--       Droplet condensation and evaporation is calculated within
561!--       advec_particles.
562          IF ( cloud_droplets  .AND.  &
563               intermediate_timestep_count == intermediate_timestep_count_max )&
564          THEN
565             CALL interaction_droplets_ptq
566          ENDIF
567
568!
569!--       Exchange of ghost points (lateral boundary conditions)
[2118]570          CALL cpu_log( log_point(26), 'exchange-horiz-progn', 'start' )
[1113]571
[2118]572          CALL exchange_horiz( u_p, nbgp )
573          CALL exchange_horiz( v_p, nbgp )
574          CALL exchange_horiz( w_p, nbgp )
575          CALL exchange_horiz( pt_p, nbgp )
576          IF ( .NOT. constant_diffusion )  CALL exchange_horiz( e_p, nbgp )
577          IF ( ocean )  THEN
578             CALL exchange_horiz( sa_p, nbgp )
579             CALL exchange_horiz( rho_ocean, nbgp )
580             CALL exchange_horiz( prho, nbgp )
581          ENDIF
582          IF ( humidity )  THEN
583             CALL exchange_horiz( q_p, nbgp )
584             IF ( cloud_physics .AND. microphysics_seifert )  THEN
585                CALL exchange_horiz( qr_p, nbgp )
586                CALL exchange_horiz( nr_p, nbgp )
[1053]587             ENDIF
[2118]588          ENDIF
589          IF ( cloud_droplets )  THEN
590             CALL exchange_horiz( ql, nbgp )
591             CALL exchange_horiz( ql_c, nbgp )
592             CALL exchange_horiz( ql_v, nbgp )
593             CALL exchange_horiz( ql_vp, nbgp )
594          ENDIF
595          IF ( wang_kernel  .OR.  collision_turbulence  .OR.                &
596               use_sgs_for_particles )  THEN
597             CALL exchange_horiz( diss, nbgp )
598          ENDIF
599          IF ( passive_scalar )  CALL exchange_horiz( s_p, nbgp )
[1]600
[2118]601          CALL cpu_log( log_point(26), 'exchange-horiz-progn', 'stop' )
[1128]602
[1]603!
604!--       Boundary conditions for the prognostic quantities (except of the
605!--       velocities at the outflow in case of a non-cyclic lateral wall)
[1113]606          CALL boundary_conds
[1]607!
[73]608!--       Swap the time levels in preparation for the next time step.
609          CALL swap_timelevel
610
[1764]611          IF ( nested_run )  THEN
[1797]612
[1764]613             CALL cpu_log( log_point(60), 'nesting', 'start' )
[1762]614!
[1933]615!--          Domain nesting. The data transfer subroutines pmci_parent_datatrans
616!--          and pmci_child_datatrans are called inside the wrapper
[1797]617!--          subroutine pmci_datatrans according to the control parameters
618!--          nesting_mode and nesting_datatransfer_mode.
619!--          TO_DO: why is nesting_mode given as a parameter here?
620             CALL pmci_datatrans( nesting_mode )
[1762]621
[2174]622             IF ( TRIM( nesting_mode ) == 'two-way' .OR.                       &
[1933]623                  nesting_mode == 'vertical' )  THEN
[1762]624!
[1933]625!--             Exchange_horiz is needed for all parent-domains after the
[1764]626!--             anterpolation
627                CALL exchange_horiz( u, nbgp )
628                CALL exchange_horiz( v, nbgp )
629                CALL exchange_horiz( w, nbgp )
[2174]630                IF ( .NOT. neutral )  CALL exchange_horiz( pt, nbgp )
631
632                IF ( humidity )  THEN
633
634                   CALL exchange_horiz( q, nbgp )
635
636                   IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
637!                        CALL exchange_horiz( qc, nbgp )
638                       CALL exchange_horiz( qr, nbgp )
639!                        CALL exchange_horiz( nc, nbgp )
640                       CALL exchange_horiz( nr, nbgp )
641                   ENDIF
642
643                ENDIF
644
645                IF ( passive_scalar )  CALL exchange_horiz( s, nbgp )
646                IF ( .NOT. constant_diffusion )  CALL exchange_horiz( e, nbgp )
[1762]647             ENDIF
648!
[1764]649!--          Correct the w top-BC in nest domains to ensure mass conservation.
[1933]650!--          This action must never be done for the root domain. Vertical
651!--          nesting implies mass conservation.
[1764]652             IF ( nest_domain )  THEN
653                CALL pmci_ensure_nest_mass_conservation
654             ENDIF
655
656             CALL cpu_log( log_point(60), 'nesting', 'stop' )
657
[1762]658          ENDIF
659
660!
[1]661!--       Temperature offset must be imposed at cyclic boundaries in x-direction
662!--       when a sloping surface is used
663          IF ( sloping_surface )  THEN
[707]664             IF ( nxl ==  0 )  pt(:,:,nxlg:nxl-1) = pt(:,:,nxlg:nxl-1) - &
665                                                    pt_slope_offset
666             IF ( nxr == nx )  pt(:,:,nxr+1:nxrg) = pt(:,:,nxr+1:nxrg) + &
667                                                    pt_slope_offset
[1]668          ENDIF
669
670!
[151]671!--       Impose a turbulent inflow using the recycling method
672          IF ( turbulent_inflow )  CALL  inflow_turbulence
673
674!
[2259]675!--       Impose a turbulent inflow using synthetic generated turbulence
676         IF ( use_synthetic_turbulence_generator ) THEN
677            CALL  stg_main
678         ENDIF
679
680!
[2050]681!--       Set values at outflow boundary using the special outflow condition
682          IF ( turbulent_outflow )  CALL  outflow_turbulence
683
684!
[1]685!--       Impose a random perturbation on the horizontal velocity field
[106]686          IF ( create_disturbances  .AND.                                      &
687               ( call_psolver_at_all_substeps  .AND.                           &
[1]688               intermediate_timestep_count == intermediate_timestep_count_max )&
[106]689          .OR. ( .NOT. call_psolver_at_all_substeps  .AND.                     &
690               intermediate_timestep_count == 1 ) )                            &
[1]691          THEN
692             time_disturb = time_disturb + dt_3d
693             IF ( time_disturb >= dt_disturb )  THEN
[1736]694                IF ( disturbance_energy_limit /= 0.0_wp  .AND.                 &
695                     hom(nzb+5,1,pr_palm,0) < disturbance_energy_limit )  THEN
[2232]696                   CALL disturb_field( 'u', tend, u )
697                   CALL disturb_field( 'v', tend, v )
[2178]698                ELSEIF ( ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )            &
699                     .AND. .NOT. nest_domain )  THEN
[1]700!
701!--                Runs with a non-cyclic lateral wall need perturbations
702!--                near the inflow throughout the whole simulation
703                   dist_range = 1
[2232]704                   CALL disturb_field( 'u', tend, u )
705                   CALL disturb_field( 'v', tend, v )
[1]706                   dist_range = 0
707                ENDIF
708                time_disturb = time_disturb - dt_disturb
709             ENDIF
710          ENDIF
711
712!
713!--       Reduce the velocity divergence via the equation for perturbation
714!--       pressure.
[106]715          IF ( intermediate_timestep_count == 1  .OR. &
716                call_psolver_at_all_substeps )  THEN
[1]717             CALL pres
718          ENDIF
719
720!
721!--       If required, compute liquid water content
[1015]722          IF ( cloud_physics )  THEN
723             CALL calc_liquid_water_content
724          ENDIF
[2174]725!
[1115]726!--       If required, compute virtual potential temperature
727          IF ( humidity )  THEN
728             CALL compute_vpt 
729          ENDIF 
[1585]730
[1]731!
732!--       Compute the diffusion quantities
733          IF ( .NOT. constant_diffusion )  THEN
734
735!
[1276]736!--          Determine surface fluxes shf and qsws and surface values
737!--          pt_surface and q_surface in dependence on data from external
738!--          file LSF_DATA respectively
739             IF ( ( large_scale_forcing .AND. lsf_surf ) .AND. &
740                 intermediate_timestep_count == intermediate_timestep_count_max )&
741             THEN
742                CALL ls_forcing_surf ( simulated_time )
743             ENDIF
744
745!
[2232]746!--          First the vertical (and horizontal) fluxes in the surface
747!--          (constant flux) layer are computed
[1691]748             IF ( constant_flux_layer )  THEN
749                CALL cpu_log( log_point(19), 'surface_layer_fluxes', 'start' )
750                CALL surface_layer_fluxes
751                CALL cpu_log( log_point(19), 'surface_layer_fluxes', 'stop' )
[1]752             ENDIF
[1241]753
[1]754!
[1691]755!--          If required, solve the energy balance for the surface and run soil
[2232]756!--          model. Call for horizontal as well as vertical surfaces
[1691]757             IF ( land_surface .AND. simulated_time > skip_time_do_lsm)  THEN
758
759                CALL cpu_log( log_point(54), 'land_surface', 'start' )
[2232]760!
761!--             Call for horizontal upward-facing surfaces
762                CALL lsm_energy_balance( .TRUE., -1 )
763                CALL lsm_soil_model( .TRUE., -1 )
764!
765!--             Call for northward-facing surfaces
766                CALL lsm_energy_balance( .FALSE., 0 )
767                CALL lsm_soil_model( .FALSE., 0 )
768!
769!--             Call for southward-facing surfaces
770                CALL lsm_energy_balance( .FALSE., 1 )
771                CALL lsm_soil_model( .FALSE., 1 )
772!
773!--             Call for eastward-facing surfaces
774                CALL lsm_energy_balance( .FALSE., 2 )
775                CALL lsm_soil_model( .FALSE., 2 )
776!
777!--             Call for westward-facing surfaces
778                CALL lsm_energy_balance( .FALSE., 3 )
779                CALL lsm_soil_model( .FALSE., 3 )
780
[1691]781                CALL cpu_log( log_point(54), 'land_surface', 'stop' )
782             ENDIF
[2007]783
[1691]784!
[2007]785!--          If required, solve the energy balance for urban surfaces and run
786!--          the material heat model
787             IF (urban_surface) THEN
788                CALL cpu_log( log_point(74), 'urban_surface', 'start' )
789                CALL usm_surface_energy_balance
790                IF ( usm_material_model )  THEN
791                   CALL usm_material_heat_model
792                ENDIF
793                CALL cpu_log( log_point(74), 'urban_surface', 'stop' )
794             ENDIF
795
796!
[1]797!--          Compute the diffusion coefficients
798             CALL cpu_log( log_point(17), 'diffusivities', 'start' )
[75]799             IF ( .NOT. humidity ) THEN
[97]800                IF ( ocean )  THEN
[388]801                   CALL diffusivities( prho, prho_reference )
[97]802                ELSE
803                   CALL diffusivities( pt, pt_reference )
804                ENDIF
[1]805             ELSE
[97]806                CALL diffusivities( vpt, pt_reference )
[1]807             ENDIF
[2232]808
[1]809             CALL cpu_log( log_point(17), 'diffusivities', 'stop' )
810
811          ENDIF
812
[1691]813!
814!--       If required, calculate radiative fluxes and heating rates
815          IF ( radiation .AND. intermediate_timestep_count                     &
816               == intermediate_timestep_count_max .AND. simulated_time >    &
817               skip_time_do_radiation )  THEN
818
819               time_radiation = time_radiation + dt_3d
820
821             IF ( time_radiation >= dt_radiation .OR. force_radiation_call )   &
822             THEN
823
824                CALL cpu_log( log_point(50), 'radiation', 'start' )
825
826                IF ( .NOT. force_radiation_call )  THEN
827                   time_radiation = time_radiation - dt_radiation
828                ENDIF
829
[1976]830                CALL radiation_control
[1691]831
832                CALL cpu_log( log_point(50), 'radiation', 'stop' )
[2007]833
834                IF (urban_surface)  THEN
835                   CALL cpu_log( log_point(75), 'usm_radiation', 'start' )
836                   CALL usm_radiation
837                   CALL cpu_log( log_point(75), 'usm_radiation', 'stop' )
838                ENDIF
839
[1691]840             ENDIF
841          ENDIF
842
[1]843       ENDDO   ! Intermediate step loop
844
845!
846!--    Increase simulation time and output times
[1111]847       nr_timesteps_this_run      = nr_timesteps_this_run + 1
[291]848       current_timestep_number    = current_timestep_number + 1
849       simulated_time             = simulated_time   + dt_3d
850       simulated_time_chr         = time_to_string( simulated_time )
851       time_since_reference_point = simulated_time - coupling_start_time
852
[1957]853
854
[1]855       IF ( simulated_time >= skip_time_data_output_av )  THEN
856          time_do_av         = time_do_av       + dt_3d
857       ENDIF
858       IF ( simulated_time >= skip_time_do2d_xy )  THEN
859          time_do2d_xy       = time_do2d_xy     + dt_3d
860       ENDIF
861       IF ( simulated_time >= skip_time_do2d_xz )  THEN
862          time_do2d_xz       = time_do2d_xz     + dt_3d
863       ENDIF
864       IF ( simulated_time >= skip_time_do2d_yz )  THEN
865          time_do2d_yz       = time_do2d_yz     + dt_3d
866       ENDIF
867       IF ( simulated_time >= skip_time_do3d    )  THEN
868          time_do3d          = time_do3d        + dt_3d
869       ENDIF
[410]870       DO  mid = 1, masks
871          IF ( simulated_time >= skip_time_domask(mid) )  THEN
872             time_domask(mid)= time_domask(mid) + dt_3d
873          ENDIF
874       ENDDO
[1]875       time_dvrp          = time_dvrp        + dt_3d
876       IF ( simulated_time >= skip_time_dosp )  THEN
877          time_dosp       = time_dosp        + dt_3d
878       ENDIF
879       time_dots          = time_dots        + dt_3d
[849]880       IF ( .NOT. first_call_lpm )  THEN
[1]881          time_dopts      = time_dopts       + dt_3d
882       ENDIF
883       IF ( simulated_time >= skip_time_dopr )  THEN
884          time_dopr       = time_dopr        + dt_3d
885       ENDIF
886       time_dopr_listing          = time_dopr_listing        + dt_3d
887       time_run_control   = time_run_control + dt_3d
888
889!
[102]890!--    Data exchange between coupled models
[291]891       IF ( coupling_mode /= 'uncoupled'  .AND.  run_coupled )  THEN
[102]892          time_coupling = time_coupling + dt_3d
[343]893
[108]894!
895!--       In case of model termination initiated by the local model
896!--       (terminate_coupled > 0), the coupler must be skipped because it would
897!--       cause an MPI intercomminucation hang.
898!--       If necessary, the coupler will be called at the beginning of the
899!--       next restart run.
900          DO WHILE ( time_coupling >= dt_coupling .AND. terminate_coupled == 0 )
[102]901             CALL surface_coupler
902             time_coupling = time_coupling - dt_coupling
903          ENDDO
904       ENDIF
905
906!
[46]907!--    Execute user-defined actions
908       CALL user_actions( 'after_integration' )
909
910!
[1]911!--    If Galilei transformation is used, determine the distance that the
912!--    model has moved so far
913       IF ( galilei_transformation )  THEN
914          advected_distance_x = advected_distance_x + u_gtrans * dt_3d
915          advected_distance_y = advected_distance_y + v_gtrans * dt_3d
916       ENDIF
917
918!
919!--    Check, if restart is necessary (because cpu-time is expiring or
920!--    because it is forced by user) and set stop flag
[108]921!--    This call is skipped if the remote model has already initiated a restart.
922       IF ( .NOT. terminate_run )  CALL check_for_restart
[1]923
924!
925!--    Carry out statistical analysis and output at the requested output times.
926!--    The MOD function is used for calculating the output time counters (like
927!--    time_dopr) in order to regard a possible decrease of the output time
928!--    interval in case of restart runs
929
930!
931!--    Set a flag indicating that so far no statistics have been created
932!--    for this time step
933       flow_statistics_called = .FALSE.
934
935!
936!--    If required, call flow_statistics for averaging in time
[1342]937       IF ( averaging_interval_pr /= 0.0_wp  .AND.  &
[1]938            ( dt_dopr - time_dopr ) <= averaging_interval_pr  .AND.  &
939            simulated_time >= skip_time_dopr )  THEN
940          time_dopr_av = time_dopr_av + dt_3d
941          IF ( time_dopr_av >= dt_averaging_input_pr )  THEN
942             do_sum = .TRUE.
943             time_dopr_av = MOD( time_dopr_av, &
944                                    MAX( dt_averaging_input_pr, dt_3d ) )
945          ENDIF
946       ENDIF
947       IF ( do_sum )  CALL flow_statistics
948
949!
[410]950!--    Sum-up 3d-arrays for later output of time-averaged 2d/3d/masked data
[1342]951       IF ( averaging_interval /= 0.0_wp  .AND.                                &
[1]952            ( dt_data_output_av - time_do_av ) <= averaging_interval  .AND. &
953            simulated_time >= skip_time_data_output_av )                    &
954       THEN
955          time_do_sla = time_do_sla + dt_3d
956          IF ( time_do_sla >= dt_averaging_input )  THEN
957             CALL sum_up_3d_data
958             average_count_3d = average_count_3d + 1
959             time_do_sla = MOD( time_do_sla, MAX( dt_averaging_input, dt_3d ) )
960          ENDIF
961       ENDIF
962
963!
964!--    Calculate spectra for time averaging
[1342]965       IF ( averaging_interval_sp /= 0.0_wp  .AND.  &
[1]966            ( dt_dosp - time_dosp ) <= averaging_interval_sp  .AND.  &
967            simulated_time >= skip_time_dosp )  THEN
968          time_dosp_av = time_dosp_av + dt_3d
969          IF ( time_dosp_av >= dt_averaging_input_pr )  THEN
970             CALL calc_spectra
971             time_dosp_av = MOD( time_dosp_av, &
972                                 MAX( dt_averaging_input_pr, dt_3d ) )
973          ENDIF
974       ENDIF
975
976!
[1957]977!--    Call flight module and output data
978       IF ( virtual_flight )  THEN
979          CALL flight_measurement
980          CALL data_output_flight
981       ENDIF
982
983!
[1]984!--    Profile output (ASCII) on file
985       IF ( time_dopr_listing >= dt_dopr_listing )  THEN
986          CALL print_1d
987          time_dopr_listing = MOD( time_dopr_listing, MAX( dt_dopr_listing, &
988                                                           dt_3d ) )
989       ENDIF
990
991!
992!--    Graphic output for PROFIL
993       IF ( time_dopr >= dt_dopr )  THEN
994          IF ( dopr_n /= 0 )  CALL data_output_profiles
995          time_dopr = MOD( time_dopr, MAX( dt_dopr, dt_3d ) )
[1342]996          time_dopr_av = 0.0_wp    ! due to averaging (see above)
[1]997       ENDIF
998
999!
1000!--    Graphic output for time series
1001       IF ( time_dots >= dt_dots )  THEN
[48]1002          CALL data_output_tseries
[1]1003          time_dots = MOD( time_dots, MAX( dt_dots, dt_3d ) )
1004       ENDIF
1005
1006!
1007!--    Output of spectra (formatted for use with PROFIL), in case of no
1008!--    time averaging, spectra has to be calculated before
1009       IF ( time_dosp >= dt_dosp )  THEN
1010          IF ( average_count_sp == 0 )  CALL calc_spectra
1011          CALL data_output_spectra
1012          time_dosp = MOD( time_dosp, MAX( dt_dosp, dt_3d ) )
1013       ENDIF
1014
1015!
1016!--    2d-data output (cross-sections)
1017       IF ( time_do2d_xy >= dt_do2d_xy )  THEN
1018          CALL data_output_2d( 'xy', 0 )
1019          time_do2d_xy = MOD( time_do2d_xy, MAX( dt_do2d_xy, dt_3d ) )
1020       ENDIF
1021       IF ( time_do2d_xz >= dt_do2d_xz )  THEN
1022          CALL data_output_2d( 'xz', 0 )
1023          time_do2d_xz = MOD( time_do2d_xz, MAX( dt_do2d_xz, dt_3d ) )
1024       ENDIF
1025       IF ( time_do2d_yz >= dt_do2d_yz )  THEN
1026          CALL data_output_2d( 'yz', 0 )
1027          time_do2d_yz = MOD( time_do2d_yz, MAX( dt_do2d_yz, dt_3d ) )
1028       ENDIF
1029
1030!
1031!--    3d-data output (volume data)
1032       IF ( time_do3d >= dt_do3d )  THEN
1033          CALL data_output_3d( 0 )
1034          time_do3d = MOD( time_do3d, MAX( dt_do3d, dt_3d ) )
1035       ENDIF
1036
1037!
[1783]1038!--    Masked data output
[410]1039       DO  mid = 1, masks
1040          IF ( time_domask(mid) >= dt_domask(mid) )  THEN
1041             CALL data_output_mask( 0 )
1042             time_domask(mid) = MOD( time_domask(mid),  &
1043                                     MAX( dt_domask(mid), dt_3d ) )
1044          ENDIF
1045       ENDDO
1046
1047!
1048!--    Output of time-averaged 2d/3d/masked data
[1]1049       IF ( time_do_av >= dt_data_output_av )  THEN
1050          CALL average_3d_data
1051          CALL data_output_2d( 'xy', 1 )
1052          CALL data_output_2d( 'xz', 1 )
1053          CALL data_output_2d( 'yz', 1 )
1054          CALL data_output_3d( 1 )
[410]1055          DO  mid = 1, masks
1056             CALL data_output_mask( 1 )
1057          ENDDO
[1]1058          time_do_av = MOD( time_do_av, MAX( dt_data_output_av, dt_3d ) )
1059       ENDIF
1060
1061!
1062!--    Output of particle time series
[253]1063       IF ( particle_advection )  THEN
1064          IF ( time_dopts >= dt_dopts  .OR. &
1065               ( simulated_time >= particle_advection_start  .AND. &
[849]1066                 first_call_lpm ) )  THEN
[253]1067             CALL data_output_ptseries
1068             time_dopts = MOD( time_dopts, MAX( dt_dopts, dt_3d ) )
1069          ENDIF
[1]1070       ENDIF
1071
1072!
1073!--    Output of dvrp-graphics (isosurface, particles, slicer)
1074#if defined( __dvrp_graphics )
1075       CALL DVRP_LOG_EVENT( -2, current_timestep_number-1 )
1076#endif
1077       IF ( time_dvrp >= dt_dvrp )  THEN
1078          CALL data_output_dvrp
1079          time_dvrp = MOD( time_dvrp, MAX( dt_dvrp, dt_3d ) )
1080       ENDIF
1081#if defined( __dvrp_graphics )
1082       CALL DVRP_LOG_EVENT( 2, current_timestep_number )
1083#endif
1084
1085!
1086!--    If required, set the heat flux for the next time step at a random value
[2232]1087       IF ( constant_heatflux  .AND.  random_heatflux )  THEN
1088          IF ( surf_def_h(0)%ns >= 1 )  CALL disturb_heatflux( surf_def_h(0) )
1089          IF ( surf_lsm_h%ns    >= 1 )  CALL disturb_heatflux( surf_lsm_h    )
1090          IF ( surf_usm_h%ns    >= 1 )  CALL disturb_heatflux( surf_usm_h    )
1091       ENDIF
[1]1092
1093!
1094!--    Execute user-defined actions
1095       CALL user_actions( 'after_timestep' )
1096
[1402]1097!
[1918]1098!--    Determine size of next time step. Save timestep dt_3d because it is
1099!--    newly calculated in routine timestep, but required further below for
1100!--    steering the run control output interval
1101       dt_3d_old = dt_3d
1102       CALL timestep
1103
1104!
[1925]1105!--    Synchronize the timestep in case of nested run.
1106       IF ( nested_run )  THEN
1107!
1108!--       Synchronize by unifying the time step.
1109!--       Global minimum of all time-steps is used for all.
1110          CALL pmci_synchronize
1111       ENDIF
1112
1113!
[1918]1114!--    Computation and output of run control parameters.
1115!--    This is also done whenever perturbations have been imposed
1116       IF ( time_run_control >= dt_run_control  .OR.                     &
1117            timestep_scheme(1:5) /= 'runge'  .OR.  disturbance_created ) &
1118       THEN
1119          CALL run_control
1120          IF ( time_run_control >= dt_run_control )  THEN
1121             time_run_control = MOD( time_run_control, &
1122                                     MAX( dt_run_control, dt_3d_old ) )
1123          ENDIF
1124       ENDIF
1125
1126!
[1402]1127!--    Output elapsed simulated time in form of a progress bar on stdout
1128       IF ( myid == 0 )  CALL output_progress_bar
1129
[1]1130       CALL cpu_log( log_point_s(10), 'timesteps', 'stop' )
1131
[667]1132
[1]1133    ENDDO   ! time loop
1134
[1402]1135    IF ( myid == 0 )  CALL finish_progress_bar
1136
[1]1137#if defined( __dvrp_graphics )
1138    CALL DVRP_LOG_EVENT( -2, current_timestep_number )
1139#endif
1140
[1402]1141    CALL location_message( 'finished time-stepping', .TRUE. )
[1384]1142
[1]1143 END SUBROUTINE time_integration
Note: See TracBrowser for help on using the repository browser.