Changeset 1019 for palm/trunk


Ignore:
Timestamp:
Sep 28, 2012 6:46:45 AM (12 years ago)
Author:
raasch
Message:

subroutine prognostic_equations_noopt has been removed

Location:
palm/trunk/SOURCE
Files:
4 edited

Legend:

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

    r1017 r1019  
    194194
    195195!
    196 !--    Arrays needed for reasons of speed optimization for cache and noopt
    197 !--    version. For the vector version the buffer arrays are not necessary,
     196!--    Arrays needed for reasons of speed optimization for cache version.
     197!--    For the vector version the buffer arrays are not necessary,
    198198!--    because the the fluxes can swapped directly inside the loops of the
    199199!--    advection routines.
  • palm/trunk/SOURCE/check_parameters.f90

    r1017 r1019  
    44! Current revisions:
    55! -----------------
    6 !
     6! non-optimized version of prognostic_equations not allowed any more
    77!
    88! Former revisions:
     
    488488    SELECT CASE ( TRIM( loop_optimization ) )
    489489
    490        CASE ( 'acc', 'cache', 'noopt', 'vector' )
     490       CASE ( 'acc', 'cache', 'vector' )
    491491          CONTINUE
    492492
     
    646646!
    647647!-- Advection schemes:
    648 !       
    649 !-- Set the LOGICALS to enhance the performance.
    650     IF ( momentum_advec == 'ws-scheme' )    ws_scheme_mom = .TRUE.
    651     IF ( scalar_advec   == 'ws-scheme'   )  ws_scheme_sca = .TRUE.
    652    
    653648    IF ( momentum_advec /= 'pw-scheme'  .AND.  momentum_advec /= 'ws-scheme' ) &
    654649    THEN
     
    673668       CALL message( 'check_parameters', 'PA0024', 1, 2, 0, 6, 0 )
    674669    ENDIF
     670    IF ( scalar_advec == 'bc-scheme'  .AND.  loop_optimization == 'cache' ) &
     671    THEN
     672       message_string = 'advection_scheme scalar_advec = "' &
     673         // TRIM( scalar_advec ) // '" not implemented for & loop_optimization = "' // &
     674         TRIM( loop_optimization ) // '"'
     675       CALL message( 'check_parameters', 'PA0026', 1, 2, 0, 6, 0 )
     676    ENDIF
    675677
    676678    IF ( use_sgs_for_particles  .AND.  .NOT. use_upstream_for_tke )  THEN
     
    686688       CALL message( 'check_parameters', 'PA0349', 1, 2, 0, 6, 0 )
    687689    ENDIF
     690
     691!
     692!-- Set LOGICAL switches to enhance performance
     693    IF ( momentum_advec == 'ws-scheme' )    ws_scheme_mom = .TRUE.
     694    IF ( scalar_advec   == 'ws-scheme'   )  ws_scheme_sca = .TRUE.
    688695
    689696!
  • palm/trunk/SOURCE/prognostic_equations.f90

    r1017 r1019  
    44! Current revisions:
    55! -----------------
    6 !
     6! non-optimized version of prognostic_equations removed
    77!
    88! Former revisions:
     
    139139
    140140    PRIVATE
    141     PUBLIC prognostic_equations_noopt, prognostic_equations_cache, &
    142            prognostic_equations_vector, prognostic_equations_acc
    143 
    144     INTERFACE prognostic_equations_noopt
    145        MODULE PROCEDURE prognostic_equations_noopt
    146     END INTERFACE prognostic_equations_noopt
     141    PUBLIC prognostic_equations_cache, prognostic_equations_vector, &
     142           prognostic_equations_acc
    147143
    148144    INTERFACE prognostic_equations_cache
     
    160156
    161157 CONTAINS
    162 
    163 
    164  SUBROUTINE prognostic_equations_noopt
    165 
    166 !------------------------------------------------------------------------------!
    167 ! Version with single loop optimization
    168 !
    169 ! (Optimized over each single prognostic equation.)
    170 !------------------------------------------------------------------------------!
    171 
    172     IMPLICIT NONE
    173 
    174     CHARACTER (LEN=9) ::  time_to_string
    175     INTEGER ::  i, i_omp_start, j, k, tn = 0
    176     REAL    ::  sbt
    177 
    178 !
    179 !-- Calculate those variables needed in the tendency terms which need
    180 !-- global communication
    181     IF ( .NOT. neutral )  CALL calc_mean_profile( pt, 4 )
    182     IF ( ocean         )  CALL calc_mean_profile( rho, 64 )
    183     IF ( humidity      )  CALL calc_mean_profile( vpt, 44 )
    184     IF ( ( ws_scheme_mom .OR. ws_scheme_sca )  .AND.  &
    185          intermediate_timestep_count == 1 )  CALL ws_statistics
    186 
    187 !
    188 !-- u-velocity component
    189     CALL cpu_log( log_point(5), 'u-equation', 'start' )
    190 
    191     i_omp_start = nxlu
    192     DO  i = nxlu, nxr
    193        DO  j = nys, nyn
    194 !
    195 !--       Tendency terms
    196           tend(:,j,i) = 0.0
    197           IF ( timestep_scheme(1:5) == 'runge' )  THEN
    198              IF ( ws_scheme_mom )  THEN
    199                  CALL advec_u_ws( i, j, i_omp_start, tn )
    200              ELSE
    201                  CALL advec_u_pw( i, j )
    202              ENDIF
    203           ELSE
    204              CALL advec_u_up( i, j )
    205           ENDIF
    206           CALL diffusion_u( i, j )
    207           CALL coriolis( i, j, 1 )
    208           IF ( sloping_surface  .AND.  .NOT. neutral )  THEN
    209              CALL buoyancy( i, j, pt, pt_reference, 1, 4 )
    210           ENDIF
    211 
    212 !
    213 !--       Drag by plant canopy
    214           IF ( plant_canopy )  CALL plant_canopy_model( i, j, 1 )
    215 
    216 !
    217 !--       External pressure gradient
    218           IF ( dp_external )  THEN
    219              DO  k = dp_level_ind_b+1, nzt
    220                 tend(k,j,i) = tend(k,j,i) - dpdxy(1) * dp_smooth_factor(k)
    221              ENDDO
    222           ENDIF
    223 
    224           CALL user_actions( i, j, 'u-tendency' )
    225 
    226 !
    227 !--       Prognostic equation for u-velocity component
    228           DO  k = nzb_u_inner(j,i)+1, nzt
    229              u_p(k,j,i) = u(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) +          &
    230                                                tsc(3) * tu_m(k,j,i) )          &
    231                                    - tsc(5) * rdf(k) * ( u(k,j,i) - ug(k) )
    232           ENDDO
    233 
    234 !
    235 !--       Calculate tendencies for the next Runge-Kutta step
    236           IF ( timestep_scheme(1:5) == 'runge' )  THEN
    237              IF ( intermediate_timestep_count == 1 )  THEN
    238                 DO  k = nzb_u_inner(j,i)+1, nzt
    239                    tu_m(k,j,i) = tend(k,j,i)
    240                 ENDDO
    241              ELSEIF ( intermediate_timestep_count < &
    242                       intermediate_timestep_count_max )  THEN
    243                 DO  k = nzb_u_inner(j,i)+1, nzt
    244                    tu_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tu_m(k,j,i)
    245                 ENDDO
    246              ENDIF
    247           ENDIF
    248 
    249        ENDDO
    250     ENDDO
    251 
    252     CALL cpu_log( log_point(5), 'u-equation', 'stop' )
    253 
    254 !
    255 !-- v-velocity component
    256     CALL cpu_log( log_point(6), 'v-equation', 'start' )
    257 
    258     i_omp_start = nxl
    259     DO  i = nxl, nxr
    260        DO  j = nysv, nyn
    261 !
    262 !--       Tendency terms
    263           tend(:,j,i) = 0.0
    264           IF ( timestep_scheme(1:5) == 'runge' )  THEN
    265              IF ( ws_scheme_mom )  THEN
    266                  CALL advec_v_ws( i, j, i_omp_start, tn )
    267              ELSE
    268                  CALL advec_v_pw( i, j )
    269              ENDIF
    270 
    271           ELSE
    272              CALL advec_v_up( i, j )
    273           ENDIF
    274           CALL diffusion_v( i, j )
    275           CALL coriolis( i, j, 2 )
    276 
    277 !
    278 !--       Drag by plant canopy
    279           IF ( plant_canopy )  CALL plant_canopy_model( i, j, 2 )     
    280 
    281 !
    282 !--       External pressure gradient
    283           IF ( dp_external )  THEN
    284              DO  k = dp_level_ind_b+1, nzt
    285                 tend(k,j,i) = tend(k,j,i) - dpdxy(2) * dp_smooth_factor(k)
    286              ENDDO
    287           ENDIF
    288 
    289           CALL user_actions( i, j, 'v-tendency' )
    290 
    291 !
    292 !--       Prognostic equation for v-velocity component
    293           DO  k = nzb_v_inner(j,i)+1, nzt
    294              v_p(k,j,i) = v(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) +          &
    295                                                tsc(3) * tv_m(k,j,i) )          &
    296                                    - tsc(5) * rdf(k) * ( v(k,j,i) - vg(k) )
    297           ENDDO
    298 
    299 !
    300 !--       Calculate tendencies for the next Runge-Kutta step
    301           IF ( timestep_scheme(1:5) == 'runge' )  THEN
    302              IF ( intermediate_timestep_count == 1 )  THEN
    303                 DO  k = nzb_v_inner(j,i)+1, nzt
    304                    tv_m(k,j,i) = tend(k,j,i)
    305                 ENDDO
    306              ELSEIF ( intermediate_timestep_count < &
    307                       intermediate_timestep_count_max )  THEN
    308                 DO  k = nzb_v_inner(j,i)+1, nzt
    309                    tv_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tv_m(k,j,i)
    310                 ENDDO
    311              ENDIF
    312           ENDIF
    313 
    314        ENDDO
    315     ENDDO
    316 
    317     CALL cpu_log( log_point(6), 'v-equation', 'stop' )
    318 
    319 !
    320 !-- w-velocity component
    321     CALL cpu_log( log_point(7), 'w-equation', 'start' )
    322 
    323     DO  i = nxl, nxr
    324        DO  j = nys, nyn
    325 !
    326 !--       Tendency terms
    327           tend(:,j,i) = 0.0
    328           IF ( timestep_scheme(1:5) == 'runge' )  THEN
    329              IF ( ws_scheme_mom )  THEN
    330                  CALL advec_w_ws( i, j, i_omp_start, tn )
    331              ELSE 
    332                  CALL advec_w_pw( i, j )
    333              ENDIF
    334 
    335           ELSE
    336              CALL advec_w_up( i, j )
    337           ENDIF
    338           CALL diffusion_w( i, j )
    339           CALL coriolis( i, j, 3 )
    340 
    341           IF ( .NOT. neutral )  THEN
    342              IF ( ocean )  THEN
    343                 CALL buoyancy( i, j, rho, rho_reference, 3, 64 )
    344              ELSE
    345                 IF ( .NOT. humidity )  THEN
    346                    CALL buoyancy( i, j, pt, pt_reference, 3, 4 )
    347                 ELSE
    348                    CALL buoyancy( i, j, vpt, pt_reference, 3, 44 )
    349                 ENDIF
    350              ENDIF
    351           ENDIF
    352 
    353 !
    354 !--       Drag by plant canopy
    355           IF ( plant_canopy )  CALL plant_canopy_model( i, j, 3 )
    356 
    357           CALL user_actions( i, j, 'w-tendency' )
    358 
    359 !
    360 !--       Prognostic equation for w-velocity component
    361           DO  k = nzb_w_inner(j,i)+1, nzt-1
    362              w_p(k,j,i) = w(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) +          &
    363                                                tsc(3) * tw_m(k,j,i) )          &
    364                                    - tsc(5) * rdf(k) * w(k,j,i)
    365           ENDDO
    366 
    367 !
    368 !--       Calculate tendencies for the next Runge-Kutta step
    369           IF ( timestep_scheme(1:5) == 'runge' )  THEN
    370              IF ( intermediate_timestep_count == 1 )  THEN
    371                 DO  k = nzb_w_inner(j,i)+1, nzt-1
    372                    tw_m(k,j,i) = tend(k,j,i)
    373                 ENDDO
    374              ELSEIF ( intermediate_timestep_count < &
    375                       intermediate_timestep_count_max )  THEN
    376                 DO  k = nzb_w_inner(j,i)+1, nzt-1
    377                    tw_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tw_m(k,j,i)
    378                 ENDDO
    379              ENDIF
    380           ENDIF
    381 
    382        ENDDO
    383     ENDDO
    384 
    385     CALL cpu_log( log_point(7), 'w-equation', 'stop' )
    386 
    387 !
    388 !-- If required, compute prognostic equation for potential temperature
    389     IF ( .NOT. neutral )  THEN
    390 
    391        CALL cpu_log( log_point(13), 'pt-equation', 'start' )
    392 
    393 !
    394 !--    pt-tendency terms with communication
    395        sbt = tsc(2)
    396        IF ( scalar_advec == 'bc-scheme' )  THEN
    397 
    398           IF ( timestep_scheme(1:5) /= 'runge' )  THEN
    399 !
    400 !--          Bott-Chlond scheme always uses Euler time step. Thus:
    401              sbt = 1.0
    402           ENDIF
    403           tend = 0.0
    404           CALL advec_s_bc( pt, 'pt' )
    405 
    406        ENDIF
    407 
    408 !
    409 !--    pt-tendency terms with no communication
    410        DO  i = nxl, nxr
    411           DO  j = nys, nyn
    412 !
    413 !--          Tendency terms
    414              IF ( scalar_advec /= 'bc-scheme' )  THEN
    415                 tend(:,j,i) = 0.0
    416                 IF ( timestep_scheme(1:5) == 'runge' )  THEN
    417                    IF ( ws_scheme_sca )  THEN
    418                       CALL advec_s_ws( i, j, pt, 'pt', flux_s_pt,       &
    419                                        diss_s_pt, flux_l_pt, diss_l_pt, &
    420                                        i_omp_start, tn )
    421                    ELSE
    422                        CALL advec_s_pw( i, j, pt )
    423                    ENDIF
    424                 ELSE
    425                    CALL advec_s_up( i, j, pt )
    426                 ENDIF
    427              ENDIF
    428 
    429              CALL diffusion_s( i, j, pt, shf, tswst, wall_heatflux )
    430 
    431 !
    432 !--          If required compute heating/cooling due to long wave radiation
    433 !--          processes
    434              IF ( radiation )  THEN
    435                 CALL calc_radiation( i, j )
    436              ENDIF
    437 
    438 !
    439 !--          If required compute impact of latent heat due to precipitation
    440              IF ( precipitation )  THEN
    441                 CALL impact_of_latent_heat( i, j )
    442              ENDIF
    443 
    444 !
    445 !--          Consideration of heat sources within the plant canopy
    446              IF ( plant_canopy .AND. ( cthf /= 0.0 ) ) THEN
    447                 CALL plant_canopy_model( i, j, 4 )
    448              ENDIF
    449 
    450 !
    451 !--          If required compute influence of large-scale subsidence/ascent
    452              IF ( large_scale_subsidence )  THEN
    453                 CALL subsidence( i, j, tend, pt, pt_init )
    454              ENDIF
    455 
    456              CALL user_actions( i, j, 'pt-tendency' )
    457 
    458 !
    459 !--          Prognostic equation for potential temperature
    460              DO  k = nzb_s_inner(j,i)+1, nzt
    461                 pt_p(k,j,i) = pt(k,j,i) + dt_3d * ( sbt * tend(k,j,i) +        &
    462                                                     tsc(3) * tpt_m(k,j,i) )    &
    463                                         - tsc(5) * ( pt(k,j,i) - pt_init(k) ) *&
    464                                           ( rdf_sc(k) + ptdf_x(i) + ptdf_y(j) )
    465              ENDDO
    466 
    467 !
    468 !--          Calculate tendencies for the next Runge-Kutta step
    469              IF ( timestep_scheme(1:5) == 'runge' )  THEN
    470                 IF ( intermediate_timestep_count == 1 )  THEN
    471                    DO  k = nzb_s_inner(j,i)+1, nzt
    472                       tpt_m(k,j,i) = tend(k,j,i)
    473                    ENDDO
    474                 ELSEIF ( intermediate_timestep_count < &
    475                          intermediate_timestep_count_max )  THEN
    476                    DO  k = nzb_s_inner(j,i)+1, nzt
    477                       tpt_m(k,j,i) = -9.5625 * tend(k,j,i) + &
    478                                      5.3125 * tpt_m(k,j,i)
    479                    ENDDO
    480                 ENDIF
    481              ENDIF
    482 
    483           ENDDO
    484        ENDDO
    485 
    486        CALL cpu_log( log_point(13), 'pt-equation', 'stop' )
    487 
    488     ENDIF
    489 
    490 !
    491 !-- If required, compute prognostic equation for salinity
    492     IF ( ocean )  THEN
    493 
    494        CALL cpu_log( log_point(37), 'sa-equation', 'start' )
    495 
    496 !
    497 !--    sa-tendency terms with communication
    498        sbt = tsc(2)
    499        IF ( scalar_advec == 'bc-scheme' )  THEN
    500 
    501           IF ( timestep_scheme(1:5) /= 'runge' )  THEN
    502 !
    503 !--          Bott-Chlond scheme always uses Euler time step. Thus:
    504              sbt = 1.0
    505           ENDIF
    506           tend = 0.0
    507           CALL advec_s_bc( sa, 'sa' )
    508 
    509        ENDIF
    510 
    511 !
    512 !--    sa terms with no communication
    513        DO  i = nxl, nxr
    514           DO  j = nys, nyn
    515 !
    516 !--          Tendency-terms
    517              IF ( scalar_advec /= 'bc-scheme' )  THEN
    518                 tend(:,j,i) = 0.0
    519                 IF ( timestep_scheme(1:5) == 'runge' ) THEN
    520                    IF ( ws_scheme_sca )  THEN
    521                        CALL advec_s_ws( i, j, sa, 'sa', flux_s_sa,  &
    522                                     diss_s_sa, flux_l_sa, diss_l_sa, i_omp_start, tn )
    523                    ELSE
    524                        CALL advec_s_pw( i, j, sa )
    525                    ENDIF
    526 
    527                 ELSE
    528                    CALL advec_s_up( i, j, sa )
    529                 ENDIF
    530              ENDIF
    531 
    532              CALL diffusion_s( i, j, sa, saswsb, saswst, wall_salinityflux )
    533        
    534              CALL user_actions( i, j, 'sa-tendency' )
    535 
    536 !
    537 !--          Prognostic equation for salinity
    538              DO  k = nzb_s_inner(j,i)+1, nzt
    539                 sa_p(k,j,i) = sa(k,j,i) + dt_3d * ( sbt * tend(k,j,i) +        &
    540                                                     tsc(3) * tsa_m(k,j,i) )    &
    541                                         - tsc(5) * rdf_sc(k) *                 &
    542                                           ( sa(k,j,i) - sa_init(k) )
    543                 IF ( sa_p(k,j,i) < 0.0 )  sa_p(k,j,i) = 0.1 * sa(k,j,i)
    544              ENDDO
    545 
    546 !
    547 !--          Calculate tendencies for the next Runge-Kutta step
    548              IF ( timestep_scheme(1:5) == 'runge' )  THEN
    549                 IF ( intermediate_timestep_count == 1 )  THEN
    550                    DO  k = nzb_s_inner(j,i)+1, nzt
    551                       tsa_m(k,j,i) = tend(k,j,i)
    552                    ENDDO
    553                 ELSEIF ( intermediate_timestep_count < &
    554                          intermediate_timestep_count_max )  THEN
    555                    DO  k = nzb_s_inner(j,i)+1, nzt
    556                       tsa_m(k,j,i) = -9.5625 * tend(k,j,i) + &
    557                                       5.3125 * tsa_m(k,j,i)
    558                    ENDDO
    559                 ENDIF
    560              ENDIF
    561 
    562 !
    563 !--          Calculate density by the equation of state for seawater
    564              CALL eqn_state_seawater( i, j )
    565 
    566           ENDDO
    567        ENDDO
    568 
    569        CALL cpu_log( log_point(37), 'sa-equation', 'stop' )
    570 
    571     ENDIF
    572 
    573 !
    574 !-- If required, compute prognostic equation for total water content / scalar
    575     IF ( humidity  .OR.  passive_scalar )  THEN
    576 
    577        CALL cpu_log( log_point(29), 'q/s-equation', 'start' )
    578 
    579 !
    580 !--    Scalar/q-tendency terms with communication
    581        sbt = tsc(2)
    582        IF ( scalar_advec == 'bc-scheme' )  THEN
    583 
    584           IF ( timestep_scheme(1:5) /= 'runge' )  THEN
    585 !
    586 !--          Bott-Chlond scheme always uses Euler time step. Thus:
    587              sbt = 1.0
    588           ENDIF
    589           tend = 0.0
    590           CALL advec_s_bc( q, 'q' )
    591 
    592        ENDIF
    593 
    594 !
    595 !--    Scalar/q-tendency terms with no communication
    596        DO  i = nxl, nxr
    597           DO  j = nys, nyn
    598 !
    599 !--          Tendency-terms
    600              IF ( scalar_advec /= 'bc-scheme' )  THEN
    601                 tend(:,j,i) = 0.0
    602                 IF ( timestep_scheme(1:5) == 'runge' ) THEN
    603                    IF ( ws_scheme_sca )  THEN
    604                        CALL advec_s_ws( i, j, q, 'q', flux_s_q, &
    605                                    diss_s_q, flux_l_q, diss_l_q, i_omp_start, tn )
    606                    ELSE
    607                        CALL advec_s_pw( i, j, q )
    608                    ENDIF
    609                 ELSE
    610                    CALL advec_s_up( i, j, q )
    611                 ENDIF
    612              ENDIF
    613 
    614              CALL diffusion_s( i, j, q, qsws, qswst, wall_qflux )
    615        
    616 !
    617 !--          If required compute decrease of total water content due to
    618 !--          precipitation
    619              IF ( precipitation )  THEN
    620                 CALL calc_precipitation( i, j )
    621              ENDIF
    622 
    623 !
    624 !--          Sink or source of scalar concentration due to canopy elements
    625              IF ( plant_canopy ) CALL plant_canopy_model( i, j, 5 )
    626              
    627 !
    628 !--          If required compute influence of large-scale subsidence/ascent
    629              IF ( large_scale_subsidence )  THEN
    630                 CALL subsidence( i, j, tend, q, q_init )
    631              ENDIF
    632 
    633              CALL user_actions( i, j, 'q-tendency' )
    634 
    635 !
    636 !--          Prognostic equation for total water content / scalar
    637              DO  k = nzb_s_inner(j,i)+1, nzt
    638                 q_p(k,j,i) = q(k,j,i) + dt_3d * ( sbt * tend(k,j,i) +          &
    639                                                   tsc(3) * tq_m(k,j,i) )       &
    640                                       - tsc(5) * rdf_sc(k) *                   &
    641                                         ( q(k,j,i) - q_init(k) )
    642                 IF ( q_p(k,j,i) < 0.0 )  q_p(k,j,i) = 0.1 * q(k,j,i)
    643              ENDDO
    644 
    645 !
    646 !--          Calculate tendencies for the next Runge-Kutta step
    647              IF ( timestep_scheme(1:5) == 'runge' )  THEN
    648                 IF ( intermediate_timestep_count == 1 )  THEN
    649                    DO  k = nzb_s_inner(j,i)+1, nzt
    650                       tq_m(k,j,i) = tend(k,j,i)
    651                    ENDDO
    652                 ELSEIF ( intermediate_timestep_count < &
    653                          intermediate_timestep_count_max )  THEN
    654                    DO  k = nzb_s_inner(j,i)+1, nzt
    655                       tq_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tq_m(k,j,i)
    656                    ENDDO
    657                 ENDIF
    658              ENDIF
    659 
    660           ENDDO
    661        ENDDO
    662 
    663        CALL cpu_log( log_point(29), 'q/s-equation', 'stop' )
    664 
    665     ENDIF
    666 
    667 !
    668 !-- If required, compute prognostic equation for turbulent kinetic
    669 !-- energy (TKE)
    670     IF ( .NOT. constant_diffusion )  THEN
    671 
    672        CALL cpu_log( log_point(16), 'tke-equation', 'start' )
    673 
    674 !
    675 !--    TKE-tendency terms with communication
    676        CALL production_e_init
    677 
    678        sbt = tsc(2)
    679        IF ( .NOT. use_upstream_for_tke )  THEN
    680           IF ( scalar_advec == 'bc-scheme' )  THEN
    681 
    682              IF ( timestep_scheme(1:5) /= 'runge' )  THEN
    683 !
    684 !--             Bott-Chlond scheme always uses Euler time step. Thus:
    685                 sbt = 1.0
    686              ENDIF
    687              tend = 0.0
    688              CALL advec_s_bc( e, 'e' )
    689           ENDIF
    690        ENDIF
    691 
    692 !
    693 !--    TKE-tendency terms with no communication
    694        DO  i = nxl, nxr
    695           DO  j = nys, nyn
    696 !
    697 !--          Tendency-terms
    698              IF ( scalar_advec /= 'bc-scheme' .OR. use_upstream_for_tke )  THEN
    699                 IF ( use_upstream_for_tke )  THEN
    700                    tend(:,j,i) = 0.0
    701                    CALL advec_s_up( i, j, e )
    702                 ELSE
    703                    tend(:,j,i) = 0.0
    704                    IF ( timestep_scheme(1:5) == 'runge' )  THEN
    705                       IF ( ws_scheme_sca )  THEN
    706                           CALL advec_s_ws( i, j, e, 'e', flux_s_e, &
    707                                       diss_s_e, flux_l_e, diss_l_e, i_omp_start, tn )
    708                       ELSE
    709                           CALL advec_s_pw( i, j, e )
    710                       ENDIF
    711                    ELSE
    712                       CALL advec_s_up( i, j, e )
    713                    ENDIF
    714                 ENDIF
    715              ENDIF
    716 
    717              IF ( .NOT. humidity )  THEN
    718                 IF ( ocean )  THEN
    719                    CALL diffusion_e( i, j, prho, prho_reference )
    720                 ELSE
    721                    CALL diffusion_e( i, j, pt, pt_reference )
    722                 ENDIF
    723              ELSE
    724                 CALL diffusion_e( i, j, vpt, pt_reference )
    725              ENDIF
    726 
    727              CALL production_e( i, j )
    728 
    729 !
    730 !--          Additional sink term for flows through plant canopies
    731              IF ( plant_canopy )  CALL plant_canopy_model( i, j, 6 )         
    732 
    733              CALL user_actions( i, j, 'e-tendency' )
    734 
    735 !
    736 !--          Prognostic equation for TKE.
    737 !--          Eliminate negative TKE values, which can occur due to numerical
    738 !--          reasons in the course of the integration. In such cases the old TKE
    739 !--          value is reduced by 90%.
    740              DO  k = nzb_s_inner(j,i)+1, nzt
    741                 e_p(k,j,i) = e(k,j,i) + dt_3d * ( sbt * tend(k,j,i) +          &
    742                                                   tsc(3) * te_m(k,j,i) )
    743                 IF ( e_p(k,j,i) < 0.0 )  e_p(k,j,i) = 0.1 * e(k,j,i)
    744              ENDDO
    745 
    746 !
    747 !--          Calculate tendencies for the next Runge-Kutta step
    748              IF ( timestep_scheme(1:5) == 'runge' )  THEN
    749                 IF ( intermediate_timestep_count == 1 )  THEN
    750                    DO  k = nzb_s_inner(j,i)+1, nzt
    751                       te_m(k,j,i) = tend(k,j,i)
    752                    ENDDO
    753                 ELSEIF ( intermediate_timestep_count < &
    754                          intermediate_timestep_count_max )  THEN
    755                    DO  k = nzb_s_inner(j,i)+1, nzt
    756                       te_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * te_m(k,j,i)
    757                    ENDDO
    758                 ENDIF
    759              ENDIF
    760 
    761           ENDDO
    762        ENDDO
    763 
    764        CALL cpu_log( log_point(16), 'tke-equation', 'stop' )
    765 
    766     ENDIF
    767 
    768 
    769  END SUBROUTINE prognostic_equations_noopt
    770158
    771159
  • palm/trunk/SOURCE/time_integration.f90

    r1017 r1019  
    44! Current revisions:
    55! -----------------
    6 !
     6! non-optimized version of prognostic_equations removed
    77!
    88! Former revisions:
     
    177177!--       in the other versions a good vectorization is prohibited due to
    178178!--       inlining problems.
    179           IF ( loop_optimization == 'vector' )  THEN
     179          IF ( loop_optimization == 'cache' )  THEN
     180             CALL prognostic_equations_cache
     181          ELSEIF ( loop_optimization == 'vector' )  THEN
    180182             CALL prognostic_equations_vector
    181183          ELSEIF ( loop_optimization == 'acc' )  THEN
    182184             CALL prognostic_equations_acc
    183           ELSE
    184              IF ( scalar_advec == 'bc-scheme' )  THEN
    185                 CALL prognostic_equations_noopt
    186              ELSE
    187                 CALL prognostic_equations_cache
    188              ENDIF
    189185          ENDIF
    190186
Note: See TracChangeset for help on using the changeset viewer.