source: palm/trunk/SOURCE/init_1d_model.f90 @ 82

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

vorlaeufige Standalone-Version fuer Linux-Cluster

  • Property svn:keywords set to Id
File size: 33.6 KB
Line 
1 SUBROUTINE init_1d_model
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
6! Preprocessor strings for different linux clusters changed to "lc",
7! routine local_flush is used for buffer flushing
8!
9! Former revisions:
10! -----------------
11! $Id: init_1d_model.f90 82 2007-04-16 15:40:52Z raasch $
12!
13! 75 2007-03-22 09:54:05Z raasch
14! Bugfix: preset of tendencies te_em, te_um, te_vm,
15! moisture renamed humidity
16!
17! RCS Log replace by Id keyword, revision history cleaned up
18!
19! Revision 1.21  2006/06/02 15:19:57  raasch
20! cpp-directives extended for lctit
21!
22! Revision 1.1  1998/03/09 16:22:10  raasch
23! Initial revision
24!
25!
26! Description:
27! ------------
28! 1D-model to initialize the 3D-arrays.
29! The temperature profile is set as steady and a corresponding steady solution
30! of the wind profile is being computed.
31! All subroutines required can be found within this file.
32!------------------------------------------------------------------------------!
33
34    USE arrays_3d
35    USE indices
36    USE model_1d
37    USE control_parameters
38
39    IMPLICIT NONE
40
41    CHARACTER (LEN=9) ::  time_to_string
42    INTEGER ::  k
43    REAL    ::  lambda
44
45!
46!-- Allocate required 1D-arrays
47    ALLOCATE( e1d(nzb:nzt+1),    e1d_m(nzb:nzt+1),   e1d_p(nzb:nzt+1), &
48              kh1d(nzb:nzt+1),   kh1d_m(nzb:nzt+1),  km1d(nzb:nzt+1),  &
49              km1d_m(nzb:nzt+1), l_black(nzb:nzt+1), l1d(nzb:nzt+1),   &
50              l1d_m(nzb:nzt+1),  rif1d(nzb:nzt+1),   te_e(nzb:nzt+1),  &
51              te_em(nzb:nzt+1),  te_u(nzb:nzt+1),    te_um(nzb:nzt+1), &
52              te_v(nzb:nzt+1),   te_vm(nzb:nzt+1),   u1d(nzb:nzt+1),   &
53              u1d_m(nzb:nzt+1),  u1d_p(nzb:nzt+1),   v1d(nzb:nzt+1),   &
54              v1d_m(nzb:nzt+1),  v1d_p(nzb:nzt+1) )
55
56!
57!-- Initialize arrays
58    IF ( constant_diffusion )  THEN
59       km1d   = km_constant
60       km1d_m = km_constant
61       kh1d   = km_constant / prandtl_number
62       kh1d_m = km_constant / prandtl_number
63    ELSE
64       e1d = 0.0; e1d_m = 0.0; e1d_p = 0.0
65       kh1d = 0.0; kh1d_m = 0.0; km1d = 0.0; km1d_m = 0.0
66       rif1d = 0.0
67!
68!--    Compute the mixing length
69       l_black(nzb) = 0.0
70
71       IF ( TRIM( mixing_length_1d ) == 'blackadar' )  THEN
72!
73!--       Blackadar mixing length
74          IF ( f /= 0.0 )  THEN
75             lambda = 2.7E-4 * SQRT( ug(nzt+1)**2 + vg(nzt+1)**2 ) / f + 1E-10
76          ELSE
77             lambda = 30.0
78          ENDIF
79
80          DO  k = nzb+1, nzt+1
81             l_black(k) = kappa * zu(k) / ( 1.0 + kappa * zu(k) / lambda )
82          ENDDO
83
84       ELSEIF ( TRIM( mixing_length_1d ) == 'as_in_3d_model' )  THEN
85!
86!--       Use the same mixing length as in 3D model
87          l_black(1:nzt) = l_grid
88          l_black(nzt+1) = l_black(nzt)
89
90       ENDIF
91
92!
93!--    Adjust mixing length to the prandtl mixing length (within the prandtl
94!--    layer)
95       IF ( adjust_mixing_length  .AND.  prandtl_layer )  THEN
96          k = nzb+1
97          l_black(k) = MIN( l_black(k), kappa * zu(k) )
98       ENDIF
99    ENDIF
100    l1d   = l_black
101    l1d_m = l_black
102    u1d   = u_init
103    u1d_m = u_init
104    u1d_p = u_init
105    v1d   = v_init
106    v1d_m = v_init
107    v1d_p = v_init
108
109!
110!-- Set initial horizontal velocities at the lowest grid levels to a very small
111!-- value in order to avoid too small time steps caused by the diffusion limit
112!-- in the initial phase of a run (at k=1, dz/2 occurs in the limiting formula!)
113    u1d(0:1)   = 0.1
114    u1d_m(0:1) = 0.1
115    u1d_p(0:1) = 0.1
116    v1d(0:1)   = 0.1
117    v1d_m(0:1) = 0.1
118    v1d_p(0:1) = 0.1
119
120!
121!-- For u*, theta* and the momentum fluxes plausible values are set
122    IF ( prandtl_layer )  THEN
123       us1d = 0.1   ! without initial friction the flow would not change
124    ELSE
125       e1d(nzb+1)  = 1.0
126       km1d(nzb+1) = 1.0
127       us1d = 0.0
128    ENDIF
129    ts1d = 0.0
130    usws1d = 0.0; usws1d_m = 0.0
131    vsws1d = 0.0; vsws1d_m = 0.0
132    z01d = roughness_length
133    IF ( humidity .OR. passive_scalar )  qs1d = 0.0
134
135!
136!-- Tendencies must be preset in order to avoid runtime errors within the
137!-- first Runge-Kutta step
138    te_em = 0.0
139    te_um = 0.0
140    te_vm = 0.0
141
142!
143!-- Set start time in hh:mm:ss - format
144    simulated_time_chr = time_to_string( simulated_time_1d )
145
146!
147!-- Integrate the 1D-model equations using the leap-frog scheme
148    CALL time_integration_1d
149
150
151 END SUBROUTINE init_1d_model
152
153
154
155 SUBROUTINE time_integration_1d
156
157!------------------------------------------------------------------------------!
158! Description:
159! ------------
160! Leap-frog time differencing scheme for the 1D-model.
161!------------------------------------------------------------------------------!
162
163    USE arrays_3d
164    USE control_parameters
165    USE indices
166    USE model_1d
167    USE pegrid
168
169    IMPLICIT NONE
170
171    CHARACTER (LEN=9) ::  time_to_string
172    INTEGER ::  k
173    REAL    ::  a, b, dissipation, dpt_dz, flux, kmzm, kmzp, l_stable, pt_0, &
174                uv_total
175
176!
177!-- Determine the time step at the start of a 1D-simulation and
178!-- determine and printout quantities used for run control
179    CALL timestep_1d
180    CALL run_control_1d
181
182!
183!-- Start of time loop
184    DO  WHILE ( simulated_time_1d < end_time_1d  .AND.  .NOT. stop_dt_1d )
185
186!
187!--    Depending on the timestep scheme, carry out one or more intermediate
188!--    timesteps
189
190       intermediate_timestep_count = 0
191       DO  WHILE ( intermediate_timestep_count < &
192                   intermediate_timestep_count_max )
193
194          intermediate_timestep_count = intermediate_timestep_count + 1
195
196          CALL timestep_scheme_steering
197
198!
199!--       Compute all tendency terms. If a Prandtl-layer is simulated, k starts
200!--       at nzb+2.
201          DO  k = nzb_diff, nzt
202
203             kmzm = 0.5 * ( km1d_m(k-1) + km1d_m(k) )
204             kmzp = 0.5 * ( km1d_m(k) + km1d_m(k+1) )
205!
206!--          u-component
207             te_u(k) =  f * ( v1d(k) - vg(k) ) + ( &
208                              kmzp * ( u1d_m(k+1) - u1d_m(k) ) * ddzu(k+1) &
209                            - kmzm * ( u1d_m(k) - u1d_m(k-1) ) * ddzu(k)   &
210                                              ) * ddzw(k)
211!
212!--          v-component
213             te_v(k) = -f * ( u1d(k) - ug(k) ) + ( &
214                              kmzp * ( v1d_m(k+1) - v1d_m(k) ) * ddzu(k+1) &
215                            - kmzm * ( v1d_m(k) - v1d_m(k-1) ) * ddzu(k)   &
216                                              ) * ddzw(k)
217          ENDDO
218          IF ( .NOT. constant_diffusion )  THEN
219             DO  k = nzb_diff, nzt
220!
221!--             TKE
222                kmzm = 0.5 * ( km1d_m(k-1) + km1d_m(k) )
223                kmzp = 0.5 * ( km1d_m(k) + km1d_m(k+1) )
224                IF ( .NOT. humidity )  THEN
225                   pt_0 = pt_init(k)
226                   flux =  ( pt_init(k+1)-pt_init(k-1) ) * dd2zu(k)
227                ELSE
228                   pt_0 = pt_init(k) * ( 1.0 + 0.61 * q_init(k) )
229                   flux = ( ( pt_init(k+1) - pt_init(k-1) ) +                 &
230                            0.61 * pt_init(k) * ( q_init(k+1) - q_init(k-1) ) &
231                          ) * dd2zu(k)
232                ENDIF
233
234                IF ( dissipation_1d == 'detering' )  THEN
235!
236!--                According to Detering, c_e=0.064
237                   dissipation = 0.064 * e1d_m(k) * SQRT( e1d_m(k) ) / l1d_m(k)
238                ELSEIF ( dissipation_1d == 'as_in_3d_model' )  THEN
239                   dissipation = ( 0.19 + 0.74 * l1d_m(k) / l_grid(k) ) &
240                                 * e1d_m(k) * SQRT( e1d_m(k) ) / l1d_m(k)
241                ENDIF
242
243                te_e(k) = km1d(k) * ( ( ( u1d(k+1) - u1d(k-1) ) * dd2zu(k) )**2&
244                                    + ( ( v1d(k+1) - v1d(k-1) ) * dd2zu(k) )**2&
245                                    )                                          &
246                                    - g / pt_0 * kh1d(k) * flux                &
247                                    +            (                             &
248                                 kmzp * ( e1d_m(k+1) - e1d_m(k) ) * ddzu(k+1)  &
249                               - kmzm * ( e1d_m(k) - e1d_m(k-1) ) * ddzu(k)    &
250                                                 ) * ddzw(k)                   &
251                               - dissipation
252             ENDDO
253          ENDIF
254
255!
256!--       Tendency terms at the top of the Prandtl-layer.
257!--       Finite differences of the momentum fluxes are computed using half the
258!--       normal grid length (2.0*ddzw(k)) for the sake of enhanced accuracy
259          IF ( prandtl_layer )  THEN
260
261             k = nzb+1
262             kmzm = 0.5 * ( km1d_m(k-1) + km1d_m(k) )
263             kmzp = 0.5 * ( km1d_m(k) + km1d_m(k+1) )
264             IF ( .NOT. humidity )  THEN
265                pt_0 = pt_init(k)
266                flux =  ( pt_init(k+1)-pt_init(k-1) ) * dd2zu(k)
267             ELSE
268                pt_0 = pt_init(k) * ( 1.0 + 0.61 * q_init(k) )
269                flux = ( ( pt_init(k+1) - pt_init(k-1) ) +                 &
270                         0.61 * pt_init(k) * ( q_init(k+1) - q_init(k-1) ) &
271                       ) * dd2zu(k)
272             ENDIF
273
274             IF ( dissipation_1d == 'detering' )  THEN
275!
276!--             According to Detering, c_e=0.064
277                dissipation = 0.064 * e1d_m(k) * SQRT( e1d_m(k) ) / l1d_m(k)
278             ELSEIF ( dissipation_1d == 'as_in_3d_model' )  THEN
279                dissipation = ( 0.19 + 0.74 * l1d_m(k) / l_grid(k) ) &
280                              * e1d_m(k) * SQRT( e1d_m(k) ) / l1d_m(k)
281             ENDIF
282
283!
284!--          u-component
285             te_u(k) = f * ( v1d(k) - vg(k) ) + ( &
286                       kmzp * ( u1d_m(k+1) - u1d_m(k) ) * ddzu(k+1) + usws1d_m &
287                                               ) * 2.0 * ddzw(k)
288!
289!--          v-component
290             te_v(k) = -f * ( u1d(k) - ug(k) ) + ( &
291                       kmzp * ( v1d_m(k+1) - v1d_m(k) ) * ddzu(k+1) + vsws1d_m &
292                                           ) * 2.0 * ddzw(k)
293!
294!--          TKE
295             te_e(k) = km1d(k) * ( ( ( u1d(k+1) - u1d(k-1) ) * dd2zu(k) )**2   &
296                                 + ( ( v1d(k+1) - v1d(k-1) ) * dd2zu(k) )**2   &
297                                 )                                             &
298                                 - g / pt_0 * kh1d(k) * flux                   &
299                                 +           (                                 &
300                              kmzp * ( e1d_m(k+1) - e1d_m(k) ) * ddzu(k+1)     &
301                            - kmzm * ( e1d_m(k) - e1d_m(k-1) ) * ddzu(k)       &
302                                              ) * ddzw(k)                      &
303                            - dissipation
304          ENDIF
305
306!
307!--       Prognostic equations for all 1D variables
308          DO  k = nzb+1, nzt
309
310             u1d_p(k) = ( 1. - tsc(1) ) * u1d_m(k) + &
311                        tsc(1) * u1d(k) + dt_1d * ( tsc(2) * te_u(k) + &
312                                                    tsc(3) * te_um(k) ) 
313             v1d_p(k) = ( 1. - tsc(1) ) * v1d_m(k) + &
314                        tsc(1) * v1d(k) + dt_1d * ( tsc(2) * te_v(k) + &
315                                                    tsc(3) * te_vm(k) ) 
316
317          ENDDO
318          IF ( .NOT. constant_diffusion )  THEN
319             DO  k = nzb+1, nzt
320
321                e1d_p(k) = ( 1. - tsc(1) ) * e1d_m(k) + &
322                           tsc(1) * e1d(k) + dt_1d * ( tsc(2) * te_e(k) + &
323                                                       tsc(3) * te_em(k) ) 
324
325             ENDDO
326!
327!--          Eliminate negative TKE values, which can result from the
328!--          integration due to numerical inaccuracies. In such cases the TKE
329!--          value is reduced to 10 percent of its old value.
330             WHERE ( e1d_p < 0.0 )  e1d_p = 0.1 * e1d
331          ENDIF
332
333!
334!--       Calculate tendencies for the next Runge-Kutta step
335          IF ( timestep_scheme(1:5) == 'runge' ) THEN
336             IF ( intermediate_timestep_count == 1 )  THEN
337
338                DO  k = nzb+1, nzt
339                   te_um(k) = te_u(k)
340                   te_vm(k) = te_v(k)
341                ENDDO
342
343                IF ( .NOT. constant_diffusion )  THEN
344                   DO k = nzb+1, nzt
345                      te_em(k) = te_e(k)
346                   ENDDO
347                ENDIF
348
349             ELSEIF ( intermediate_timestep_count < &
350                         intermediate_timestep_count_max )  THEN
351
352                DO  k = nzb+1, nzt
353                   te_um(k) = -9.5625 * te_u(k) + 5.3125 * te_um(k)
354                   te_vm(k) = -9.5625 * te_v(k) + 5.3125 * te_vm(k)
355                ENDDO
356
357                IF ( .NOT. constant_diffusion )  THEN
358                   DO k = nzb+1, nzt
359                      te_em(k) = -9.5625 * te_e(k) + 5.3125 * te_em(k)
360                   ENDDO
361                ENDIF
362
363             ENDIF
364          ENDIF
365
366
367!
368!--       Boundary conditions for the prognostic variables.
369!--       At the top boundary (nzt+1) u,v and e keep their initial values
370!--       (ug(nzt+1), vg(nzt+1), 0), at the bottom boundary the mirror
371!--       boundary condition applies to u and v.
372!--       The boundary condition for e is set further below ( (u*/cm)**2 ).
373          u1d_p(nzb) = -u1d_p(nzb+1)
374          v1d_p(nzb) = -v1d_p(nzb+1)
375
376!
377!--       If necessary, apply the time filter
378          IF ( asselin_filter_factor /= 0.0  .AND. &
379               timestep_scheme(1:5) /= 'runge' )  THEN
380
381             u1d = u1d + asselin_filter_factor * ( u1d_p - 2.0 * u1d + u1d_m )
382             v1d = v1d + asselin_filter_factor * ( v1d_p - 2.0 * v1d + v1d_m )
383
384             IF ( .NOT. constant_diffusion )  THEN
385                e1d = e1d + asselin_filter_factor * &
386                            ( e1d_p - 2.0 * e1d + e1d_m )
387             ENDIF
388
389          ENDIF
390
391!
392!--       Swap the time levels in preparation for the next time step.
393          IF ( timestep_scheme(1:4) == 'leap' )  THEN
394             u1d_m  = u1d
395             v1d_m  = v1d
396             IF ( .NOT. constant_diffusion )  THEN
397                e1d_m  = e1d
398                kh1d_m = kh1d  ! The old diffusion quantities are required for
399                km1d_m = km1d  ! explicit diffusion in the leap-frog scheme.
400                l1d_m  = l1d
401                IF ( prandtl_layer )  THEN
402                   usws1d_m = usws1d
403                   vsws1d_m = vsws1d
404                ENDIF
405             ENDIF
406          ENDIF
407          u1d  = u1d_p
408          v1d  = v1d_p
409          IF ( .NOT. constant_diffusion )  THEN
410             e1d  = e1d_p
411          ENDIF
412
413!
414!--       Compute diffusion quantities
415          IF ( .NOT. constant_diffusion )  THEN
416
417!
418!--          First compute the vertical fluxes in the Prandtl-layer
419             IF ( prandtl_layer )  THEN
420!
421!--             Compute theta* using Rif numbers of the previous time step
422                IF ( rif1d(1) >= 0.0 )  THEN
423!
424!--                Stable stratification
425                   ts1d = kappa * ( pt_init(nzb+1) - pt_init(nzb) ) /      &
426                          ( LOG( zu(nzb+1) / z01d ) + 5.0 * rif1d(nzb+1) * &
427                                          ( zu(nzb+1) - z01d ) / zu(nzb+1) &
428                          )
429                ELSE
430!
431!--                Unstable stratification
432                   a = SQRT( 1.0 - 16.0 * rif1d(nzb+1) )
433                   b = SQRT( 1.0 - 16.0 * rif1d(nzb+1) / zu(nzb+1) * z01d )
434!
435!--                In the borderline case the formula for stable stratification
436!--                must be applied, because otherwise a zero division would
437!--                occur in the argument of the logarithm.
438                   IF ( a == 0.0  .OR.  b == 0.0 )  THEN
439                      ts1d = kappa * ( pt_init(nzb+1) - pt_init(nzb) ) /      &
440                             ( LOG( zu(nzb+1) / z01d ) + 5.0 * rif1d(nzb+1) * &
441                                             ( zu(nzb+1) - z01d ) / zu(nzb+1) &
442                             )
443                   ELSE
444                      ts1d = kappa * ( pt_init(nzb+1) - pt_init(nzb) ) / &
445                             LOG( (a-1.0) / (a+1.0) * (b+1.0) / (b-1.0) )
446                   ENDIF
447                ENDIF
448
449             ENDIF    ! prandtl_layer
450
451!
452!--          Compute the Richardson-flux numbers,
453!--          first at the top of the Prandtl-layer using u* of the previous
454!--          time step (+1E-30, if u* = 0), then in the remaining area. There
455!--          the rif-numbers of the previous time step are used.
456
457             IF ( prandtl_layer )  THEN
458                IF ( .NOT. humidity )  THEN
459                   pt_0 = pt_init(nzb+1)
460                   flux = ts1d
461                ELSE
462                   pt_0 = pt_init(nzb+1) * ( 1.0 + 0.61 * q_init(nzb+1) )
463                   flux = ts1d + 0.61 * pt_init(k) * qs1d
464                ENDIF
465                rif1d(nzb+1) = zu(nzb+1) * kappa * g * flux / &
466                               ( pt_0 * ( us1d**2 + 1E-30 ) )
467             ENDIF
468
469             DO  k = nzb_diff, nzt
470                IF ( .NOT. humidity )  THEN
471                   pt_0 = pt_init(k)
472                   flux = ( pt_init(k+1) - pt_init(k-1) ) * dd2zu(k)
473                ELSE
474                   pt_0 = pt_init(k) * ( 1.0 + 0.61 * q_init(k) )
475                   flux = ( ( pt_init(k+1) - pt_init(k-1) )                    &
476                            + 0.61 * pt_init(k) * ( q_init(k+1) - q_init(k-1) )&
477                          ) * dd2zu(k)
478                ENDIF
479                IF ( rif1d(k) >= 0.0 )  THEN
480                   rif1d(k) = g / pt_0 * flux /                              &
481                              (  ( ( u1d(k+1) - u1d(k-1) ) * dd2zu(k) )**2   &
482                               + ( ( v1d(k+1) - v1d(k-1) ) * dd2zu(k) )**2   &
483                               + 1E-30                                       &
484                              )
485                ELSE
486                   rif1d(k) = g / pt_0 * flux /                              &
487                              (  ( ( u1d(k+1) - u1d(k-1) ) * dd2zu(k) )**2   &
488                               + ( ( v1d(k+1) - v1d(k-1) ) * dd2zu(k) )**2   &
489                               + 1E-30                                       &
490                              ) * ( 1.0 - 16.0 * rif1d(k) )**0.25
491                ENDIF
492             ENDDO
493!
494!--          Richardson-numbers must remain restricted to a realistic value
495!--          range. It is exceeded excessively for very small velocities
496!--          (u,v --> 0).
497             WHERE ( rif1d < rif_min )  rif1d = rif_min
498             WHERE ( rif1d > rif_max )  rif1d = rif_max
499
500!
501!--          Compute u* from the absolute velocity value
502             IF ( prandtl_layer )  THEN
503                uv_total = SQRT( u1d(nzb+1)**2 + v1d(nzb+1)**2 )
504
505                IF ( rif1d(nzb+1) >= 0.0 )  THEN
506!
507!--                Stable stratification
508                   us1d = kappa * uv_total / (                                 &
509                             LOG( zu(nzb+1) / z01d ) + 5.0 * rif1d(nzb+1) *    &
510                                              ( zu(nzb+1) - z01d ) / zu(nzb+1) &
511                                             )
512                ELSE
513!
514!--                Unstable stratification
515                   a = 1.0 / SQRT( SQRT( 1.0 - 16.0 * rif1d(nzb+1) ) )
516                   b = 1.0 / SQRT( SQRT( 1.0 - 16.0 * rif1d(nzb+1) / zu(nzb+1) &
517                                                    * z01d ) )
518!
519!--                In the borderline case the formula for stable stratification
520!--                must be applied, because otherwise a zero division would
521!--                occur in the argument of the logarithm.
522                   IF ( a == 1.0  .OR.  b == 1.0 )  THEN
523                      us1d = kappa * uv_total / (                            &
524                             LOG( zu(nzb+1) / z01d ) +                       &
525                             5.0 * rif1d(nzb+1) * ( zu(nzb+1) - z01d ) /     &
526                                                  zu(nzb+1) )
527                   ELSE
528                      us1d = kappa * uv_total / (                              &
529                                 LOG( (1.0+b) / (1.0-b) * (1.0-a) / (1.0+a) ) +&
530                                 2.0 * ( ATAN( b ) - ATAN( a ) )               &
531                                                )
532                   ENDIF
533                ENDIF
534
535!
536!--             Compute the momentum fluxes for the diffusion terms
537                usws1d  = - u1d(nzb+1) / uv_total * us1d**2
538                vsws1d  = - v1d(nzb+1) / uv_total * us1d**2
539
540!
541!--             Boundary condition for the turbulent kinetic energy at the top
542!--             of the Prandtl-layer. c_m = 0.4 according to Detering.
543!--             Additional Neumann condition de/dz = 0 at nzb is set to ensure
544!--             compatibility with the 3D model.
545                IF ( ibc_e_b == 2 )  THEN
546                   e1d(nzb+1) = ( us1d / 0.1 )**2
547!                  e1d(nzb+1) = ( us1d / 0.4 )**2  !not used so far, see also
548                                                   !prandtl_fluxes
549                ENDIF
550                e1d(nzb) = e1d(nzb+1)
551
552                IF ( humidity .OR. passive_scalar ) THEN
553!
554!--                Compute q*
555                   IF ( rif1d(1) >= 0.0 )  THEN
556!
557!--                Stable stratification
558                   qs1d = kappa * ( q_init(nzb+1) - q_init(nzb) ) /        &
559                          ( LOG( zu(nzb+1) / z01d ) + 5.0 * rif1d(nzb+1) * &
560                                          ( zu(nzb+1) - z01d ) / zu(nzb+1) &
561                          )
562                ELSE
563!
564!--                Unstable stratification
565                   a = SQRT( 1.0 - 16.0 * rif1d(nzb+1) )
566                   b = SQRT( 1.0 - 16.0 * rif1d(nzb+1) / zu(nzb+1) * z01d )
567!
568!--                In the borderline case the formula for stable stratification
569!--                must be applied, because otherwise a zero division would
570!--                occur in the argument of the logarithm.
571                   IF ( a == 1.0  .OR.  b == 1.0 )  THEN
572                      qs1d = kappa * ( q_init(nzb+1) - q_init(nzb) ) /        &
573                             ( LOG( zu(nzb+1) / z01d ) + 5.0 * rif1d(nzb+1) * &
574                                             ( zu(nzb+1) - z01d ) / zu(nzb+1) &
575                             )
576                   ELSE
577                      qs1d = kappa * ( q_init(nzb+1) - q_init(nzb) ) / &
578                             LOG( (a-1.0) / (a+1.0) * (b+1.0) / (b-1.0) )
579                   ENDIF
580                ENDIF               
581                ELSE
582                   qs1d = 0.0
583                ENDIF             
584
585             ENDIF   !  prandtl_layer
586
587!
588!--          Compute the diabatic mixing length
589             IF ( mixing_length_1d == 'blackadar' )  THEN
590                DO  k = nzb+1, nzt
591                   IF ( rif1d(k) >= 0.0 )  THEN
592                      l1d(k) = l_black(k) / ( 1.0 + 5.0 * rif1d(k) )
593                   ELSE
594                      l1d(k) = l_black(k) * ( 1.0 - 16.0 * rif1d(k) )**0.25
595                   ENDIF
596                   l1d(k) = l_black(k)
597                ENDDO
598
599             ELSEIF ( mixing_length_1d == 'as_in_3d_model' )  THEN
600                DO  k = nzb+1, nzt
601                   dpt_dz = ( pt_init(k+1) - pt_init(k-1) ) * dd2zu(k)
602                   IF ( dpt_dz > 0.0 )  THEN
603                      l_stable = 0.76 * SQRT( e1d(k) ) / &
604                                     SQRT( g / pt_init(k) * dpt_dz ) + 1E-5
605                   ELSE
606                      l_stable = l_grid(k)
607                   ENDIF
608                   l1d(k) = MIN( l_grid(k), l_stable )
609                ENDDO
610             ENDIF
611
612!
613!--          Adjust mixing length to the prandtl mixing length
614             IF ( adjust_mixing_length  .AND.  prandtl_layer )  THEN
615                k = nzb+1
616                IF ( rif1d(k) >= 0.0 )  THEN
617                   l1d(k) = MIN( l1d(k), kappa * zu(k) / ( 1.0 + 5.0 * &
618                                                           rif1d(k) ) )
619                ELSE
620                   l1d(k) = MIN( l1d(k), kappa * zu(k) *          &
621                                  SQRT( SQRT( 1.0 - 16.0 * rif1d(k) ) ) )
622                ENDIF
623             ENDIF
624
625!
626!--          Compute the diffusion coefficients for momentum via the
627!--          corresponding Prandtl-layer relationship and according to
628!--          Prandtl-Kolmogorov, respectively. The unstable stratification is
629!--          computed via the adiabatic mixing length, for the unstability has
630!--          already been taken account of via the TKE (cf. also Diss.).
631             IF ( prandtl_layer )  THEN
632                IF ( rif1d(nzb+1) >= 0.0 )  THEN
633                   km1d(nzb+1) = us1d * kappa * zu(nzb+1) / &
634                                 ( 1.0 + 5.0 * rif1d(nzb+1) )
635                ELSE
636                   km1d(nzb+1) = us1d * kappa * zu(nzb+1) * &
637                                 ( 1.0 - 16.0 * rif1d(nzb+1) )**0.25
638                ENDIF
639             ENDIF
640             DO  k = nzb_diff, nzt
641!                km1d(k) = 0.4 * SQRT( e1d(k) ) !changed: adjustment to 3D-model
642                km1d(k) = 0.1 * SQRT( e1d(k) )
643                IF ( rif1d(k) >= 0.0 )  THEN
644                   km1d(k) = km1d(k) * l1d(k)
645                ELSE
646                   km1d(k) = km1d(k) * l_black(k)
647                ENDIF
648             ENDDO
649
650!
651!--          Add damping layer
652             DO  k = damp_level_ind_1d+1, nzt+1
653                km1d(k) = 1.1 * km1d(k-1)
654                km1d(k) = MIN( km1d(k), 10.0 )
655             ENDDO
656
657!
658!--          Compute the diffusion coefficient for heat via the relationship
659!--          kh = phim / phih * km
660             DO  k = nzb+1, nzt
661                IF ( rif1d(k) >= 0.0 )  THEN
662                   kh1d(k) = km1d(k)
663                ELSE
664                   kh1d(k) = km1d(k) * ( 1.0 - 16.0 * rif1d(k) )**0.25
665                ENDIF
666             ENDDO
667
668          ENDIF   ! .NOT. constant_diffusion
669
670!
671!--       The Runge-Kutta scheme needs the recent diffusion quantities
672          IF ( timestep_scheme(1:5) == 'runge' )  THEN
673             u1d_m  = u1d
674             v1d_m  = v1d
675             IF ( .NOT. constant_diffusion )  THEN
676                e1d_m  = e1d
677                kh1d_m = kh1d
678                km1d_m = km1d
679                l1d_m  = l1d
680                IF ( prandtl_layer )  THEN
681                   usws1d_m = usws1d
682                   vsws1d_m = vsws1d
683                ENDIF
684             ENDIF
685          ENDIF
686
687
688       ENDDO   ! intermediate step loop
689
690!
691!--    Increment simulated time and output times
692       current_timestep_number_1d = current_timestep_number_1d + 1
693       simulated_time_1d          = simulated_time_1d + dt_1d
694       simulated_time_chr         = time_to_string( simulated_time_1d )
695       time_pr_1d                 = time_pr_1d          + dt_1d
696       time_run_control_1d        = time_run_control_1d + dt_1d
697
698!
699!--    Determine and print out quantities for run control
700       IF ( time_run_control_1d >= dt_run_control_1d )  THEN
701          CALL run_control_1d
702          time_run_control_1d = time_run_control_1d - dt_run_control_1d
703       ENDIF
704
705!
706!--    Profile output on file
707       IF ( time_pr_1d >= dt_pr_1d )  THEN
708          CALL print_1d_model
709          time_pr_1d = time_pr_1d - dt_pr_1d
710       ENDIF
711
712!
713!--    Determine size of next time step
714       CALL timestep_1d
715
716    ENDDO   ! time loop
717
718
719 END SUBROUTINE time_integration_1d
720
721
722 SUBROUTINE run_control_1d
723
724!------------------------------------------------------------------------------!
725! Description:
726! ------------
727! Compute and print out quantities for run control of the 1D model.
728!------------------------------------------------------------------------------!
729
730    USE constants
731    USE indices
732    USE model_1d
733    USE pegrid
734    USE control_parameters
735
736    IMPLICIT NONE
737
738    INTEGER ::  k
739    REAL    ::  alpha, energy, umax, uv_total, vmax
740
741!
742!-- Output
743    IF ( myid == 0 )  THEN
744!
745!--    If necessary, write header
746       IF ( .NOT. run_control_header_1d )  THEN
747          WRITE ( 15, 100 )
748          run_control_header_1d = .TRUE.
749       ENDIF
750
751!
752!--    Compute control quantities
753!--    grid level nzp is excluded due to mirror boundary condition
754       umax = 0.0; vmax = 0.0; energy = 0.0
755       DO  k = nzb+1, nzt+1
756          umax = MAX( ABS( umax ), ABS( u1d(k) ) )
757          vmax = MAX( ABS( vmax ), ABS( v1d(k) ) )
758          energy = energy + 0.5 * ( u1d(k)**2 + v1d(k)**2 )
759       ENDDO
760       energy = energy / REAL( nzt - nzb + 1 )
761
762       uv_total = SQRT( u1d(nzb+1)**2 + v1d(nzb+1)**2 )
763       IF ( ABS( v1d(nzb+1) ) .LT. 1.0E-5 )  THEN
764          alpha = ACOS( SIGN( 1.0 , u1d(nzb+1) ) )
765       ELSE
766          alpha = ACOS( u1d(nzb+1) / uv_total )
767          IF ( v1d(nzb+1) <= 0.0 )  alpha = 2.0 * pi - alpha
768       ENDIF
769       alpha = alpha / ( 2.0 * pi ) * 360.0
770
771       WRITE ( 15, 101 )  current_timestep_number_1d, simulated_time_chr, &
772                          dt_1d, umax, vmax, us1d, alpha, energy
773!
774!--    Write buffer contents to disc immediately
775       CALL local_flush( 15 )
776
777    ENDIF
778
779!
780!-- formats
781100 FORMAT (///'1D-Zeitschrittkontrollausgaben:'/ &
782              &'------------------------------'// &
783           &'ITER.  HH:MM:SS    DT      UMAX   VMAX    U*   ALPHA   ENERG.'/ &
784           &'-------------------------------------------------------------')
785101 FORMAT (I5,2X,A9,1X,F6.2,2X,F6.2,1X,F6.2,2X,F5.3,2X,F5.1,2X,F7.2)
786
787
788 END SUBROUTINE run_control_1d
789
790
791
792 SUBROUTINE timestep_1d
793
794!------------------------------------------------------------------------------!
795! Description:
796! ------------
797! Compute the time step w.r.t. the diffusion criterion
798!------------------------------------------------------------------------------!
799
800    USE arrays_3d
801    USE indices
802    USE model_1d
803    USE pegrid
804    USE control_parameters
805
806    IMPLICIT NONE
807
808    INTEGER ::  k
809    REAL    ::  div, dt_diff, fac, percent_change, value
810
811
812!
813!-- Compute the currently feasible time step according to the diffusion
814!-- criterion. At nzb+1 the half grid length is used.
815    IF ( timestep_scheme(1:4) == 'leap' )  THEN
816       fac = 0.25
817    ELSE
818       fac = 0.35
819    ENDIF
820    dt_diff = dt_max_1d
821    DO  k = nzb+2, nzt
822       value   = fac * dzu(k) * dzu(k) / ( km1d(k) + 1E-20 )
823       dt_diff = MIN( value, dt_diff )
824    ENDDO
825    value   = fac * zu(nzb+1) * zu(nzb+1) / ( km1d(nzb+1) + 1E-20 )
826    dt_1d = MIN( value, dt_diff )
827
828!
829!-- Set flag when the time step becomes too small
830    IF ( dt_1d < ( 0.00001 * dt_max_1d ) )  THEN
831       stop_dt_1d = .TRUE.
832       IF ( myid == 0 )  THEN
833          PRINT*,'+++ timestep_1d: timestep has exceeded the lower limit'
834          PRINT*,'                 dt_1d = ',dt_1d,' s   simulation stopped!'
835       ENDIF
836       CALL local_stop
837    ENDIF
838
839    IF ( timestep_scheme(1:4) == 'leap' )  THEN
840
841!
842!--    The current time step will only be changed if the new time step exceeds
843!--    its previous value by 5 % or falls 2 % below. After a time step
844!--    reduction at least 30 iterations must be done with this value before a
845!--    new reduction will be allowed again.
846!--    The control parameters for application of Euler- or leap-frog schemes are
847!--    set accordingly.
848       percent_change = dt_1d / old_dt_1d - 1.0
849       IF ( percent_change > 0.05  .OR.  percent_change < -0.02 )  THEN
850
851!
852!--       Each time step increase is by at most 2 %
853          IF ( percent_change > 0.0  .AND.  simulated_time_1d /= 0.0 )  THEN
854             dt_1d = 1.02 * old_dt_1d
855          ENDIF
856
857!
858!--       A more or less simple new time step value is obtained taking only the
859!--       first two significant digits
860          div = 1000.0
861          DO  WHILE ( dt_1d < div )
862             div = div / 10.0
863          ENDDO
864          dt_1d = NINT( dt_1d * 100.0 / div ) * div / 100.0
865
866!
867!--       Now the time step can be changed.
868          IF ( percent_change < 0.0 )  THEN
869!
870!--          Time step reduction
871             old_dt_1d = dt_1d
872             last_dt_change_1d = current_timestep_number_1d
873          ELSE
874!
875!--          Time step will only be increased if at least 30 iterations have
876!--          been done since the previous time step change, and of course at
877!--          simulation start, respectively.
878             IF ( current_timestep_number_1d >= last_dt_change_1d + 30  .OR. &
879                     simulated_time_1d == 0.0 )  THEN
880                old_dt_1d = dt_1d
881                last_dt_change_1d = current_timestep_number_1d
882             ELSE
883                dt_1d = old_dt_1d
884             ENDIF
885          ENDIF
886       ELSE
887!
888!--       No time step change since the difference is too small
889          dt_1d = old_dt_1d
890       ENDIF
891
892    ELSE    ! Runge-Kutta
893
894!--    A more or less simple new time step value is obtained taking only the
895!--    first two significant digits
896       div = 1000.0
897       DO  WHILE ( dt_1d < div )
898          div = div / 10.0
899       ENDDO
900       dt_1d = NINT( dt_1d * 100.0 / div ) * div / 100.0
901
902       old_dt_1d = dt_1d
903       last_dt_change_1d = current_timestep_number_1d
904
905    ENDIF
906
907 END SUBROUTINE timestep_1d
908
909
910
911 SUBROUTINE print_1d_model
912
913!------------------------------------------------------------------------------!
914! Description:
915! ------------
916! List output of profiles from the 1D-model
917!------------------------------------------------------------------------------!
918
919    USE arrays_3d
920    USE indices
921    USE model_1d
922    USE pegrid
923    USE control_parameters
924
925    IMPLICIT NONE
926
927
928    INTEGER ::  k
929
930
931    IF ( myid == 0 )  THEN
932!
933!--    Open list output file for profiles from the 1D-model
934       CALL check_open( 17 )
935
936!
937!--    Write Header
938       WRITE ( 17, 100 )  TRIM( run_description_header ), &
939                          TRIM( simulated_time_chr )
940       WRITE ( 17, 101 )
941
942!
943!--    Write the values
944       WRITE ( 17, 102 )
945       WRITE ( 17, 101 )
946       DO  k = nzt+1, nzb, -1
947          WRITE ( 17, 103)  k, zu(k), u1d(k), v1d(k), pt_init(k), e1d(k), &
948                            rif1d(k), km1d(k), kh1d(k), l1d(k), zu(k), k
949       ENDDO
950       WRITE ( 17, 101 )
951       WRITE ( 17, 102 )
952       WRITE ( 17, 101 )
953
954!
955!--    Write buffer contents to disc immediately
956       CALL local_flush( 17 )
957
958    ENDIF
959
960!
961!-- Formats
962100 FORMAT (//1X,A/1X,10('-')/' 1d-model profiles'/ &
963            ' Time: ',A)
964101 FORMAT (1X,79('-'))
965102 FORMAT ('   k     zu      u      v     pt      e    rif    Km    Kh     ', &
966            'l      zu      k')
967103 FORMAT (1X,I4,1X,F7.1,1X,F6.2,1X,F6.2,1X,F6.2,1X,F6.2,1X,F5.2,1X,F5.2, &
968            1X,F5.2,1X,F6.2,1X,F7.1,2X,I4)
969
970
971 END SUBROUTINE print_1d_model
Note: See TracBrowser for help on using the repository browser.