source: palm/trunk/SOURCE/lpm_droplet_condensation.f90 @ 1354

Last change on this file since 1354 was 1347, checked in by heinze, 10 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 21.8 KB
Line 
1 SUBROUTINE lpm_droplet_condensation
2
3!--------------------------------------------------------------------------------!
4! This file is part of PALM.
5!
6! PALM is free software: you can redistribute it and/or modify it under the terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2014 Leibniz Universitaet Hannover
18!--------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: lpm_droplet_condensation.f90 1347 2014-03-27 13:23:00Z heinze $
27!
28! 1346 2014-03-27 13:18:20Z heinze
29! Bugfix: REAL constants provided with KIND-attribute especially in call of
30! intrinsic function like MAX, MIN, SIGN
31!
32! 1322 2014-03-20 16:38:49Z raasch
33! REAL constants defined as wp-kind
34!
35! 1320 2014-03-20 08:40:49Z raasch
36! ONLY-attribute added to USE-statements,
37! kind-parameters added to all INTEGER and REAL declaration statements,
38! kinds are defined in new module kinds,
39! comment fields (!:) to be used for variable explanations added to
40! all variable declaration statements
41!
42! 1318 2014-03-17 13:35:16Z raasch
43! module interfaces removed
44!
45! 1092 2013-02-02 11:24:22Z raasch
46! unused variables removed
47!
48! 1071 2012-11-29 16:54:55Z franke
49! Ventilation effect for evaporation of large droplets included
50! Check for unreasonable results included in calculation of Rosenbrock method
51! since physically unlikely results were observed and for the same
52! reason the first internal time step in Rosenbrock method should be < 1.0E02 in
53! case of evaporation
54! Unnecessary calculation of ql_int removed
55! Unnecessary calculations in Rosenbrock method (d2rdt2, drdt_m, dt_ros_last)
56! removed
57! Bugfix: factor in calculation of surface tension changed from 0.00155 to
58! 0.000155
59!
60! 1036 2012-10-22 13:43:42Z raasch
61! code put under GPL (PALM 3.9)
62!
63! 849 2012-03-15 10:35:09Z raasch
64! initial revision (former part of advec_particles)
65!
66!
67! Description:
68! ------------
69! Calculates change in droplet radius by condensation/evaporation, using
70! either an analytic formula or by numerically integrating the radius growth
71! equation including curvature and solution effects using Rosenbrocks method
72! (see Numerical recipes in FORTRAN, 2nd edition, p. 731).
73! The analytical formula and growth equation follow those given in
74! Rogers and Yau (A short course in cloud physics, 3rd edition, p. 102/103).
75!------------------------------------------------------------------------------!
76
77    USE arrays_3d,                                                             &
78        ONLY:  hyp, pt, q,  ql_c, ql_v, zu
79
80    USE cloud_parameters,                                                      &
81        ONLY:  bfactor, curvature_solution_effects, diff_coeff_l,              &
82               eps_ros, l_d_rv, l_v, rho_l,  r_v, thermal_conductivity_l
83
84    USE constants,                                                             &
85        ONLY:  pi
86
87    USE control_parameters,                                                    &
88        ONLY:  atmos_ocean_sign, dt_3d, dz, message_string,                    &
89               molecular_viscosity, rho_surface
90    USE cpulog,                                                                &
91        ONLY:  cpu_log, log_point_s
92
93    USE grid_variables,                                                        &
94        ONLY:  dx, ddx, dy, ddy
95
96    USE lpm_collision_kernels_mod,                                             &
97        ONLY:  rclass_lbound, rclass_ubound
98
99    USE kinds
100
101    USE particle_attributes,                                                   &
102        ONLY:  hall_kernel, number_of_particles, offset_ocean_nzt,             &
103               offset_ocean_nzt_m1, particles, radius_classes,                 &
104               use_kernel_tables, wang_kernel
105
106
107    IMPLICIT NONE
108
109    INTEGER(iwp) :: i                          !:
110    INTEGER(iwp) :: internal_timestep_count    !:
111    INTEGER(iwp) :: j                          !:
112    INTEGER(iwp) :: jtry                       !:
113    INTEGER(iwp) :: k                          !:
114    INTEGER(iwp) :: n                          !:
115    INTEGER(iwp) :: ros_count                  !:
116 
117    INTEGER(iwp), PARAMETER ::  maxtry = 40    !:
118
119    LOGICAL ::  repeat                         !:
120
121    REAL(wp) ::  aa                            !:
122    REAL(wp) ::  afactor                       !:
123    REAL(wp) ::  arg                           !:
124    REAL(wp) ::  bb                            !:
125    REAL(wp) ::  cc                            !:
126    REAL(wp) ::  dd                            !:
127    REAL(wp) ::  ddenom                        !:
128    REAL(wp) ::  delta_r                       !:
129    REAL(wp) ::  drdt                          !:
130    REAL(wp) ::  drdt_ini                      !:
131    REAL(wp) ::  dt_ros                        !:
132    REAL(wp) ::  dt_ros_next                   !:
133    REAL(wp) ::  dt_ros_sum                    !:
134    REAL(wp) ::  dt_ros_sum_ini                !:
135    REAL(wp) ::  d2rdtdr                       !:
136    REAL(wp) ::  errmax                        !:
137    REAL(wp) ::  err_ros                       !:
138    REAL(wp) ::  g1                            !:
139    REAL(wp) ::  g2                            !:
140    REAL(wp) ::  g3                            !:
141    REAL(wp) ::  g4                            !:
142    REAL(wp) ::  e_a                           !:
143    REAL(wp) ::  e_s                           !:
144    REAL(wp) ::  gg                            !:
145    REAL(wp) ::  new_r                         !:
146    REAL(wp) ::  p_int                         !:
147    REAL(wp) ::  pt_int                        !:
148    REAL(wp) ::  pt_int_l                      !:
149    REAL(wp) ::  pt_int_u                      !:
150    REAL(wp) ::  q_int                         !:
151    REAL(wp) ::  q_int_l                       !:
152    REAL(wp) ::  q_int_u                       !:
153    REAL(wp) ::  r_ros                         !:
154    REAL(wp) ::  r_ros_ini                     !:
155    REAL(wp) ::  sigma                         !:
156    REAL(wp) ::  t_int                         !:
157    REAL(wp) ::  x                             !:
158    REAL(wp) ::  y                             !:
159    REAL(wp) ::  re_p                          !:
160 
161!-- Parameters for Rosenbrock method
162    REAL(wp), PARAMETER ::  a21 = 2.0             !:
163    REAL(wp), PARAMETER ::  a31 = 48.0/25.0_wp    !:
164    REAL(wp), PARAMETER ::  a32 = 6.0/25.0_wp     !:
165    REAL(wp), PARAMETER ::  b1 = 19.0/9.0_wp      !:
166    REAL(wp), PARAMETER ::  b2 = 0.5              !:
167    REAL(wp), PARAMETER ::  b3 = 25.0/108.0_wp    !:
168    REAL(wp), PARAMETER ::  b4 = 125.0/108.0_wp   !:
169    REAL(wp), PARAMETER ::  c21 = -8.0            !:
170    REAL(wp), PARAMETER ::  c31 = 372.0/25.0_wp   !:
171    REAL(wp), PARAMETER ::  c32 = 12.0/5.0_wp     !:
172    REAL(wp), PARAMETER ::  c41 = -112.0/125.0_wp !:
173    REAL(wp), PARAMETER ::  c42 = -54.0/125.0_wp  !:
174    REAL(wp), PARAMETER ::  c43 = -2.0/5.0_wp     !:
175    REAL(wp), PARAMETER ::  errcon = 0.1296       !:
176    REAL(wp), PARAMETER ::  e1 = 17.0/54.0_wp     !:
177    REAL(wp), PARAMETER ::  e2 = 7.0/36.0_wp      !:
178    REAL(wp), PARAMETER ::  e3 = 0.0              !:
179    REAL(wp), PARAMETER ::  e4 = 125.0/108.0_wp   !:
180    REAL(wp), PARAMETER ::  gam = 0.5             !:
181    REAL(wp), PARAMETER ::  grow = 1.5            !:
182    REAL(wp), PARAMETER ::  pgrow = -0.25         !:
183    REAL(wp), PARAMETER ::  pshrnk = -1.0/3.0_wp  !:
184    REAL(wp), PARAMETER ::  shrnk = 0.5           !:
185
186
187
188    CALL cpu_log( log_point_s(42), 'lpm_droplet_condens', 'start' )
189
190    DO  n = 1, number_of_particles
191!
192!--    Interpolate temperature and humidity.
193!--    First determine left, south, and bottom index of the arrays.
194       i = particles(n)%x * ddx
195       j = particles(n)%y * ddy
196       k = ( particles(n)%z + 0.5 * dz * atmos_ocean_sign ) / dz  &
197           + offset_ocean_nzt                     ! only exact if equidistant
198
199       x  = particles(n)%x - i * dx
200       y  = particles(n)%y - j * dy
201       aa = x**2          + y**2
202       bb = ( dx - x )**2 + y**2
203       cc = x**2          + ( dy - y )**2
204       dd = ( dx - x )**2 + ( dy - y )**2
205       gg = aa + bb + cc + dd
206
207       pt_int_l = ( ( gg - aa ) * pt(k,j,i)   + ( gg - bb ) * pt(k,j,i+1)   &
208                  + ( gg - cc ) * pt(k,j+1,i) + ( gg - dd ) * pt(k,j+1,i+1) &
209                  ) / ( 3.0 * gg )
210
211       pt_int_u = ( ( gg-aa ) * pt(k+1,j,i)   + ( gg-bb ) * pt(k+1,j,i+1)   &
212                  + ( gg-cc ) * pt(k+1,j+1,i) + ( gg-dd ) * pt(k+1,j+1,i+1) &
213                  ) / ( 3.0 * gg )
214
215       pt_int = pt_int_l + ( particles(n)%z - zu(k) ) / dz * &
216                           ( pt_int_u - pt_int_l )
217
218       q_int_l = ( ( gg - aa ) * q(k,j,i)   + ( gg - bb ) * q(k,j,i+1)   &
219                 + ( gg - cc ) * q(k,j+1,i) + ( gg - dd ) * q(k,j+1,i+1) &
220                 ) / ( 3.0 * gg )
221
222       q_int_u = ( ( gg-aa ) * q(k+1,j,i)   + ( gg-bb ) * q(k+1,j,i+1)   &
223                 + ( gg-cc ) * q(k+1,j+1,i) + ( gg-dd ) * q(k+1,j+1,i+1) &
224                 ) / ( 3.0 * gg )
225
226       q_int = q_int_l + ( particles(n)%z - zu(k) ) / dz * &
227                           ( q_int_u - q_int_l )
228
229!
230!--    Calculate real temperature and saturation vapor pressure
231       p_int = hyp(k) + ( particles(n)%z - zu(k) ) / dz * ( hyp(k+1)-hyp(k) )
232       t_int = pt_int * ( p_int / 100000.0 )**0.286
233
234       e_s = 611.0 * EXP( l_d_rv * ( 3.6609E-3 - 1.0 / t_int ) )
235
236!
237!--    Current vapor pressure
238       e_a = q_int * p_int / ( 0.378 * q_int + 0.622 )
239
240!
241!--    Thermal conductivity for water (from Rogers and Yau, Table 7.1),
242!--    diffusivity for water vapor (after Hall und Pruppacher, 1976)
243       thermal_conductivity_l = 7.94048E-05 * t_int + 0.00227011 
244       diff_coeff_l           = 0.211E-4 * ( t_int / 273.15 )**1.94 * &
245                                ( 101325.0 / p_int)
246!
247!--    Change in radius by condensation/evaporation
248       IF ( particles(n)%radius >= 4.0E-5  .AND.  e_a/e_s < 1.0 )  THEN
249!
250!--       Approximation for large radii, where curvature and solution effects
251!--       can be neglected but ventilation effect has to be included in case of
252!--       evaporation.
253!--       First calculate the droplet's Reynolds number
254          re_p = 2.0 * particles(n)%radius * ABS( particles(n)%speed_z ) / &
255                 molecular_viscosity
256!
257!--       Ventilation coefficient after Rogers and Yau, 1989
258          IF ( re_p > 2.5 )  THEN
259             afactor = 0.78 + 0.28 * SQRT( re_p )
260          ELSE
261             afactor = 1.0 + 0.09 * re_p
262          ENDIF
263
264          arg = particles(n)%radius**2 + 2.0 * dt_3d * afactor *              &
265                           ( e_a / e_s - 1.0 ) /                              &
266                           ( ( l_d_rv / t_int - 1.0 ) * l_v * rho_l / t_int / &
267                             thermal_conductivity_l +                         &
268                             rho_l * r_v * t_int / diff_coeff_l / e_s )
269
270          new_r = SQRT( arg )
271
272       ELSEIF ( particles(n)%radius >= 1.0E-6  .OR.  &
273                .NOT. curvature_solution_effects )  THEN
274!
275!--       Approximation for larger radii in case that curvature and solution
276!--       effects are neglected and ventilation effects does not play a role
277          arg = particles(n)%radius**2 + 2.0 * dt_3d *                        &
278                           ( e_a / e_s - 1.0 ) /                              &
279                           ( ( l_d_rv / t_int - 1.0 ) * l_v * rho_l / t_int / &
280                             thermal_conductivity_l +                         &
281                             rho_l * r_v * t_int / diff_coeff_l / e_s )
282          IF ( arg < 1.0E-16 )  THEN
283             new_r = 1.0E-8
284          ELSE
285             new_r = SQRT( arg )
286          ENDIF
287       ENDIF
288
289       IF ( curvature_solution_effects  .AND.                              &
290            ( ( particles(n)%radius < 1.0E-6 ) .OR. ( new_r < 1.0E-6 ) ) ) &
291       THEN
292!
293!--       Curvature and solutions effects are included in growth equation.
294!--       Change in Radius is calculated with 4th-order Rosenbrock method
295!--       for stiff o.d.e's with monitoring local truncation error to adjust
296!--       stepsize (see Numerical recipes in FORTRAN, 2nd edition, p. 731).
297!--       For larger radii the simple analytic method (see ELSE) gives
298!--       almost the same results.
299
300          ros_count = 0
301          repeat = .TRUE.
302!
303!--       Carry out the Rosenbrock algorithm. In case of unreasonable results
304!--       the switch "repeat" will be set true and the algorithm will be carried
305!--       out again with the internal time step set to its initial (small) value.
306!--       Unreasonable results may occur if the external conditions, especially the
307!--       supersaturation, has significantly changed compared to the last PALM
308!--       timestep.
309          DO WHILE ( repeat )
310
311             repeat = .FALSE.
312!
313!--          Surface tension after (Straka, 2009)
314             sigma = 0.0761 - 0.000155 * ( t_int - 273.15 )
315
316             r_ros = particles(n)%radius
317             dt_ros_sum  = 0.0      ! internal integrated time (s)
318             internal_timestep_count = 0
319
320             ddenom  = 1.0 / ( rho_l * r_v * t_int / ( e_s * diff_coeff_l ) +  &
321                               ( l_v / ( r_v * t_int ) - 1.0 ) *               &
322                               rho_l * l_v / ( thermal_conductivity_l * t_int )&
323                             )
324
325             afactor = 2.0 * sigma / ( rho_l * r_v * t_int )
326
327!
328!--          Take internal time step values from the end of last PALM time step
329             dt_ros_next = particles(n)%rvar1
330
331!
332!--          Internal time step should not be > 1.0E-2 in case of evaporation
333!--          because larger values may lead to secondary solutions which are
334!--          physically unlikely
335             IF ( dt_ros_next > 1.0E-2  .AND.  e_a/e_s < 1.0 )  THEN
336                dt_ros_next = 1.0E-3
337             ENDIF
338!
339!--          If calculation of Rosenbrock method is repeated due to unreasonalble
340!--          results during previous try the initial internal time step has to be
341!--          reduced
342             IF ( ros_count > 1 )  THEN
343                dt_ros_next = dt_ros_next - ( 0.2 * dt_ros_next )
344             ELSEIF ( ros_count > 5 )  THEN
345!
346!--             Prevent creation of infinite loop
347                message_string = 'ros_count > 5 in Rosenbrock method'
348                CALL message( 'lpm_droplet_condensation', 'PA0018', 2, 2, &
349                               0, 6, 0 )
350             ENDIF
351
352!
353!--          Internal time step must not be larger than PALM time step
354             dt_ros = MIN( dt_ros_next, dt_3d )
355!
356!--          Integrate growth equation in time unless PALM time step is reached
357             DO WHILE ( dt_ros_sum < dt_3d )
358
359                internal_timestep_count = internal_timestep_count + 1
360
361!
362!--             Derivative at starting value
363                drdt = ddenom / r_ros * ( e_a / e_s - 1.0 - afactor / r_ros + &
364                                          bfactor / r_ros**3 )
365                drdt_ini       = drdt
366                dt_ros_sum_ini = dt_ros_sum
367                r_ros_ini      = r_ros
368
369!
370!--             Calculate radial derivative of dr/dt
371                d2rdtdr = ddenom * ( ( 1.0 - e_a/e_s ) / r_ros**2 + &
372                                     2.0 * afactor / r_ros**3 -     &
373                                     4.0 * bfactor / r_ros**5 )
374!
375!--             Adjust stepsize unless required accuracy is reached
376                DO  jtry = 1, maxtry+1
377
378                   IF ( jtry == maxtry+1 )  THEN
379                      message_string = 'maxtry > 40 in Rosenbrock method'
380                      CALL message( 'lpm_droplet_condensation', 'PA0347', 2, 2, &
381                                    0, 6, 0 )
382                   ENDIF
383
384                   aa    = 1.0 / ( gam * dt_ros ) - d2rdtdr
385                   g1    = drdt_ini / aa
386                   r_ros = r_ros_ini + a21 * g1
387                   drdt  = ddenom / r_ros * ( e_a / e_s - 1.0 - &
388                                              afactor / r_ros + &
389                                              bfactor / r_ros**3 )
390
391                   g2    = ( drdt + c21 * g1 / dt_ros )&
392                           / aa
393                   r_ros = r_ros_ini + a31 * g1 + a32 * g2
394                   drdt  = ddenom / r_ros * ( e_a / e_s - 1.0 - &
395                                              afactor / r_ros + &
396                                              bfactor / r_ros**3 )
397
398                   g3    = ( drdt +  &
399                             ( c31 * g1 + c32 * g2 ) / dt_ros ) / aa
400                   g4    = ( drdt +  &
401                             ( c41 * g1 + c42 * g2 + c43 * g3 ) / dt_ros ) / aa
402                   r_ros = r_ros_ini + b1 * g1 + b2 * g2 + b3 * g3 + b4 * g4
403
404                   dt_ros_sum = dt_ros_sum_ini + dt_ros
405
406                   IF ( dt_ros_sum == dt_ros_sum_ini )  THEN
407                      message_string = 'zero stepsize in Rosenbrock method'
408                      CALL message( 'lpm_droplet_condensation', 'PA0348', 2, 2, &
409                                    0, 6, 0 )
410                   ENDIF
411!
412!--                Calculate error
413                   err_ros = e1*g1 + e2*g2 + e3*g3 + e4*g4
414                   errmax  = 0.0
415                   errmax  = MAX( errmax, ABS( err_ros / r_ros_ini ) ) / eps_ros
416!
417!--                Leave loop if accuracy is sufficient, otherwise try again
418!--                with a reduced stepsize
419                   IF ( errmax <= 1.0 )  THEN
420                      EXIT
421                   ELSE
422                      dt_ros = SIGN( MAX( ABS( 0.9 * dt_ros * errmax**pshrnk ), &
423                                          shrnk * ABS( dt_ros ) ), dt_ros )
424                   ENDIF
425
426                ENDDO  ! loop for stepsize adjustment
427
428!
429!--             Calculate next internal time step
430                IF ( errmax > errcon )  THEN
431                   dt_ros_next = 0.9 * dt_ros * errmax**pgrow
432                ELSE
433                   dt_ros_next = grow * dt_ros
434                ENDIF
435
436!
437!--             Estimated time step is reduced if the PALM time step is exceeded
438                IF ( ( dt_ros_next + dt_ros_sum ) >= dt_3d )  THEN
439                   dt_ros = dt_3d - dt_ros_sum
440                ELSE
441                   dt_ros = dt_ros_next
442                ENDIF
443
444             ENDDO
445!
446!--          Store internal time step value for next PALM step
447             particles(n)%rvar1 = dt_ros_next
448
449             new_r = r_ros
450!
451!--          Radius should not fall below 1E-8 because Rosenbrock method may
452!--          lead to errors otherwise
453             new_r = MAX( new_r, 1.0E-8_wp )
454!
455!--          Check if calculated droplet radius change is reasonable since in
456!--          case of droplet evaporation the Rosenbrock method may lead to
457!--          secondary solutions which are physically unlikely.
458!--          Due to the solution effect the droplets may grow for relative
459!--          humidities below 100%, but change of radius should not be too large.
460!--          In case of unreasonable droplet growth the Rosenbrock method is
461!--          recalculated using a smaller initial time step.
462!--          Limiting values are tested for droplets down to 1.0E-7
463             IF ( new_r - particles(n)%radius >= 3.0E-7  .AND.  &
464                  e_a/e_s < 0.97 )  THEN
465                ros_count = ros_count + 1
466                repeat = .TRUE.
467             ENDIF
468
469          ENDDO    ! Rosenbrock method
470
471       ENDIF
472
473       delta_r = new_r - particles(n)%radius
474
475!
476!--    Sum up the change in volume of liquid water for the respective grid
477!--    volume (this is needed later in lpm_calc_liquid_water_content for
478!--    calculating the release of latent heat)
479       i = ( particles(n)%x + 0.5 * dx ) * ddx
480       j = ( particles(n)%y + 0.5 * dy ) * ddy
481       k = particles(n)%z / dz + 1 + offset_ocean_nzt_m1
482           ! only exact if equidistant
483
484       ql_c(k,j,i) = ql_c(k,j,i) + particles(n)%weight_factor *            &
485                                   rho_l * 1.33333333 * pi *               &
486                                   ( new_r**3 - particles(n)%radius**3 ) / &
487                                   ( rho_surface * dx * dy * dz )
488       IF ( ql_c(k,j,i) > 100.0 )  THEN
489          WRITE( message_string, * ) 'k=',k,' j=',j,' i=',i,      &
490                       ' ql_c=',ql_c(k,j,i), ' &part(',n,')%wf=', &
491                       particles(n)%weight_factor,' delta_r=',delta_r
492          CALL message( 'lpm_droplet_condensation', 'PA0143', 2, 2, -1, 6, 1 )
493       ENDIF
494
495!
496!--    Change the droplet radius
497       IF ( ( new_r - particles(n)%radius ) < 0.0  .AND.  new_r < 0.0 ) &
498       THEN
499          WRITE( message_string, * ) '#1 k=',k,' j=',j,' i=',i,   &
500                       ' e_s=',e_s, ' e_a=',e_a,' t_int=',t_int,  &
501                       ' &delta_r=',delta_r,                      &
502                       ' particle_radius=',particles(n)%radius
503          CALL message( 'lpm_droplet_condensation', 'PA0144', 2, 2, -1, 6, 1 )
504       ENDIF
505
506!
507!--    Sum up the total volume of liquid water (needed below for
508!--    re-calculating the weighting factors)
509       ql_v(k,j,i) = ql_v(k,j,i) + particles(n)%weight_factor * new_r**3
510
511       particles(n)%radius = new_r
512
513!
514!--    Determine radius class of the particle needed for collision
515       IF ( ( hall_kernel  .OR.  wang_kernel )  .AND.  use_kernel_tables ) &
516       THEN
517          particles(n)%class = ( LOG( new_r ) - rclass_lbound ) /  &
518                               ( rclass_ubound - rclass_lbound ) * &
519                               radius_classes
520          particles(n)%class = MIN( particles(n)%class, radius_classes )
521          particles(n)%class = MAX( particles(n)%class, 1 )
522       ENDIF
523
524    ENDDO
525
526    CALL cpu_log( log_point_s(42), 'lpm_droplet_condens', 'stop' )
527
528
529 END SUBROUTINE lpm_droplet_condensation
Note: See TracBrowser for help on using the repository browser.