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

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

preliminary uncomplete changes for ocean version

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