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

Last change on this file since 138 was 135, checked in by raasch, 16 years ago

bug concerning mixing length calculation in case of negative f removed; inconsistencies concerning variable grid removed

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