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

Last change on this file since 978 was 978, checked in by fricke, 12 years ago

merge fricke branch back into trunk

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