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

Last change on this file since 48 was 48, checked in by raasch, 17 years ago

preliminary version, several changes to be explained later

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