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

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

preliminary version for coupled runs

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