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

Last change on this file since 2297 was 2292, checked in by schwenkel, 7 years ago

implementation of new bulk microphysics scheme

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