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

Last change on this file since 46 was 46, checked in by raasch, 18 years ago

updating parts of Marcus changes

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