source: palm/trunk/SOURCE/init_3d_model.f90 @ 120

Last change on this file since 120 was 110, checked in by raasch, 16 years ago

New:
---
Allows runs for a coupled atmosphere-ocean LES,
coupling frequency is controlled by new d3par-parameter dt_coupling,
the coupling mode (atmosphere_to_ocean or ocean_to_atmosphere) for the
respective processes is read from environment variable coupling_mode,
which is set by the mpiexec-command,
communication between the two models is done using the intercommunicator
comm_inter,
local files opened by the ocean model get the additional suffic "_O".
Assume saturation at k=nzb_s_inner(j,i) for atmosphere coupled to ocean.

A momentum flux can be set as top boundary condition using the new
inipar parameter top_momentumflux_u|v.

Non-cyclic boundary conditions can be used along all horizontal directions.

Quantities w*p* and w"e can be output as vertical profiles.

Initial profiles are reset to constant profiles in case that initializing_actions /= 'set_constant_profiles'. (init_rankine)

Optionally calculate km and kh from initial TKE e_init.

Changed:


Remaining variables iran changed to iran_part (advec_particles, init_particles).

In case that the presure solver is not called for every Runge-Kutta substep
(call_psolver_at_all_substeps = .F.), it is called after the first substep
instead of the last. In that case, random perturbations are also added to the
velocity field after the first substep.

Initialization of km,kh = 0.00001 for ocean = .T. (for ocean = .F. it remains 0.01).

Allow data_output_pr= q, wq, w"q", w*q* for humidity = .T. (instead of cloud_physics = .T.).

Errors:


Bugs from code parts for non-cyclic boundary conditions are removed: loops for
u and v are starting from index nxlu, nysv, respectively. The radiation boundary
condition is used for every Runge-Kutta substep. Velocity phase speeds for
the radiation boundary conditions are calculated for the first Runge-Kutta
substep only and reused for the further substeps. New arrays c_u, c_v, and c_w
are defined for this purpose. Several index errors are removed from the
radiation boundary condition code parts. Upper bounds for calculating
u_0 and v_0 (in production_e) are nxr+1 and nyn+1 because otherwise these
values are not available in case of non-cyclic boundary conditions.

+dots_num_palm in module user, +module netcdf_control in user_init (both in user_interface)

Bugfix: wrong sign removed from the buoyancy production term in the case use_reference = .T. (production_e)

Bugfix: Error message concerning output of particle concentration (pc) modified (check_parameters).

Bugfix: Rayleigh damping for ocean fixed.

  • Property svn:keywords set to Id
