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

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

Initial repository layout and content

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