Ignore:
Timestamp:
Oct 30, 2007 12:12:24 PM (16 years ago)
Author:
letzel
Message:

prognostic_equations include the respective wall_*flux in the parameter list of
calls of diffusion_s. Same as before, only the values of wall_heatflux(0:4)
can be assigned. At present, wall_humidityflux, wall_qflux, wall_salinityflux,
wall_scalarflux are kept zero. diffusion_s uses the respective wall_*flux
instead of wall_heatflux. This update serves two purposes:

  • it avoids errors in calculations with humidity/scalar/salinity and prescribed

non-zero wall_heatflux,

  • it prepares PALM for a possible assignment of wall fluxes of

humidity/scalar/salinity in a future release.

Bugfix: assignment of fluxes at walls

Updates to documentation:
chapter_4.2.html#mode_dvrp
chapter_3.5.4.html#time_series

Default for mrun options -q and -n is "sla3" for lctit. Queues bes1 and bes2
removed. DOC/misc/Tsubame.html updated.

Modified default paths (/work/...) for lctit in .mrun.config.default

File:
1 edited

Legend:

Unmodified
Added
Removed
  • palm/trunk/SOURCE/diffusion_s.f90

    r39 r129  
    44! Actual revisions:
    55! -----------------
    6 !
     6! replace wall_heatflux by wall_s_flux that is now included in the parameter
     7! list, bugfix for assignment of fluxes at walls
    78!
    89! Former revisions:
     
    4546! Call for all grid points
    4647!------------------------------------------------------------------------------!
    47     SUBROUTINE diffusion_s( ddzu, ddzw, kh, s, s_flux_b, s_flux_t, tend )
     48    SUBROUTINE diffusion_s( ddzu, ddzw, kh, s, s_flux_b, s_flux_t, &
     49                            wall_s_flux, tend )
    4850
    4951       USE control_parameters
     
    5759       REAL    ::  ddzu(1:nzt+1), ddzw(1:nzt+1)
    5860       REAL    ::  tend(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1)
     61       REAL    ::  wall_s_flux(0:4)
    5962       REAL, DIMENSION(:,:),   POINTER ::  s_flux_b, s_flux_t
    6063       REAL, DIMENSION(:,:,:), POINTER ::  kh, s
     
    8689                                          + 0.5 * ( fwxp(j,i) *               &
    8790                        ( kh(k,j,i) + kh(k,j,i+1) ) * ( s(k,j,i+1)-s(k,j,i) ) &
    88                         - ( 1.0 - fwxp(j,i) ) * wall_heatflux(1)              &
     91                        + ( 1.0 - fwxp(j,i) ) * wall_s_flux(1)                &
    8992                                                   -fwxm(j,i) *               &
    9093                        ( kh(k,j,i) + kh(k,j,i-1) ) * ( s(k,j,i)-s(k,j,i-1) ) &
    91                         + ( 1.0 - fwxm(j,i) ) * wall_heatflux(3)              &
     94                        + ( 1.0 - fwxm(j,i) ) * wall_s_flux(2)                &
    9295                                                  ) * ddx2                    &
    9396                                          + 0.5 * ( fwyp(j,i) *               &
    9497                        ( kh(k,j,i) + kh(k,j+1,i) ) * ( s(k,j+1,i)-s(k,j,i) ) &
    95                         - ( 1.0 - fwyp(j,i) ) * wall_heatflux(2)              &
     98                        + ( 1.0 - fwyp(j,i) ) * wall_s_flux(3)                &
    9699                                                   -fwym(j,i) *               &
    97100                        ( kh(k,j,i) + kh(k,j-1,i) ) * ( s(k,j,i)-s(k,j-1,i) ) &
    98                         + ( 1.0 - fwym(j,i) ) * wall_heatflux(4)              &
     101                        + ( 1.0 - fwym(j,i) ) * wall_s_flux(4)                &
    99102                                                  ) * ddy2
    100103                ENDDO
     
    156159!------------------------------------------------------------------------------!
    157160    SUBROUTINE diffusion_s_ij( i, j, ddzu, ddzw, kh, s, s_flux_b, s_flux_t, &
    158                                tend )
     161                               wall_s_flux, tend )
    159162
    160163       USE control_parameters
     
    168171       REAL    ::  ddzu(1:nzt+1), ddzw(1:nzt+1)
    169172       REAL    ::  tend(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1)
     173       REAL    ::  wall_s_flux(0:4)
    170174       REAL, DIMENSION(:,:),   POINTER ::  s_flux_b, s_flux_t
    171175       REAL, DIMENSION(:,:,:), POINTER ::  kh, s
     
    195199                                          + 0.5 * ( fwxp(j,i) *               &
    196200                        ( kh(k,j,i) + kh(k,j,i+1) ) * ( s(k,j,i+1)-s(k,j,i) ) &
    197                         - ( 1.0 - fwxp(j,i) ) * wall_heatflux(1)              &
     201                        + ( 1.0 - fwxp(j,i) ) * wall_s_flux(1)                &
    198202                                                   -fwxm(j,i) *               &
    199203                        ( kh(k,j,i) + kh(k,j,i-1) ) * ( s(k,j,i)-s(k,j,i-1) ) &
    200                         + ( 1.0 - fwxm(j,i) ) * wall_heatflux(3)              &
     204                        + ( 1.0 - fwxm(j,i) ) * wall_s_flux(2)                &
    201205                                                  ) * ddx2                    &
    202206                                          + 0.5 * ( fwyp(j,i) *               &
    203207                        ( kh(k,j,i) + kh(k,j+1,i) ) * ( s(k,j+1,i)-s(k,j,i) ) &
    204                         - ( 1.0 - fwyp(j,i) ) * wall_heatflux(2)              &
     208                        + ( 1.0 - fwyp(j,i) ) * wall_s_flux(3)                &
    205209                                                   -fwym(j,i) *               &
    206210                        ( kh(k,j,i) + kh(k,j-1,i) ) * ( s(k,j,i)-s(k,j-1,i) ) &
    207                         + ( 1.0 - fwym(j,i) ) * wall_heatflux(4)              &
     211                        + ( 1.0 - fwym(j,i) ) * wall_s_flux(4)                &
    208212                                                  ) * ddy2
    209213          ENDDO
Note: See TracChangeset for help on using the changeset viewer.