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

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

preliminary version, several changes to be explained later

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