File size: 37.4 KB
Line 
1#if defined( __ibmy_special )
2@PROCESS NOOPTimize
3#endif
4 SUBROUTINE init_3d_model
5
6!------------------------------------------------------------------------------!
7! Actual revisions:
8! -----------------
9!
10!
11! Former revisions:
12! -----------------
13! $Id: init_3d_model.f90 110 2007-10-05 05:13:14Z raasch $
14!
15! 108 2007-08-24 15:10:38Z letzel
16! Flux initialization in case of coupled runs, +momentum fluxes at top boundary,
17! +arrays for phase speed c_u, c_v, c_w, indices for u|v|w_m_l|r changed
18! +qswst_remote in case of atmosphere model with humidity coupled to ocean
19! Rayleigh damping for ocean, optionally calculate km and kh from initial
20! TKE e_init
21!
22! 97 2007-06-21 08:23:15Z raasch
23! Initialization of salinity, call of init_ocean
24!
25! 87 2007-05-22 15:46:47Z raasch
26! var_hom and var_sum renamed pr_palm
27!
28! 75 2007-03-22 09:54:05Z raasch
29! Arrays for radiation boundary conditions are allocated (u_m_l, u_m_r, etc.),
30! bugfix for cases with the outflow damping layer extending over more than one
31! subdomain, moisture renamed humidity,
32! new initializing action "by_user" calls user_init_3d_model,
33! precipitation_amount/rate, ts_value are allocated, +module netcdf_control,
34! initial velocities at nzb+1 are regarded for volume
35! flow control in case they have been set zero before (to avoid small timesteps)
36! -uvmean_outflow, uxrp, vynp eliminated
37!
38! 19 2007-02-23 04:53:48Z raasch
39! +handling of top fluxes
40!
41! RCS Log replace by Id keyword, revision history cleaned up
42!
43! Revision 1.49  2006/08/22 15:59:07  raasch
44! No optimization of this file on the ibmy (Yonsei Univ.)
45!
46! Revision 1.1  1998/03/09 16:22:22  raasch
47! Initial revision
48!
49!
50! Description:
51! ------------
52! Allocation of arrays and initialization of the 3D model via
53! a) pre-run the 1D model
54! or
55! b) pre-set constant linear profiles
56! or
57! c) read values of a previous run
58!------------------------------------------------------------------------------!
59
60    USE arrays_3d
61    USE averaging
62    USE cloud_parameters
63    USE constants
64    USE control_parameters
65    USE cpulog
66    USE indices
67    USE interfaces
68    USE model_1d
69    USE netcdf_control
70    USE particle_attributes
71    USE pegrid
72    USE profil_parameter
73    USE random_function_mod
74    USE statistics
75
76    IMPLICIT NONE
77
78    INTEGER ::  i, j, k, sr
79
80    INTEGER, DIMENSION(:), ALLOCATABLE ::  ngp_2dh_l, ngp_3d_inner_l
81
82    INTEGER, DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_outer_l
83
84    REAL, DIMENSION(1:2) ::  volume_flow_area_l, volume_flow_initial_l
85
86
87!
88!-- Allocate arrays
89    ALLOCATE( ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions), &
90              ngp_3d(0:statistic_regions),                                  &
91              ngp_3d_inner(0:statistic_regions),                            &
92              ngp_3d_inner_l(0:statistic_regions),                          &
93              sums_divnew_l(0:statistic_regions),                           &
94              sums_divold_l(0:statistic_regions) )
95    ALLOCATE( rdf(nzb+1:nzt) )
96    ALLOCATE( hom_sum(nzb:nzt+1,pr_palm+max_pr_user,0:statistic_regions),   &
97              ngp_2dh_outer(nzb:nzt+1,0:statistic_regions),                 &
98              ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions),               &
99              rmask(nys-1:nyn+1,nxl-1:nxr+1,0:statistic_regions),           &
100              sums(nzb:nzt+1,pr_palm+max_pr_user),                          &
101              sums_l(nzb:nzt+1,pr_palm+max_pr_user,0:threads_per_task-1),   &
102              sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1), &
103              sums_up_fraction_l(10,3,0:statistic_regions),                 &
104              sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions),                &
105              ts_value(var_ts,0:statistic_regions) )
106    ALLOCATE( km_damp_x(nxl-1:nxr+1), km_damp_y(nys-1:nyn+1) )
107
108    ALLOCATE( rif_1(nys-1:nyn+1,nxl-1:nxr+1), shf_1(nys-1:nyn+1,nxl-1:nxr+1), &
109              ts(nys-1:nyn+1,nxl-1:nxr+1), tswst_1(nys-1:nyn+1,nxl-1:nxr+1),  &
110              us(nys-1:nyn+1,nxl-1:nxr+1), usws_1(nys-1:nyn+1,nxl-1:nxr+1),   &
111              uswst_1(nys-1:nyn+1,nxl-1:nxr+1),                               &
112              vsws_1(nys-1:nyn+1,nxl-1:nxr+1),                                &
113              vswst_1(nys-1:nyn+1,nxl-1:nxr+1), z0(nys-1:nyn+1,nxl-1:nxr+1) )
114
115    IF ( timestep_scheme(1:5) /= 'runge' )  THEN
116!
117!--    Leapfrog scheme needs two timelevels of diffusion quantities
118       ALLOCATE( rif_2(nys-1:nyn+1,nxl-1:nxr+1),   &
119                 shf_2(nys-1:nyn+1,nxl-1:nxr+1),   &
120                 tswst_2(nys-1:nyn+1,nxl-1:nxr+1), &
121                 usws_2(nys-1:nyn+1,nxl-1:nxr+1),  &
122                 uswst_2(nys-1:nyn+1,nxl-1:nxr+1), &
123                 vswst_2(nys-1:nyn+1,nxl-1:nxr+1), &
124                 vsws_2(nys-1:nyn+1,nxl-1:nxr+1) )
125    ENDIF
126
127    ALLOCATE( d(nzb+1:nzta,nys:nyna,nxl:nxra),         &
128              e_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
129              e_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
130              e_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
131              kh_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
132              km_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
133              p(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),    &
134              pt_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
135              pt_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
136              pt_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
137              tend(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
138              u_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
139              u_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
140              u_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
141              v_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
142              v_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
143              v_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
144              w_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
145              w_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
146              w_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
147
148    IF ( timestep_scheme(1:5) /= 'runge' )  THEN
149       ALLOCATE( kh_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
150                 km_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
151    ENDIF
152
153    IF ( humidity  .OR.  passive_scalar ) THEN
154!
155!--    2D-humidity/scalar arrays
156       ALLOCATE ( qs(nys-1:nyn+1,nxl-1:nxr+1),     &
157                  qsws_1(nys-1:nyn+1,nxl-1:nxr+1), &
158                  qswst_1(nys-1:nyn+1,nxl-1:nxr+1) )
159
160       IF ( timestep_scheme(1:5) /= 'runge' )  THEN
161          ALLOCATE( qsws_2(nys-1:nyn+1,nxl-1:nxr+1), &
162                    qswst_2(nys-1:nyn+1,nxl-1:nxr+1) )
163       ENDIF
164!
165!--    3D-humidity/scalar arrays
166       ALLOCATE( q_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
167                 q_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
168                 q_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
169
170!
171!--    3D-arrays needed for humidity only
172       IF ( humidity )  THEN
173          ALLOCATE( vpt_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
174
175          IF ( timestep_scheme(1:5) /= 'runge' )  THEN
176             ALLOCATE( vpt_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
177          ENDIF
178
179          IF ( cloud_physics ) THEN
180!
181!--          Liquid water content
182             ALLOCATE ( ql_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
183!
184!--          Precipitation amount and rate (only needed if output is switched)
185             ALLOCATE( precipitation_amount(nys-1:nyn+1,nxl-1:nxr+1), &
186                       precipitation_rate(nys-1:nyn+1,nxl-1:nxr+1) )
187          ENDIF
188
189          IF ( cloud_droplets )  THEN
190!
191!--          Liquid water content, change in liquid water content,
192!--          real volume of particles (with weighting), volume of particles
193             ALLOCATE ( ql_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
194                        ql_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
195                        ql_v(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
196                        ql_vp(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
197          ENDIF
198
199       ENDIF
200
201    ENDIF
202
203    IF ( ocean )  THEN
204       ALLOCATE( saswsb_1(nys-1:nyn+1,nxl-1:nxr+1), &
205                 saswst_1(nys-1:nyn+1,nxl-1:nxr+1) )
206       ALLOCATE( rho_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
207                 sa_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
208                 sa_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
209                 sa_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
210       rho => rho_1  ! routine calc_mean_profile requires density to be a
211                     ! pointer
212       IF ( humidity_remote )  THEN
213          ALLOCATE( qswst_remote(nys-1:nyn+1,nxl-1:nxr+1) )
214          qswst_remote = 0.0
215       ENDIF
216    ENDIF
217
218!
219!-- 3D-array for storing the dissipation, needed for calculating the sgs
220!-- particle velocities
221    IF ( use_sgs_for_particles )  THEN
222       ALLOCATE ( diss(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
223    ENDIF
224
225    IF ( dt_dosp /= 9999999.9 )  THEN
226       ALLOCATE( spectrum_x( 1:nx/2, 1:10, 1:10 ), &
227                 spectrum_y( 1:ny/2, 1:10, 1:10 ) )
228    ENDIF
229
230!
231!-- 4D-array for storing the Rif-values at vertical walls
232    IF ( topography /= 'flat' )  THEN
233       ALLOCATE( rif_wall(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1,1:4) )
234       rif_wall = 0.0
235    ENDIF
236
237!
238!-- Velocities at nzb+1 needed for volume flow control
239    IF ( conserve_volume_flow )  THEN
240       ALLOCATE( u_nzb_p1_for_vfc(nys:nyn), v_nzb_p1_for_vfc(nxl:nxr) )
241       u_nzb_p1_for_vfc = 0.0
242       v_nzb_p1_for_vfc = 0.0
243    ENDIF
244
245!
246!-- Arrays to store velocity data from t-dt and the phase speeds which
247!-- are needed for radiation boundary conditions
248    IF ( outflow_l )  THEN
249       ALLOCATE( u_m_l(nzb:nzt+1,nys-1:nyn+1,1:2), &
250                 v_m_l(nzb:nzt+1,nys-1:nyn+1,0:1), &
251                 w_m_l(nzb:nzt+1,nys-1:nyn+1,0:1) )
252    ENDIF
253    IF ( outflow_r )  THEN
254       ALLOCATE( u_m_r(nzb:nzt+1,nys-1:nyn+1,nx-1:nx), &
255                 v_m_r(nzb:nzt+1,nys-1:nyn+1,nx-1:nx), &
256                 w_m_r(nzb:nzt+1,nys-1:nyn+1,nx-1:nx) )
257    ENDIF
258    IF ( outflow_l  .OR.  outflow_r )  THEN
259       ALLOCATE( c_u(nzb:nzt+1,nys-1:nyn+1), c_v(nzb:nzt+1,nys-1:nyn+1), &
260                 c_w(nzb:nzt+1,nys-1:nyn+1) )
261    ENDIF
262    IF ( outflow_s )  THEN
263       ALLOCATE( u_m_s(nzb:nzt+1,0:1,nxl-1:nxr+1), &
264                 v_m_s(nzb:nzt+1,1:2,nxl-1:nxr+1), &
265                 w_m_s(nzb:nzt+1,0:1,nxl-1:nxr+1) )
266    ENDIF
267    IF ( outflow_n )  THEN
268       ALLOCATE( u_m_n(nzb:nzt+1,ny-1:ny,nxl-1:nxr+1), &
269                 v_m_n(nzb:nzt+1,ny-1:ny,nxl-1:nxr+1), &
270                 w_m_n(nzb:nzt+1,ny-1:ny,nxl-1:nxr+1) )
271    ENDIF
272    IF ( outflow_s  .OR.  outflow_n )  THEN
273       ALLOCATE( c_u(nzb:nzt+1,nxl-1:nxr+1), c_v(nzb:nzt+1,nxl-1:nxr+1), &
274                 c_w(nzb:nzt+1,nxl-1:nxr+1) )
275    ENDIF
276
277!
278!-- Initial assignment of the pointers
279    IF ( timestep_scheme(1:5) /= 'runge' )  THEN
280
281       rif_m   => rif_1;    rif   => rif_2
282       shf_m   => shf_1;    shf   => shf_2
283       tswst_m => tswst_1;  tswst => tswst_2
284       usws_m  => usws_1;   usws  => usws_2
285       uswst_m => uswst_1;  uswst => uswst_2
286       vsws_m  => vsws_1;   vsws  => vsws_2
287       vswst_m => vswst_1;  vswst => vswst_2
288       e_m  => e_1;   e  => e_2;   e_p  => e_3;   te_m  => e_3
289       kh_m => kh_1;  kh => kh_2
290       km_m => km_1;  km => km_2
291       pt_m => pt_1;  pt => pt_2;  pt_p => pt_3;  tpt_m => pt_3
292       u_m  => u_1;   u  => u_2;   u_p  => u_3;   tu_m  => u_3
293       v_m  => v_1;   v  => v_2;   v_p  => v_3;   tv_m  => v_3
294       w_m  => w_1;   w  => w_2;   w_p  => w_3;   tw_m  => w_3
295
296       IF ( humidity  .OR.  passive_scalar )  THEN
297          qsws_m  => qsws_1;   qsws  => qsws_2
298          qswst_m => qswst_1;  qswst => qswst_2
299          q_m    => q_1;     q    => q_2;     q_p => q_3;     tq_m => q_3
300          IF ( humidity )        vpt_m  => vpt_1;   vpt  => vpt_2
301          IF ( cloud_physics )   ql   => ql_1
302          IF ( cloud_droplets )  THEN
303             ql   => ql_1
304             ql_c => ql_2
305          ENDIF
306       ENDIF
307
308    ELSE
309
310       rif   => rif_1
311       shf   => shf_1
312       tswst => tswst_1
313       usws  => usws_1
314       uswst => uswst_1
315       vsws  => vsws_1
316       vswst => vswst_1
317       e     => e_1;   e_p  => e_2;   te_m  => e_3;   e_m  => e_3
318       kh    => kh_1
319       km    => km_1
320       pt    => pt_1;  pt_p => pt_2;  tpt_m => pt_3;  pt_m => pt_3
321       u     => u_1;   u_p  => u_2;   tu_m  => u_3;   u_m  => u_3
322       v     => v_1;   v_p  => v_2;   tv_m  => v_3;   v_m  => v_3
323       w     => w_1;   w_p  => w_2;   tw_m  => w_3;   w_m  => w_3
324
325       IF ( humidity  .OR.  passive_scalar )  THEN
326          qsws   => qsws_1
327          qswst  => qswst_1
328          q      => q_1;     q_p  => q_2;     tq_m  => q_3;    q_m => q_3
329          IF ( humidity )        vpt  => vpt_1
330          IF ( cloud_physics )   ql   => ql_1
331          IF ( cloud_droplets )  THEN
332             ql   => ql_1
333             ql_c => ql_2
334          ENDIF
335       ENDIF
336
337       IF ( ocean )  THEN
338          saswsb => saswsb_1
339          saswst => saswst_1
340          sa     => sa_1;    sa_p => sa_2;    tsa_m => sa_3
341       ENDIF
342
343    ENDIF
344
345!
346!-- Initialize model variables
347    IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
348!
349!--    First model run of a possible job queue.
350!--    Initial profiles of the variables must be computes.
351       IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
352!
353!--       Use solutions of the 1D model as initial profiles,
354!--       start 1D model
355          CALL init_1d_model
356!
357!--       Transfer initial profiles to the arrays of the 3D model
358          DO  i = nxl-1, nxr+1
359             DO  j = nys-1, nyn+1
360                e(:,j,i)  = e1d
361                kh(:,j,i) = kh1d
362                km(:,j,i) = km1d
363                pt(:,j,i) = pt_init
364                u(:,j,i)  = u1d
365                v(:,j,i)  = v1d
366             ENDDO
367          ENDDO
368
369          IF ( humidity  .OR.  passive_scalar )  THEN
370             DO  i = nxl-1, nxr+1
371                DO  j = nys-1, nyn+1
372                   q(:,j,i) = q_init
373                ENDDO
374             ENDDO
375          ENDIF
376
377          IF ( .NOT. constant_diffusion )  THEN
378             DO  i = nxl-1, nxr+1
379                DO  j = nys-1, nyn+1
380                   e(:,j,i)  = e1d
381                ENDDO
382             ENDDO
383!
384!--          Store initial profiles for output purposes etc.
385             hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 )
386
387             IF ( prandtl_layer )  THEN
388                rif  = rif1d(nzb+1)
389                ts   = 0.0  ! could actually be computed more accurately in the
390                            ! 1D model. Update when opportunity arises.
391                us   = us1d
392                usws = usws1d
393                vsws = vsws1d
394             ELSE
395                ts   = 0.0  ! must be set, because used in
396                rif  = 0.0  ! flowste
397                us   = 0.0
398                usws = 0.0
399                vsws = 0.0
400             ENDIF
401
402          ELSE
403             e    = 0.0  ! must be set, because used in
404             rif  = 0.0  ! flowste
405             ts   = 0.0
406             us   = 0.0
407             usws = 0.0
408             vsws = 0.0
409          ENDIF
410          uswst = top_momentumflux_u
411          vswst = top_momentumflux_v
412
413!
414!--       In every case qs = 0.0 (see also pt)
415!--       This could actually be computed more accurately in the 1D model.
416!--       Update when opportunity arises!
417          IF ( humidity  .OR.  passive_scalar )  qs = 0.0
418
419!
420!--       inside buildings set velocities back to zero
421          IF ( topography /= 'flat' )  THEN
422             DO  i = nxl-1, nxr+1
423                DO  j = nys-1, nyn+1
424                   u(nzb:nzb_u_inner(j,i),j,i) = 0.0
425                   v(nzb:nzb_v_inner(j,i),j,i) = 0.0
426                ENDDO
427             ENDDO
428!
429!--          WARNING: The extra boundary conditions set after running the
430!--          -------  1D model impose an error on the divergence one layer
431!--                   below the topography; need to correct later
432!--          ATTENTION: Provisional correction for Piacsek & Williams
433!--          ---------  advection scheme: keep u and v zero one layer below
434!--                     the topography.
435             IF ( ibc_uv_b == 0 )  THEN
436!
437!--             Satisfying the Dirichlet condition with an extra layer below
438!--             the surface where the u and v component change their sign.
439                DO  i = nxl-1, nxr+1
440                   DO  j = nys-1, nyn+1
441                      IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = -u(1,j,i)
442                      IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = -v(1,j,i)
443                   ENDDO
444                ENDDO
445
446             ELSE
447!
448!--             Neumann condition
449                DO  i = nxl-1, nxr+1
450                   DO  j = nys-1, nyn+1
451                      IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = u(1,j,i)
452                      IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = v(1,j,i)
453                   ENDDO
454                ENDDO
455
456             ENDIF
457
458          ENDIF
459
460       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 ) &
461       THEN
462!
463!--       Use constructed initial profiles (velocity constant with height,
464!--       temperature profile with constant gradient)
465          DO  i = nxl-1, nxr+1
466             DO  j = nys-1, nyn+1
467                pt(:,j,i) = pt_init
468                u(:,j,i)  = u_init
469                v(:,j,i)  = v_init
470             ENDDO
471          ENDDO
472
473!
474!--       Set initial horizontal velocities at the lowest computational grid levels
475!--       to zero in order to avoid too small time steps caused by the diffusion
476!--       limit in the initial phase of a run (at k=1, dz/2 occurs in the
477!--       limiting formula!). The original values are stored to be later used for
478!--       volume flow control.
479          DO  i = nxl-1, nxr+1
480             DO  j = nys-1, nyn+1
481                u(nzb:nzb_u_inner(j,i)+1,j,i) = 0.0
482                v(nzb:nzb_v_inner(j,i)+1,j,i) = 0.0
483             ENDDO
484          ENDDO
485          IF ( conserve_volume_flow )  THEN
486             IF ( nxr == nx )  THEN
487                DO  j = nys, nyn
488                   k = nzb_u_inner(j,nx) + 1
489                   u_nzb_p1_for_vfc(j) = u_init(k) * dzu(k)
490                ENDDO
491             ENDIF
492             IF ( nyn == ny )  THEN
493                DO  i = nxl, nxr
494                   k = nzb_v_inner(ny,i) + 1
495                   v_nzb_p1_for_vfc(i) = v_init(k) * dzu(k)
496                ENDDO
497             ENDIF
498          ENDIF
499
500          IF ( humidity  .OR.  passive_scalar )  THEN
501             DO  i = nxl-1, nxr+1
502                DO  j = nys-1, nyn+1
503                   q(:,j,i) = q_init
504                ENDDO
505             ENDDO
506          ENDIF
507
508          IF ( ocean )  THEN
509             DO  i = nxl-1, nxr+1
510                DO  j = nys-1, nyn+1
511                   sa(:,j,i) = sa_init
512                ENDDO
513             ENDDO
514          ENDIF
515         
516          IF ( constant_diffusion )  THEN
517             km   = km_constant
518             kh   = km / prandtl_number
519             e    = 0.0
520          ELSEIF ( e_init > 0.0 )  THEN
521             DO  k = nzb+1, nzt
522                km(k,:,:) = 0.1 * l_grid(k) * SQRT( e_init )
523             ENDDO
524             km(nzb,:,:)   = km(nzb+1,:,:)
525             km(nzt+1,:,:) = km(nzt,:,:)
526             kh   = km / prandtl_number
527             e    = e_init
528          ELSE
529             IF ( .NOT. ocean )  THEN
530                kh   = 0.01   ! there must exist an initial diffusion, because
531                km   = 0.01   ! otherwise no TKE would be produced by the
532                              ! production terms, as long as not yet
533                              ! e = (u*/cm)**2 at k=nzb+1
534             ELSE
535                kh   = 0.00001
536                km   = 0.00001
537             ENDIF
538             e    = 0.0
539          ENDIF
540          rif   = 0.0
541          ts    = 0.0
542          us    = 0.0
543          usws  = 0.0
544          uswst = top_momentumflux_u
545          vsws  = 0.0
546          vswst = top_momentumflux_v
547          IF ( humidity  .OR.  passive_scalar ) qs = 0.0
548
549!
550!--       Compute initial temperature field and other constants used in case
551!--       of a sloping surface
552          IF ( sloping_surface )  CALL init_slope
553
554       ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 ) &
555       THEN
556!
557!--       Initialization will completely be done by the user
558          CALL user_init_3d_model
559
560       ENDIF
561
562!
563!--    Calculate virtual potential temperature
564       IF ( humidity ) vpt = pt * ( 1.0 + 0.61 * q )
565
566!
567!--    Store initial profiles for output purposes etc.
568       hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 )
569       hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 )
570       IF ( ibc_uv_b == 0 )  THEN
571          hom(nzb,1,5,:) = -hom(nzb+1,1,5,:)  ! due to satisfying the Dirichlet
572          hom(nzb,1,6,:) = -hom(nzb+1,1,6,:)  ! condition with an extra layer
573              ! below the surface where the u and v component change their sign
574       ENDIF
575       hom(:,1,7,:)  = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
576       hom(:,1,23,:) = SPREAD( km(:,nys,nxl), 2, statistic_regions+1 )
577       hom(:,1,24,:) = SPREAD( kh(:,nys,nxl), 2, statistic_regions+1 )
578
579       IF ( ocean )  THEN
580!
581!--       Store initial salinity profile
582          hom(:,1,26,:)  = SPREAD( sa(:,nys,nxl), 2, statistic_regions+1 )
583       ENDIF
584
585       IF ( humidity )  THEN
586!
587!--       Store initial profile of total water content, virtual potential
588!--       temperature
589          hom(:,1,26,:) = SPREAD(   q(:,nys,nxl), 2, statistic_regions+1 )
590          hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 )
591          IF ( cloud_physics  .OR.  cloud_droplets ) THEN
592!
593!--          Store initial profile of specific humidity and potential
594!--          temperature
595             hom(:,1,27,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
596             hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
597          ENDIF
598       ENDIF
599
600       IF ( passive_scalar )  THEN
601!
602!--       Store initial scalar profile
603          hom(:,1,26,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
604       ENDIF
605
606!
607!--    Initialize fluxes at bottom surface
608       IF ( use_surface_fluxes )  THEN
609
610          IF ( constant_heatflux )  THEN
611!
612!--          Heat flux is prescribed
613             IF ( random_heatflux )  THEN
614                CALL disturb_heatflux
615             ELSE
616                shf = surface_heatflux
617!
618!--             Over topography surface_heatflux is replaced by wall_heatflux(0)
619                IF ( TRIM( topography ) /= 'flat' )  THEN
620                   DO  i = nxl-1, nxr+1
621                      DO  j = nys-1, nyn+1
622                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
623                            shf(j,i) = wall_heatflux(0)
624                         ENDIF
625                      ENDDO
626                   ENDDO
627                ENDIF
628             ENDIF
629             IF ( ASSOCIATED( shf_m ) )  shf_m = shf
630          ENDIF
631
632!
633!--       Determine the near-surface water flux
634          IF ( humidity  .OR.  passive_scalar )  THEN
635             IF ( constant_waterflux )  THEN
636                qsws   = surface_waterflux
637                IF ( ASSOCIATED( qsws_m ) )  qsws_m = qsws
638             ENDIF
639          ENDIF
640
641       ENDIF
642
643!
644!--    Initialize fluxes at top surface
645!--    Currently, only the heatflux and salinity flux can be prescribed.
646!--    The latent flux is zero in this case!
647       IF ( use_top_fluxes )  THEN
648
649          IF ( constant_top_heatflux )  THEN
650!
651!--          Heat flux is prescribed
652             tswst = top_heatflux
653             IF ( ASSOCIATED( tswst_m ) )  tswst_m = tswst
654
655             IF ( humidity  .OR.  passive_scalar )  THEN
656                qswst = 0.0
657                IF ( ASSOCIATED( qswst_m ) )  qswst_m = qswst
658             ENDIF
659
660             IF ( ocean )  THEN
661                saswsb = bottom_salinityflux
662                saswst = top_salinityflux
663             ENDIF
664          ENDIF
665
666!
667!--       Initialization in case of a coupled model run
668          IF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
669             tswst = 0.0
670             IF ( ASSOCIATED( tswst_m ) )  tswst_m = tswst
671          ENDIF
672
673       ENDIF
674
675!
676!--    Initialize Prandtl layer quantities
677       IF ( prandtl_layer )  THEN
678
679          z0 = roughness_length
680
681          IF ( .NOT. constant_heatflux )  THEN 
682!
683!--          Surface temperature is prescribed. Here the heat flux cannot be
684!--          simply estimated, because therefore rif, u* and theta* would have
685!--          to be computed by iteration. This is why the heat flux is assumed
686!--          to be zero before the first time step. It approaches its correct
687!--          value in the course of the first few time steps.
688             shf   = 0.0
689             IF ( ASSOCIATED( shf_m ) )  shf_m = 0.0
690          ENDIF
691
692          IF ( humidity  .OR.  passive_scalar )  THEN
693             IF ( .NOT. constant_waterflux )  THEN
694                qsws   = 0.0
695                IF ( ASSOCIATED( qsws_m ) )   qsws_m = 0.0
696             ENDIF
697          ENDIF
698
699       ENDIF
700
701!
702!--    Calculate the initial volume flow at the right and north boundary
703       IF ( conserve_volume_flow )  THEN
704
705          volume_flow_initial_l = 0.0
706          volume_flow_area_l    = 0.0
707 
708          IF ( nxr == nx )  THEN
709             DO  j = nys, nyn
710                DO  k = nzb_2d(j,nx) + 1, nzt
711                   volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
712                                              u(k,j,nx) * dzu(k)
713                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzu(k)
714                ENDDO
715!
716!--             Correction if velocity at nzb+1 has been set zero further above
717                volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
718                                           u_nzb_p1_for_vfc(j)
719             ENDDO
720          ENDIF
721
722          IF ( nyn == ny )  THEN
723             DO  i = nxl, nxr
724                DO  k = nzb_2d(ny,i) + 1, nzt 
725                   volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
726                                              v(k,ny,i) * dzu(k)
727                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzu(k)
728                ENDDO
729!
730!--             Correction if velocity at nzb+1 has been set zero further above
731                volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
732                                           v_nzb_p1_for_vfc(i)
733             ENDDO
734          ENDIF
735
736#if defined( __parallel )
737          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
738                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
739          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
740                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
741#else
742          volume_flow_initial = volume_flow_initial_l
743          volume_flow_area    = volume_flow_area_l
744#endif 
745       ENDIF
746
747!
748!--    For the moment, perturbation pressure and vertical velocity are zero
749       p = 0.0; w = 0.0
750
751!
752!--    Initialize array sums (must be defined in first call of pres)
753       sums = 0.0
754
755!
756!--    Treating cloud physics, liquid water content and precipitation amount
757!--    are zero at beginning of the simulation
758       IF ( cloud_physics )  THEN
759          ql = 0.0
760          IF ( precipitation )  precipitation_amount = 0.0
761       ENDIF
762
763!
764!--    Initialize spectra
765       IF ( dt_dosp /= 9999999.9 )  THEN
766          spectrum_x = 0.0
767          spectrum_y = 0.0
768       ENDIF
769
770!
771!--    Impose vortex with vertical axis on the initial velocity profile
772       IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 )  THEN
773          CALL init_rankine
774       ENDIF
775
776!
777!--    Impose temperature anomaly (advection test only)
778       IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0 )  THEN
779          CALL init_pt_anomaly
780       ENDIF
781
782!
783!--    If required, change the surface temperature at the start of the 3D run
784       IF ( pt_surface_initial_change /= 0.0 )  THEN
785          pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change
786       ENDIF
787
788!
789!--    If required, change the surface humidity/scalar at the start of the 3D
790!--    run
791       IF ( ( humidity .OR. passive_scalar ) .AND. &
792            q_surface_initial_change /= 0.0 )  THEN
793          q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change
794       ENDIF
795
796!
797!--    Initialize the random number generator (from numerical recipes)
798       CALL random_function_ini
799
800!
801!--    Impose random perturbation on the horizontal velocity field and then
802!--    remove the divergences from the velocity field
803       IF ( create_disturbances )  THEN
804          CALL disturb_field( nzb_u_inner, tend, u )
805          CALL disturb_field( nzb_v_inner, tend, v )
806          n_sor = nsor_ini
807          CALL pres
808          n_sor = nsor
809       ENDIF
810
811!
812!--    Once again set the perturbation pressure explicitly to zero in order to
813!--    assure that it does not generate any divergences in the first time step.
814!--    At t=0 the velocity field is free of divergence (as constructed above).
815!--    Divergences being created during a time step are not yet known and thus
816!--    cannot be corrected during the time step yet.
817       p = 0.0
818
819!
820!--    Initialize old and new time levels.
821       IF ( timestep_scheme(1:5) /= 'runge' )  THEN
822          e_m = e; pt_m = pt; u_m = u; v_m = v; w_m = w; kh_m = kh; km_m = km
823       ELSE
824          te_m = 0.0; tpt_m = 0.0; tu_m = 0.0; tv_m = 0.0; tw_m = 0.0
825       ENDIF
826       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
827
828       IF ( humidity  .OR.  passive_scalar )  THEN
829          IF ( ASSOCIATED( q_m ) )  q_m = q
830          IF ( timestep_scheme(1:5) == 'runge' )  tq_m = 0.0
831          q_p = q
832          IF ( humidity  .AND.  ASSOCIATED( vpt_m ) )  vpt_m = vpt
833       ENDIF
834
835       IF ( ocean )  THEN
836          tsa_m = 0.0
837          sa_p  = sa
838       ENDIF
839
840!
841!--    Initialize old timelevels needed for radiation boundary conditions
842       IF ( outflow_l )  THEN
843          u_m_l(:,:,:) = u(:,:,1:2)
844          v_m_l(:,:,:) = v(:,:,0:1)
845          w_m_l(:,:,:) = w(:,:,0:1)
846       ENDIF
847       IF ( outflow_r )  THEN
848          u_m_r(:,:,:) = u(:,:,nx-1:nx)
849          v_m_r(:,:,:) = v(:,:,nx-1:nx)
850          w_m_r(:,:,:) = w(:,:,nx-1:nx)
851       ENDIF
852       IF ( outflow_s )  THEN
853          u_m_s(:,:,:) = u(:,0:1,:)
854          v_m_s(:,:,:) = v(:,1:2,:)
855          w_m_s(:,:,:) = w(:,0:1,:)
856       ENDIF
857       IF ( outflow_n )  THEN
858          u_m_n(:,:,:) = u(:,ny-1:ny,:)
859          v_m_n(:,:,:) = v(:,ny-1:ny,:)
860          w_m_n(:,:,:) = w(:,ny-1:ny,:)
861       ENDIF
862
863    ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data' ) &
864    THEN
865!
866!--    Read binary data from restart file
867       CALL read_3d_binary
868
869!
870!--    Calculate initial temperature field and other constants used in case
871!--    of a sloping surface
872       IF ( sloping_surface )  CALL init_slope
873
874!
875!--    Initialize new time levels (only done in order to set boundary values
876!--    including ghost points)
877       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
878       IF ( humidity  .OR.  passive_scalar )  q_p = q
879       IF ( ocean )  sa_p = sa
880
881    ELSE
882!
883!--    Actually this part of the programm should not be reached
884       IF ( myid == 0 )  PRINT*,'+++ init_3d_model: unknown initializing ', &
885                                                    'problem'
886       CALL local_stop
887    ENDIF
888
889!
890!-- If required, initialize dvrp-software
891    IF ( dt_dvrp /= 9999999.9 )  CALL init_dvrp
892
893    IF ( ocean )  THEN
894!
895!--    Initialize quantities needed for the ocean model
896       CALL init_ocean
897    ELSE
898!
899!--    Initialize quantities for handling cloud physics
900!--    This routine must be called before init_particles, because
901!--    otherwise, array pt_d_t, needed in data_output_dvrp (called by
902!--    init_particles) is not defined.
903       CALL init_cloud_physics
904    ENDIF
905
906!
907!-- If required, initialize particles
908    IF ( particle_advection )  CALL init_particles
909
910!
911!-- Initialize quantities for special advections schemes
912    CALL init_advec
913
914!
915!-- Initialize Rayleigh damping factors
916    rdf = 0.0
917    IF ( rayleigh_damping_factor /= 0.0 )  THEN
918       IF ( .NOT. ocean )  THEN
919          DO  k = nzb+1, nzt
920             IF ( zu(k) >= rayleigh_damping_height )  THEN
921                rdf(k) = rayleigh_damping_factor * &
922                      ( SIN( pi * 0.5 * ( zu(k) - rayleigh_damping_height )    &
923                                      / ( zu(nzt) - rayleigh_damping_height ) )&
924                      )**2
925             ENDIF
926          ENDDO
927       ELSE
928          DO  k = nzt, nzb+1, -1
929             IF ( zu(k) <= rayleigh_damping_height )  THEN
930                rdf(k) = rayleigh_damping_factor * &
931                      ( SIN( pi * 0.5 * ( rayleigh_damping_height - zu(k) )    &
932                                      / ( rayleigh_damping_height - zu(nzb+1)))&
933                      )**2
934             ENDIF
935          ENDDO
936       ENDIF
937    ENDIF
938
939!
940!-- Initialize diffusivities used within the outflow damping layer in case of
941!-- non-cyclic lateral boundaries. A linear increase is assumed over the first
942!-- half of the width of the damping layer
943    IF ( bc_lr == 'dirichlet/radiation' )  THEN
944
945       DO  i = nxl-1, nxr+1
946          IF ( i >= nx - outflow_damping_width )  THEN
947             km_damp_x(i) = km_damp_max * MIN( 1.0,                    &
948                            ( i - ( nx - outflow_damping_width ) ) /   &
949                            REAL( outflow_damping_width/2 )            &
950                                             )
951          ELSE
952             km_damp_x(i) = 0.0
953          ENDIF
954       ENDDO
955
956    ELSEIF ( bc_lr == 'radiation/dirichlet' )  THEN
957
958       DO  i = nxl-1, nxr+1
959          IF ( i <= outflow_damping_width )  THEN
960             km_damp_x(i) = km_damp_max * MIN( 1.0,                    &
961                            ( outflow_damping_width - i ) /            &
962                            REAL( outflow_damping_width/2 )            &
963                                             )
964          ELSE
965             km_damp_x(i) = 0.0
966          ENDIF
967       ENDDO
968
969    ENDIF
970
971    IF ( bc_ns == 'dirichlet/radiation' )  THEN
972
973       DO  j = nys-1, nyn+1
974          IF ( j >= ny - outflow_damping_width )  THEN
975             km_damp_y(j) = km_damp_max * MIN( 1.0,                    &
976                            ( j - ( ny - outflow_damping_width ) ) /   &
977                            REAL( outflow_damping_width/2 )            &
978                                             )
979          ELSE
980             km_damp_y(j) = 0.0
981          ENDIF
982       ENDDO
983
984    ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
985
986       DO  j = nys-1, nyn+1
987          IF ( j <= outflow_damping_width )  THEN
988             km_damp_y(j) = km_damp_max * MIN( 1.0,                    &
989                            ( outflow_damping_width - j ) /            &
990                            REAL( outflow_damping_width/2 )            &
991                                             )
992          ELSE
993             km_damp_y(j) = 0.0
994          ENDIF
995       ENDDO
996
997    ENDIF
998
999!
1000!-- Initialize local summation arrays for UP flow_statistics. This is necessary
1001!-- because they may not yet have been initialized when they are called from
1002!-- flow_statistics (or - depending on the chosen model run - are never
1003!-- initialized)
1004    sums_divnew_l      = 0.0
1005    sums_divold_l      = 0.0
1006    sums_l_l           = 0.0
1007    sums_up_fraction_l = 0.0
1008    sums_wsts_bc_l     = 0.0
1009
1010!
1011!-- Pre-set masks for regional statistics. Default is the total model domain.
1012    rmask = 1.0
1013
1014!
1015!-- User-defined initializing actions. Check afterwards, if maximum number
1016!-- of allowed timeseries is not exceeded
1017    CALL user_init
1018
1019    IF ( dots_num > dots_max )  THEN
1020       IF ( myid == 0 )  THEN
1021          PRINT*, '+++ user_init: number of time series quantities exceeds', &
1022                  ' its maximum of dots_max = ', dots_max
1023          PRINT*, '    Please increase dots_max in modules.f90.'
1024       ENDIF
1025       CALL local_stop
1026    ENDIF
1027
1028!
1029!-- Input binary data file is not needed anymore. This line must be placed
1030!-- after call of user_init!
1031    CALL close_file( 13 )
1032
1033!
1034!-- Compute total sum of active mask grid points
1035!-- ngp_2dh: number of grid points of a horizontal cross section through the
1036!--          total domain
1037!-- ngp_3d:  number of grid points of the total domain
1038!-- Note: The lower vertical index nzb_s_outer imposes a small error on the 2D
1039!-- ----  averages of staggered variables such as u and v due to the topography
1040!--       arrangement on the staggered grid. Maybe revise later.
1041    ngp_2dh_outer_l = 0
1042    ngp_2dh_outer   = 0
1043    ngp_2dh_l       = 0
1044    ngp_2dh         = 0
1045    ngp_3d_inner_l  = 0
1046    ngp_3d_inner    = 0
1047    ngp_3d          = 0
1048    ngp_sums        = ( nz + 2 ) * ( pr_palm + max_pr_user )
1049
1050    DO  sr = 0, statistic_regions
1051       DO  i = nxl, nxr
1052          DO  j = nys, nyn
1053             IF ( rmask(j,i,sr) == 1.0 )  THEN
1054!
1055!--             All xy-grid points
1056                ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1
1057!
1058!--             xy-grid points above topography
1059                DO  k = nzb_s_outer(j,i), nz + 1
1060                   ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr) + 1
1061                ENDDO
1062!
1063!--             All grid points of the total domain above topography
1064                ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr) + &
1065                                     ( nz - nzb_s_inner(j,i) + 2 )
1066             ENDIF
1067          ENDDO
1068       ENDDO
1069    ENDDO
1070
1071    sr = statistic_regions + 1
1072#if defined( __parallel )
1073    CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM,  &
1074                        comm2d, ierr )
1075    CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr, &
1076                        MPI_INTEGER, MPI_SUM, comm2d, ierr )
1077    CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner(0), sr, MPI_INTEGER, &
1078                        MPI_SUM, comm2d, ierr )
1079#else
1080    ngp_2dh       = ngp_2dh_l
1081    ngp_2dh_outer = ngp_2dh_outer_l
1082    ngp_3d_inner  = ngp_3d_inner_l
1083#endif
1084
1085    ngp_3d = ngp_2dh * ( nz + 2 )
1086
1087!
1088!-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics,
1089!-- buoyancy, etc. A zero value will occur for cases where all grid points of
1090!-- the respective subdomain lie below the surface topography
1091    ngp_2dh_outer = MAX( 1, ngp_2dh_outer(:,:) ) 
1092    ngp_3d_inner  = MAX( 1, ngp_3d_inner(:)    )
1093
1094    DEALLOCATE( ngp_2dh_l, ngp_2dh_outer_l, ngp_3d_inner_l )
1095
1096
1097 END SUBROUTINE init_3d_model
Note: See TracBrowser for help on using the repository browser.