New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 497 for trunk/NEMO/OFF_SRC/LDF/ldfslp.F90 – NEMO

Ignore:
Timestamp:
2006-09-12T13:03:53+02:00 (18 years ago)
Author:
opalod
Message:

nemo_v1_update_064:CE:re-organization of coordinate definition and scale factors

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OFF_SRC/LDF/ldfslp.F90

    r343 r497  
    2020   USE lbclnk          ! ocean lateral boundary conditions (or mpp link) 
    2121   USE in_out_manager  ! I/O manager 
     22   USE prtctl          ! Print control 
    2223 
    2324   IMPLICIT NONE 
     
    4647#  include "vectopt_loop_substitute.h90" 
    4748   !!---------------------------------------------------------------------- 
    48    !!   OPA 9.0 , LOCEAN-IPSL  (2005) 
    49    !!   $Header$ 
    50    !!   This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt 
     49   !!   OPA 9.0 , LOCEAN-IPSL (2005)  
    5150   !!---------------------------------------------------------------------- 
    5251 
     
    6968      !!      of 10cm/s) 
    7069      !!        A horizontal shapiro filter is applied to the slopes 
    71       !!        'key_s_coord' defined: add to the previously computed slopes 
     70      !!        ln_sco=T, s-coordinate, add to the previously computed slopes 
    7271      !!      the slope of the model level surface. 
    7372      !!        macro-tasked on horizontal slab (jk-loop)  (2, jpk-1) 
    7473      !!      [slopes already set to zero at level 1, and to zero or the ocean 
    75       !!      bottom slope ('key_s_coord' defined) at level jpk in inildf] 
     74      !!      bottom slope (ln_sco=T) at level jpk in inildf] 
    7675      !! 
    7776      !! ** Action : - uslp, wslpi, and vslp, wslpj, the i- and  j-slopes  
     
    8382      !!   8.1  !  99-10  (A. Jouzeau)  NEW profile 
    8483      !!   8.5  !  99-10  (G. Madec)  Free form, F90 
     84      !!   9.0  !  05-10  (A. Beckmann)  correction for s-coordinates 
    8585      !!---------------------------------------------------------------------- 
    8686      !! * Modules used 
     
    9797      !! * Local declarations 
    9898      INTEGER  ::   ji, jj, jk              ! dummy loop indices 
    99       INTEGER  ::   ii0, ii1, ij0, ij1      ! temporary integer 
    100 #if defined key_partial_steps 
    101       INTEGER  ::   iku, ikv  ! temporary integers 
    102 #endif 
     99      INTEGER  ::   ii0, ii1, ij0, ij1,  &  ! temporary integer 
     100         &          iku, ikv                !    "          " 
    103101      REAL(wp) ::   & 
    104          zeps, zmg, zm05g, zcoef1, zcoef2,   &  ! temporary scalars 
    105          zau, zbu, zav, zbv,                 & 
    106          zai, zbi, zaj, zbj,                & 
    107          zcofu, zcofv, zcofw,                & 
    108          z1u, z1v, z1wu, z1wv,               & 
     102         zeps, zmg, zm05g,               &  ! temporary scalars 
     103         zcoef1, zcoef2, zcoef3,         &  ! 
     104         zau, zbu, zav, zbv,             & 
     105         zai, zbi, zaj, zbj,             & 
     106         zcofu, zcofv, zcofw,            & 
     107         z1u, z1v, z1wu, z1wv,           & 
    109108         zalpha 
    110109      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zww 
     
    138137      END DO 
    139138 
    140 #if defined key_partial_steps 
    141       ! partial steps correction at the bottom ocean level (zps_hde routine) 
    142 # if defined key_vectopt_loop   &&   ! defined key_autotasking 
    143       jj = 1 
    144       DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling) 
     139      IF( ln_zps ) THEN      ! partial steps correction at the bottom ocean level (zps_hde routine) 
     140# if defined key_vectopt_loop   &&   ! defined key_mpp_omp 
     141         jj = 1 
     142         DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling) 
    145143# else 
    146       DO jj = 1, jpjm1 
    147          DO ji = 1, jpim1 
    148 # endif 
    149             ! last ocean level 
    150             iku = MIN ( mbathy(ji,jj), mbathy(ji+1,jj) ) - 1 
    151             ikv = MIN ( mbathy(ji,jj), mbathy(ji,jj+1) ) - 1 
    152             zgru(ji,jj,iku) = gru(ji,jj)  
    153             zgrv(ji,jj,ikv) = grv(ji,jj)                
    154 # if ! defined key_vectopt_loop   ||   defined key_autotasking 
    155          END DO 
    156 # endif 
    157       END DO 
    158 #endif 
     144         DO jj = 1, jpjm1 
     145            DO ji = 1, jpim1 
     146# endif 
     147               ! last ocean level 
     148               iku = MIN ( mbathy(ji,jj), mbathy(ji+1,jj) ) - 1 
     149               ikv = MIN ( mbathy(ji,jj), mbathy(ji,jj+1) ) - 1 
     150               zgru(ji,jj,iku) = gru(ji,jj)  
     151               zgrv(ji,jj,ikv) = grv(ji,jj)                
     152# if ! defined key_vectopt_loop   ||   defined key_mpp_omp 
     153            END DO 
     154# endif 
     155         END DO 
     156      ENDIF 
    159157 
    160158      ! Slopes of isopycnal surfaces just below the mixed layer 
     
    203201               ! uslp and vslp output in zwz and zww, resp. 
    204202               zalpha = MAX( omlmask(ji,jj,jk), omlmask(ji+1,jj,jk) ) 
    205 #if defined key_s_coord  
    206203               zwz (ji,jj,jk) = ( zau / ( zbu - zeps ) * ( 1. - zalpha)   & 
    207                   &        + zalpha * uslpml(ji,jj)   & 
    208                   &        * ( fsdepu(ji,jj,jk) - .5*fse3u(ji,jj,1) )   & 
    209                   &        / MAX( hmlpt(ji,jj), hmlpt(ji+1,jj), 5. ) )   & 
    210                   &        * umask(ji,jj,jk) 
     204                  &           + zalpha * uslpml(ji,jj)                    & 
     205                  &                    * 0.5 * ( fsdept(ji+1,jj,jk)+fsdept(ji,jj,jk)-fse3u(ji,jj,1) )   & 
     206                  &                          / MAX( hmlpt(ji,jj), hmlpt(ji+1,jj), 5. ) ) * umask(ji,jj,jk) 
    211207               zalpha = MAX( omlmask(ji,jj,jk), omlmask(ji,jj+1,jk) ) 
    212208               zww (ji,jj,jk) = ( zav / ( zbv - zeps ) * ( 1. - zalpha)   & 
    213                   &        + zalpha * vslpml(ji,jj)   & 
    214                   &        * ( fsdepv(ji,jj,jk) - .5*fse3v(ji,jj,1) )   & 
    215                   &         / MAX( hmlpt(ji,jj), hmlpt(ji,jj+1), 5. ) )   & 
    216                   &        * vmask(ji,jj,jk) 
    217 #else 
    218                ! z-coord and partial steps slope computed in the same way 
    219                zwz (ji,jj,jk) = ( zau / ( zbu - zeps ) * ( 1. - zalpha)    & 
    220                   &        + zalpha * uslpml(ji,jj)    & 
    221                   &        * ( fsdept(ji,jj,jk) - .5*fse3u(ji,jj,1))    & 
    222                   &        / MAX (hmlpt(ji,jj),hmlpt(ji+1,jj),5.) )    & 
    223                   &        * umask (ji,jj,jk) 
    224                zalpha = MAX(omlmask(ji,jj,jk),omlmask(ji,jj+1,jk)) 
    225                zww (ji,jj,jk) = ( zav / ( zbv - zeps ) * ( 1. - zalpha)    & 
    226                   &        + zalpha * vslpml(ji,jj)    & 
    227                   &        * ( fsdept(ji,jj,jk) - .5*fse3v(ji,jj,1))    & 
    228                   &        / MAX(hmlpt(ji,jj),hmlpt(ji,jj+1),5.) )    & 
    229                   &        * vmask (ji,jj,jk) 
    230 #endif 
     209                  &           + zalpha * vslpml(ji,jj)                    & 
     210                  &                    * 0.5 * ( fsdept(ji,jj+1,jk)+fsdept(ji,jj,jk)-fse3v(ji,jj,1) )   & 
     211                  &                          / MAX( hmlpt(ji,jj), hmlpt(ji,jj+1), 5. ) ) * vmask(ji,jj,jk) 
    231212            END DO 
    232213         END DO 
     
    292273            END DO 
    293274         END DO 
    294  
    295  
    296          IF( lk_sco ) THEN 
    297             ! Add the slope of level surfaces 
    298             ! ----------------------------------- 
    299             ! 'key_s_coord' defined but not 'key_traldfiso' the computation is done 
    300             ! in inildf, ldfslp never called 
    301             ! 'key_s_coord' and 'key_traldfiso' defined, the slope of level surfaces 
    302             ! is added to the slope of isopycnal surfaces. 
    303             ! c a u t i o n : minus sign as fsdep has positive value  
    304           
    305             DO jj = 2, jpjm1 
    306                DO ji = fs_2, fs_jpim1   ! vector opt. 
    307                   uslp(ji,jj,jk) = uslp(ji,jj,jk) - 1. / e1u(ji,jj)   & 
    308                      &           * ( fsdept(ji+1,jj,jk) - fsdept(ji,jj,jk) ) 
    309                   vslp(ji,jj,jk) = vslp(ji,jj,jk) - 1. / e2v(ji,jj)   & 
    310                      &           * ( fsdept(ji,jj+1,jk) - fsdept(ji,jj,jk) ) 
    311                END DO 
    312             END DO 
    313          ENDIF 
    314275 
    315276 
     
    354315               zbj = MIN( zwy (ji,jj,jk), -100.*ABS(zaj), -7.e+3/fse3w(ji,jj,jk)*ABS(zaj) ) 
    355316               ! wslpi and wslpj output in zwz and zww, resp. 
    356                zalpha = MAX(omlmask(ji,jj,jk),omlmask(ji,jj,jk-1)) 
    357                zwz(ji,jj,jk) = ( zai / ( zbi - zeps) * ( 1. - zalpha )   & 
    358                   &            + zalpha * wslpiml(ji,jj)   & 
    359                   &            * fsdepw(ji,jj,jk) / MAX( hmlp(ji,jj),10. ) )   & 
    360                   &            * tmask (ji,jj,jk) 
    361                zww(ji,jj,jk) = ( zaj / ( zbj - zeps) * ( 1. - zalpha )   & 
    362                   &            + zalpha * wslpjml(ji,jj)   & 
    363                   &            * fsdepw(ji,jj,jk) / MAX( hmlp(ji,jj),10. ) )   & 
    364                   &            * tmask (ji,jj,jk) 
     317               zalpha = MAX( omlmask(ji,jj,jk), omlmask(ji,jj,jk-1) ) 
     318               zcoef3 = fsdepw(ji,jj,jk) / MAX( hmlp(ji,jj), 10. ) 
     319               zwz(ji,jj,jk) = (     zai / ( zbi - zeps)  * ( 1. - zalpha )   & 
     320                  &             + zcoef3 * wslpiml(ji,jj) *        zalpha   ) * tmask (ji,jj,jk) 
     321               zww(ji,jj,jk) = (     zaj / ( zbj - zeps)  * ( 1. - zalpha )   & 
     322                  &             + zcoef3 * wslpjml(ji,jj) *        zalpha   ) * tmask (ji,jj,jk) 
    365323            END DO 
    366324         END DO 
     
    424382         END DO 
    425383          
    426          IF( lk_sco ) THEN 
    427           
    428             ! Slope of level surfaces 
    429             ! ----------------------- 
    430             ! 'key_s_coord' defined but not 'key_traldfiso' the computation is done 
    431             ! in inildf, ldfslp never called 
    432             ! 'key_s_coord' and 'key_traldfiso' defined, the slope of level surfaces 
    433             ! is added to the slope of isopycnal surfaces. 
    434           
    435             DO jj = 2, jpjm1 
    436                DO ji = fs_2, fs_jpim1   ! vector opt. 
    437                   wslpi(ji,jj,jk) = wslpi(ji,jj,jk) - 1. / e1t(ji,jj)   & 
    438                      &                                   * ( fsdepuw(ji+1,jj,jk) - fsdepuw(ji,jj,jk) ) 
    439                   wslpj(ji,jj,jk) = wslpj(ji,jj,jk) - 1. / e2t(ji,jj)   & 
    440                      &                                   * ( fsdepvw(ji,jj+1,jk) - fsdepvw(ji,jj,jk) ) 
    441                END DO 
    442             END DO 
    443          ENDIF 
    444384          
    445385         ! III. Specific grid points 
     
    474414      ! III Lateral boundary conditions on all slopes (uslp , vslp,  
    475415      ! -------------------------------                wslpi, wslpj ) 
    476       CALL lbc_lnk( uslp , 'U', -1. ) 
    477       CALL lbc_lnk( vslp , 'V', -1. ) 
    478       CALL lbc_lnk( wslpi, 'W', -1. ) 
    479       CALL lbc_lnk( wslpj, 'W', -1. ) 
     416      CALL lbc_lnk( uslp , 'U', -1. )      ;      CALL lbc_lnk( vslp , 'V', -1. ) 
     417      CALL lbc_lnk( wslpi, 'W', -1. )      ;      CALL lbc_lnk( wslpj, 'W', -1. ) 
     418 
     419      IF(ln_ctl) THEN 
     420         CALL prt_ctl(tab3d_1=uslp , clinfo1=' slp  - u : ', tab3d_2=vslp,  clinfo2=' v : ', kdim=jpk) 
     421         CALL prt_ctl(tab3d_1=wslpi, clinfo1=' slp  - wi: ', tab3d_2=wslpj, clinfo2=' wj: ', kdim=jpk) 
     422      ENDIF 
    480423 
    481424   END SUBROUTINE ldf_slp 
     
    546489      ! mask for mixed layer 
    547490      DO jk = 1, jpk 
    548 # if defined key_vectopt_loop   &&   ! defined key_autotasking 
     491# if defined key_vectopt_loop   &&   ! defined key_mpp_omp 
    549492         jj = 1 
    550493         DO ji = 1, jpij   ! vector opt. (forced unrolling) 
     
    560503                  omlmask(ji,jj,jk) = 0.e0 
    561504               ENDIF 
    562 # if ! defined key_vectopt_loop   ||   defined key_autotasking 
     505# if ! defined key_vectopt_loop   ||   defined key_mpp_omp 
    563506            END DO 
    564507# endif 
     
    578521      zwy(:,jpj) = 0.e0 
    579522      zwy(jpi,:) = 0.e0 
    580 # if defined key_vectopt_loop   &&   ! defined key_autotasking 
     523# if defined key_vectopt_loop   &&   ! defined key_mpp_omp 
    581524      jj = 1 
    582525      DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling) 
     
    591534               &             * ( pn2(ji,jj,ik) + pn2(ji,jj,ik+1) )   & 
    592535               &             / MAX( tmask(ji,jj,ik) + tmask (ji,jj,ik+1), 1. ) 
    593 # if ! defined key_vectopt_loop   ||   defined key_autotasking 
     536# if ! defined key_vectopt_loop   ||   defined key_mpp_omp 
    594537         END DO 
    595538# endif 
     
    599542 
    600543      ! Slope at u points 
    601 # if defined key_vectopt_loop   &&   ! defined key_autotasking 
     544# if defined key_vectopt_loop   &&   ! defined key_mpp_omp 
    602545      jj = 1 
    603546      DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling) 
     
    616559            ! uslpml 
    617560            uslpml (ji,jj) = zau / ( zbu - zeps ) * umask (ji,jj,ik) 
    618 # if ! defined key_vectopt_loop   ||   defined key_autotasking 
     561# if ! defined key_vectopt_loop   ||   defined key_mpp_omp 
    619562         END DO 
    620563# endif 
     
    628571      zwy ( :, jpj) = 0.e0 
    629572      zwy ( jpi, :) = 0.e0 
    630 # if defined key_vectopt_loop   &&   ! defined key_autotasking 
     573# if defined key_vectopt_loop   &&   ! defined key_mpp_omp 
    631574      jj = 1 
    632575      DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling) 
     
    640583               &             * ( pn2(ji,jj,ik) + pn2(ji,jj,ik+1) )   & 
    641584               &             / MAX( tmask(ji,jj,ik) + tmask (ji,jj,ik+1), 1. ) 
    642 # if ! defined key_vectopt_loop   ||   defined key_autotasking 
     585# if ! defined key_vectopt_loop   ||   defined key_mpp_omp 
    643586         END DO 
    644587# endif 
     
    649592 
    650593      ! Slope at v points 
    651 # if defined key_vectopt_loop   &&   ! defined key_autotasking 
     594# if defined key_vectopt_loop   &&   ! defined key_mpp_omp 
    652595      jj = 1 
    653596      DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling) 
     
    666609            ! vslpml 
    667610            vslpml (ji,jj) = zav / ( zbv - zeps ) * vmask (ji,jj,ik) 
    668 # if ! defined key_vectopt_loop   ||   defined key_autotasking 
     611# if ! defined key_vectopt_loop   ||   defined key_mpp_omp 
    669612         END DO 
    670613# endif 
     
    680623      ! Local vertical density gradient evaluated from N^2 
    681624      ! zwy = d/dz(prd)= - mk ( prd ) / grav * pn2 -- at w point 
    682 # if defined key_vectopt_loop   &&   ! defined key_autotasking 
     625# if defined key_vectopt_loop   &&   ! defined key_mpp_omp 
    683626      jj = 1 
    684627      DO ji = 1, jpij   ! vector opt. (forced unrolling) 
     
    692635            zwy (ji,jj) = zm05g * pn2 (ji,jj,ik) *     & 
    693636               &             ( prd (ji,jj,ik) + prd (ji,jj,ikm1) + 2. ) 
    694 # if ! defined key_vectopt_loop   ||   defined key_autotasking 
     637# if ! defined key_vectopt_loop   ||   defined key_mpp_omp 
    695638         END DO 
    696639# endif 
     
    698641 
    699642      ! Slope at w point 
    700 # if defined key_vectopt_loop   &&   ! defined key_autotasking 
     643# if defined key_vectopt_loop   &&   ! defined key_mpp_omp 
    701644      jj = 1 
    702645      DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling) 
     
    728671            wslpiml (ji,jj) = zai / ( zbi - zeps) * tmask (ji,jj,ik) 
    729672            wslpjml (ji,jj) = zaj / ( zbj - zeps) * tmask (ji,jj,ik) 
    730 # if ! defined key_vectopt_loop   ||   defined key_autotasking 
     673# if ! defined key_vectopt_loop   ||   defined key_mpp_omp 
    731674         END DO 
    732675# endif 
     
    780723 
    781724      IF( ln_traldf_hor ) THEN 
     725         IF(lwp) THEN 
     726            WRITE(numout,*) '          Horizontal mixing in s-coordinate: slope = slope of s-surfaces' 
     727         ENDIF 
    782728 
    783729         ! geopotential diffusion in s-coordinates on tracers and/or momentum 
     
    790736            DO jj = 2, jpjm1 
    791737               DO ji = fs_2, fs_jpim1   ! vector opt. 
    792                   uslp (ji,jj,jk) = -1. / e1u(ji,jj) * umask(ji,jj,jk)   & 
    793                      &                               * ( fsdept(ji+1,jj,jk) - fsdept(ji,jj,jk) ) 
    794                   vslp (ji,jj,jk) = -1. / e2v(ji,jj) * vmask(ji,jj,jk)   & 
    795                      &                               * ( fsdept(ji,jj+1,jk) - fsdept(ji,jj,jk) ) 
    796                   wslpi(ji,jj,jk) = -1. / e1t(ji,jj) * tmask(ji,jj,jk)   & 
    797                      &                               * ( fsdepuw(ji+1,jj,jk) - fsdepuw(ji,jj,jk) ) 
    798                   wslpj(ji,jj,jk) = -1. / e2t(ji,jj) * tmask(ji,jj,jk)   & 
    799                      &                               * ( fsdepvw(ji,jj+1,jk) - fsdepvw(ji,jj,jk) ) 
     738                  uslp (ji,jj,jk) = -1./e1u(ji,jj) * ( fsdept(ji+1,jj,jk) - fsdept(ji ,jj ,jk) ) * umask(ji,jj,jk) 
     739                  vslp (ji,jj,jk) = -1./e2v(ji,jj) * ( fsdept(ji,jj+1,jk) - fsdept(ji ,jj ,jk) ) * vmask(ji,jj,jk) 
     740                  wslpi(ji,jj,jk) = -1./e1t(ji,jj) * ( fsdepw(ji+1,jj,jk) - fsdepw(ji-1,jj,jk) ) * tmask(ji,jj,jk) * 0.5 
     741                  wslpj(ji,jj,jk) = -1./e2t(ji,jj) * ( fsdepw(ji,jj+1,jk) - fsdepw(ji,jj-1,jk) ) * tmask(ji,jj,jk) * 0.5 
    800742               END DO 
    801743            END DO 
     
    803745 
    804746         ! Lateral boundary conditions on the slopes 
    805          CALL lbc_lnk( uslp , 'U', -1. ) 
    806          CALL lbc_lnk( vslp , 'V', -1. ) 
    807          CALL lbc_lnk( wslpi, 'W', -1. ) 
    808          CALL lbc_lnk( wslpj, 'W', -1. ) 
     747         CALL lbc_lnk( uslp , 'U', -1. )      ;      CALL lbc_lnk( vslp , 'V', -1. ) 
     748         CALL lbc_lnk( wslpi, 'W', -1. )      ;      CALL lbc_lnk( wslpj, 'W', -1. ) 
    809749      ENDIF 
    810750 
Note: See TracChangeset for help on using the changeset viewer.