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

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

further preliminary uncomplete changes for ocean version

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