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 14017 for NEMO/branches/2020/dev_r13648_ASINTER-04_laurent_bulk_ice/src/OCE/TRA/eosbn2.F90 – NEMO

Ignore:
Timestamp:
2020-12-02T16:32:24+01:00 (4 years ago)
Author:
laurent
Message:

Keep up with trunk revision 13999

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2020/dev_r13648_ASINTER-04_laurent_bulk_ice/src/OCE/TRA/eosbn2.F90

    r13497 r14017  
    3131   !!   bn2           : compute the Brunt-Vaisala frequency 
    3232   !!   eos_pt_from_ct: compute the potential temperature from the Conservative Temperature 
    33    !!   eos_rab       : generic interface of in situ thermal/haline expansion ratio  
     33   !!   eos_rab       : generic interface of in situ thermal/haline expansion ratio 
    3434   !!   eos_rab_3d    : compute in situ thermal/haline expansion ratio 
    3535   !!   eos_rab_2d    : compute in situ thermal/haline expansion ratio for 2d fields 
     
    3939   !!---------------------------------------------------------------------- 
    4040   USE dom_oce        ! ocean space and time domain 
     41   USE domutl, ONLY : is_tile 
    4142   USE phycst         ! physical constants 
    4243   USE stopar         ! Stochastic T/S fluctuations 
     
    4546   USE in_out_manager ! I/O manager 
    4647   USE lib_mpp        ! MPP library 
    47    USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)   
     48   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
    4849   USE prtctl         ! Print control 
    4950   USE lbclnk         ! ocean lateral boundary conditions 
     
    6263   END INTERFACE 
    6364   ! 
    64    INTERFACE eos_fzp  
     65   INTERFACE eos_fzp 
    6566      MODULE PROCEDURE eos_fzp_2d, eos_fzp_0d 
    6667   END INTERFACE 
     
    8889 
    8990   !                               !!!  simplified eos coefficients (default value: Vallis 2006) 
    90    REAL(wp) ::   rn_a0      = 1.6550e-1_wp     ! thermal expansion coeff.  
    91    REAL(wp) ::   rn_b0      = 7.6554e-1_wp     ! saline  expansion coeff.  
    92    REAL(wp) ::   rn_lambda1 = 5.9520e-2_wp     ! cabbeling coeff. in T^2         
    93    REAL(wp) ::   rn_lambda2 = 5.4914e-4_wp     ! cabbeling coeff. in S^2         
    94    REAL(wp) ::   rn_mu1     = 1.4970e-4_wp     ! thermobaric coeff. in T   
    95    REAL(wp) ::   rn_mu2     = 1.1090e-5_wp     ! thermobaric coeff. in S   
    96    REAL(wp) ::   rn_nu      = 2.4341e-3_wp     ! cabbeling coeff. in theta*salt   
    97     
     91   REAL(wp) ::   rn_a0      = 1.6550e-1_wp     ! thermal expansion coeff. 
     92   REAL(wp) ::   rn_b0      = 7.6554e-1_wp     ! saline  expansion coeff. 
     93   REAL(wp) ::   rn_lambda1 = 5.9520e-2_wp     ! cabbeling coeff. in T^2 
     94   REAL(wp) ::   rn_lambda2 = 5.4914e-4_wp     ! cabbeling coeff. in S^2 
     95   REAL(wp) ::   rn_mu1     = 1.4970e-4_wp     ! thermobaric coeff. in T 
     96   REAL(wp) ::   rn_mu2     = 1.1090e-5_wp     ! thermobaric coeff. in S 
     97   REAL(wp) ::   rn_nu      = 2.4341e-3_wp     ! cabbeling coeff. in theta*salt 
     98 
    9899   ! TEOS10/EOS80 parameters 
    99100   REAL(wp) ::   r1_S0, r1_T0, r1_Z0, rdeltaS 
    100     
     101 
    101102   ! EOS parameters 
    102103   REAL(wp) ::   EOS000 , EOS100 , EOS200 , EOS300 , EOS400 , EOS500 , EOS600 
     
    116117   REAL(wp) ::   EOS022 
    117118   REAL(wp) ::   EOS003 , EOS103 
    118    REAL(wp) ::   EOS013  
    119     
     119   REAL(wp) ::   EOS013 
     120 
    120121   ! ALPHA parameters 
    121122   REAL(wp) ::   ALP000 , ALP100 , ALP200 , ALP300 , ALP400 , ALP500 
     
    132133   REAL(wp) ::   ALP012 
    133134   REAL(wp) ::   ALP003 
    134     
     135 
    135136   ! BETA parameters 
    136137   REAL(wp) ::   BET000 , BET100 , BET200 , BET300 , BET400 , BET500 
     
    159160   REAL(wp) ::   PEN002 , PEN102 
    160161   REAL(wp) ::   PEN012 
    161     
     162 
    162163   ! ALPHA_PEN parameters 
    163164   REAL(wp) ::   APE000 , APE100 , APE200 , APE300 
     
    189190 
    190191   SUBROUTINE eos_insitu( pts, prd, pdep ) 
     192      !! 
     193      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celsius] 
     194      !                                                      ! 2 : salinity               [psu] 
     195      REAL(wp), DIMENSION(:,:,:)  , INTENT(  out) ::   prd   ! in situ density            [-] 
     196      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pdep  ! depth                      [m] 
     197      !! 
     198      CALL eos_insitu_t( pts, is_tile(pts), prd, is_tile(prd), pdep, is_tile(pdep) ) 
     199   END SUBROUTINE eos_insitu 
     200 
     201   SUBROUTINE eos_insitu_t( pts, ktts, prd, ktrd, pdep, ktdep ) 
    191202      !!---------------------------------------------------------------------- 
    192203      !!                   ***  ROUTINE eos_insitu  *** 
     
    222233      !!                TEOS-10 Manual, 2010 
    223234      !!---------------------------------------------------------------------- 
    224       REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celsius] 
    225       !                                                               ! 2 : salinity               [psu] 
    226       REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(  out) ::   prd   ! in situ density            [-] 
    227       REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pdep  ! depth                      [m] 
     235      INTEGER                                 , INTENT(in   ) ::   ktts, ktrd, ktdep 
     236      REAL(wp), DIMENSION(A2D_T(ktts) ,JPK,JPTS), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celsius] 
     237      !                                                                  ! 2 : salinity               [psu] 
     238      REAL(wp), DIMENSION(A2D_T(ktrd) ,JPK     ), INTENT(  out) ::   prd   ! in situ density            [-] 
     239      REAL(wp), DIMENSION(A2D_T(ktdep),JPK     ), INTENT(in   ) ::   pdep  ! depth                      [m] 
    228240      ! 
    229241      INTEGER  ::   ji, jj, jk                ! dummy loop indices 
     
    238250      CASE( np_teos10, np_eos80 )                !==  polynomial TEOS-10 / EOS-80 ==! 
    239251         ! 
    240          DO_3D( 1, 1, 1, 1, 1, jpkm1 ) 
     252         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 ) 
    241253            ! 
    242254            zh  = pdep(ji,jj,jk) * r1_Z0                                  ! depth 
     
    274286      CASE( np_seos )                !==  simplified EOS  ==! 
    275287         ! 
    276          DO_3D( 1, 1, 1, 1, 1, jpkm1 ) 
     288         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 ) 
    277289            zt  = pts  (ji,jj,jk,jp_tem) - 10._wp 
    278290            zs  = pts  (ji,jj,jk,jp_sal) - 35._wp 
     
    283295               &  + rn_b0 * ( 1._wp - 0.5_wp*rn_lambda2*zs - rn_mu2*zh ) * zs   & 
    284296               &  - rn_nu * zt * zs 
    285                !                                  
     297               ! 
    286298            prd(ji,jj,jk) = zn * r1_rho0 * ztm                ! density anomaly (masked) 
    287299         END_3D 
     
    293305      IF( ln_timing )   CALL timing_stop('eos-insitu') 
    294306      ! 
    295    END SUBROUTINE eos_insitu 
     307   END SUBROUTINE eos_insitu_t 
    296308 
    297309 
    298310   SUBROUTINE eos_insitu_pot( pts, prd, prhop, pdep ) 
     311      !! 
     312      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pts    ! 1 : potential temperature  [Celsius] 
     313      !                                                       ! 2 : salinity               [psu] 
     314      REAL(wp), DIMENSION(:,:,:)  , INTENT(  out) ::   prd    ! in situ density            [-] 
     315      REAL(wp), DIMENSION(:,:,:)  , INTENT(  out) ::   prhop  ! potential density (surface referenced) 
     316      REAL(wp), DIMENSION(:,:,:)  , INTENT(in   ) ::   pdep   ! depth                      [m] 
     317      !! 
     318      CALL eos_insitu_pot_t( pts, is_tile(pts), prd, is_tile(prd), prhop, is_tile(prhop), pdep, is_tile(pdep) ) 
     319   END SUBROUTINE eos_insitu_pot 
     320 
     321 
     322   SUBROUTINE eos_insitu_pot_t( pts, ktts, prd, ktrd, prhop, ktrhop, pdep, ktdep ) 
    299323      !!---------------------------------------------------------------------- 
    300324      !!                  ***  ROUTINE eos_insitu_pot  *** 
     
    309333      !! 
    310334      !!---------------------------------------------------------------------- 
    311       REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in   ) ::   pts    ! 1 : potential temperature  [Celsius] 
    312       !                                                                ! 2 : salinity               [psu] 
    313       REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(  out) ::   prd    ! in situ density            [-] 
    314       REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(  out) ::   prhop  ! potential density (surface referenced) 
    315       REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pdep   ! depth                      [m] 
     335      INTEGER                                  , INTENT(in   ) ::   ktts, ktrd, ktrhop, ktdep 
     336      REAL(wp), DIMENSION(A2D_T(ktts)  ,JPK,JPTS), INTENT(in   ) ::   pts    ! 1 : potential temperature  [Celsius] 
     337      !                                                                    ! 2 : salinity               [psu] 
     338      REAL(wp), DIMENSION(A2D_T(ktrd)  ,JPK     ), INTENT(  out) ::   prd    ! in situ density            [-] 
     339      REAL(wp), DIMENSION(A2D_T(ktrhop),JPK     ), INTENT(  out) ::   prhop  ! potential density (surface referenced) 
     340      REAL(wp), DIMENSION(A2D_T(ktdep) ,JPK     ), INTENT(in   ) ::   pdep   ! depth                      [m] 
    316341      ! 
    317342      INTEGER  ::   ji, jj, jk, jsmp             ! dummy loop indices 
     
    338363            END DO 
    339364            ! 
    340             DO_3D( 1, 1, 1, 1, 1, jpkm1 ) 
     365            DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 ) 
    341366               ! 
    342367               ! compute density (2*nn_sto_eos) times: 
     
    388413         ! Non-stochastic equation of state 
    389414         ELSE 
    390             DO_3D( 1, 1, 1, 1, 1, jpkm1 ) 
     415            DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 ) 
    391416               ! 
    392417               zh  = pdep(ji,jj,jk) * r1_Z0                                  ! depth 
     
    423448            END_3D 
    424449         ENDIF 
    425           
     450 
    426451      CASE( np_seos )                !==  simplified EOS  ==! 
    427452         ! 
    428          DO_3D( 1, 1, 1, 1, 1, jpkm1 ) 
     453         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 ) 
    429454            zt  = pts  (ji,jj,jk,jp_tem) - 10._wp 
    430455            zs  = pts  (ji,jj,jk,jp_sal) - 35._wp 
     
    444469      END SELECT 
    445470      ! 
    446       IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab3d_1=prd, clinfo1=' eos-pot: ', tab3d_2=prhop, clinfo2=' pot : ', kdim=jpk ) 
     471      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab3d_1=prd, clinfo1=' eos-pot: ', & 
     472         &                                  tab3d_2=prhop, clinfo2=' pot : ', kdim=jpk ) 
    447473      ! 
    448474      IF( ln_timing )   CALL timing_stop('eos-pot') 
    449475      ! 
    450    END SUBROUTINE eos_insitu_pot 
     476   END SUBROUTINE eos_insitu_pot_t 
    451477 
    452478 
    453479   SUBROUTINE eos_insitu_2d( pts, pdep, prd ) 
     480      !! 
     481      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celsius] 
     482      !                                                    ! 2 : salinity               [psu] 
     483      REAL(wp), DIMENSION(:,:)  , INTENT(in   ) ::   pdep  ! depth                      [m] 
     484      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   prd   ! in situ density 
     485      !! 
     486      CALL eos_insitu_2d_t( pts, is_tile(pts), pdep, is_tile(pdep), prd, is_tile(prd) ) 
     487   END SUBROUTINE eos_insitu_2d 
     488 
     489 
     490   SUBROUTINE eos_insitu_2d_t( pts, ktts, pdep, ktdep, prd, ktrd ) 
    454491      !!---------------------------------------------------------------------- 
    455492      !!                  ***  ROUTINE eos_insitu_2d  *** 
     
    462499      !! 
    463500      !!---------------------------------------------------------------------- 
    464       REAL(wp), DIMENSION(jpi,jpj,jpts), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celsius] 
    465       !                                                           ! 2 : salinity               [psu] 
    466       REAL(wp), DIMENSION(jpi,jpj)     , INTENT(in   ) ::   pdep  ! depth                      [m] 
    467       REAL(wp), DIMENSION(jpi,jpj)     , INTENT(  out) ::   prd   ! in situ density 
     501      INTEGER                            , INTENT(in   ) ::   ktts, ktdep, ktrd 
     502      REAL(wp), DIMENSION(A2D_T(ktts),JPTS), INTENT(in   ) ::   pts   ! 1 : potential temperature  [Celsius] 
     503      !                                                             ! 2 : salinity               [psu] 
     504      REAL(wp), DIMENSION(A2D_T(ktdep)    ), INTENT(in   ) ::   pdep  ! depth                      [m] 
     505      REAL(wp), DIMENSION(A2D_T(ktrd)     ), INTENT(  out) ::   prd   ! in situ density 
    468506      ! 
    469507      INTEGER  ::   ji, jj, jk                ! dummy loop indices 
     
    480518      CASE( np_teos10, np_eos80 )                !==  polynomial TEOS-10 / EOS-80 ==! 
    481519         ! 
    482          DO_2D( 1, 1, 1, 1 ) 
     520         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) 
    483521            ! 
    484522            zh  = pdep(ji,jj) * r1_Z0                                  ! depth 
     
    515553      CASE( np_seos )                !==  simplified EOS  ==! 
    516554         ! 
    517          DO_2D( 1, 1, 1, 1 ) 
     555         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) 
    518556            ! 
    519557            zt    = pts  (ji,jj,jp_tem)  - 10._wp 
     
    535573      IF( ln_timing )   CALL timing_stop('eos2d') 
    536574      ! 
    537    END SUBROUTINE eos_insitu_2d 
     575   END SUBROUTINE eos_insitu_2d_t 
    538576 
    539577 
    540578   SUBROUTINE rab_3d( pts, pab, Kmm ) 
     579      !! 
     580      INTEGER                     , INTENT(in   ) ::   Kmm   ! time level index 
     581      REAL(wp), DIMENSION(:,:,:,:), INTENT(in   ) ::   pts   ! pot. temperature & salinity 
     582      REAL(wp), DIMENSION(:,:,:,:), INTENT(  out) ::   pab   ! thermal/haline expansion ratio 
     583      !! 
     584      CALL rab_3d_t( pts, is_tile(pts), pab, is_tile(pab), Kmm ) 
     585   END SUBROUTINE rab_3d 
     586 
     587 
     588   SUBROUTINE rab_3d_t( pts, ktts, pab, ktab, Kmm ) 
    541589      !!---------------------------------------------------------------------- 
    542590      !!                 ***  ROUTINE rab_3d  *** 
     
    548596      !! ** Action  : - pab     : thermal/haline expansion ratio at T-points 
    549597      !!---------------------------------------------------------------------- 
    550       INTEGER                              , INTENT(in   ) ::   Kmm   ! time level index 
    551       REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in   ) ::   pts   ! pot. temperature & salinity 
    552       REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(  out) ::   pab   ! thermal/haline expansion ratio 
     598      INTEGER                                , INTENT(in   ) ::   Kmm   ! time level index 
     599      INTEGER                                , INTENT(in   ) ::   ktts, ktab 
     600      REAL(wp), DIMENSION(A2D_T(ktts),JPK,JPTS), INTENT(in   ) ::   pts   ! pot. temperature & salinity 
     601      REAL(wp), DIMENSION(A2D_T(ktab),JPK,JPTS), INTENT(  out) ::   pab   ! thermal/haline expansion ratio 
    553602      ! 
    554603      INTEGER  ::   ji, jj, jk                ! dummy loop indices 
     
    563612      CASE( np_teos10, np_eos80 )                !==  polynomial TEOS-10 / EOS-80 ==! 
    564613         ! 
    565          DO_3D( 1, 1, 1, 1, 1, jpkm1 ) 
     614         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 ) 
    566615            ! 
    567616            zh  = gdept(ji,jj,jk,Kmm) * r1_Z0                                ! depth 
     
    616665      CASE( np_seos )                  !==  simplified EOS  ==! 
    617666         ! 
    618          DO_3D( 1, 1, 1, 1, 1, jpkm1 ) 
     667         DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1 ) 
    619668            zt  = pts (ji,jj,jk,jp_tem) - 10._wp   ! pot. temperature anomaly (t-T0) 
    620669            zs  = pts (ji,jj,jk,jp_sal) - 35._wp   ! abs. salinity anomaly (s-S0) 
     
    641690      IF( ln_timing )   CALL timing_stop('rab_3d') 
    642691      ! 
    643    END SUBROUTINE rab_3d 
     692   END SUBROUTINE rab_3d_t 
    644693 
    645694 
    646695   SUBROUTINE rab_2d( pts, pdep, pab, Kmm ) 
     696      !! 
     697      INTEGER                   , INTENT(in   ) ::   Kmm   ! time level index 
     698      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pts    ! pot. temperature & salinity 
     699      REAL(wp), DIMENSION(:,:)  , INTENT(in   ) ::   pdep   ! depth                  [m] 
     700      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   pab    ! thermal/haline expansion ratio 
     701      !! 
     702      CALL rab_2d_t(pts, is_tile(pts), pdep, is_tile(pdep), pab, is_tile(pab), Kmm) 
     703   END SUBROUTINE rab_2d 
     704 
     705 
     706   SUBROUTINE rab_2d_t( pts, ktts, pdep, ktdep, pab, ktab, Kmm ) 
    647707      !!---------------------------------------------------------------------- 
    648708      !!                 ***  ROUTINE rab_2d  *** 
     
    652712      !! ** Action  : - pab     : thermal/haline expansion ratio at T-points 
    653713      !!---------------------------------------------------------------------- 
    654       INTEGER                              , INTENT(in   ) ::   Kmm   ! time level index 
    655       REAL(wp), DIMENSION(jpi,jpj,jpts)    , INTENT(in   ) ::   pts    ! pot. temperature & salinity 
    656       REAL(wp), DIMENSION(jpi,jpj)         , INTENT(in   ) ::   pdep   ! depth                  [m] 
    657       REAL(wp), DIMENSION(jpi,jpj,jpts)    , INTENT(  out) ::   pab    ! thermal/haline expansion ratio 
     714      INTEGER                            , INTENT(in   ) ::   Kmm   ! time level index 
     715      INTEGER                            , INTENT(in   ) ::   ktts, ktdep, ktab 
     716      REAL(wp), DIMENSION(A2D_T(ktts),JPTS), INTENT(in   ) ::   pts    ! pot. temperature & salinity 
     717      REAL(wp), DIMENSION(A2D_T(ktdep)    ), INTENT(in   ) ::   pdep   ! depth                  [m] 
     718      REAL(wp), DIMENSION(A2D_T(ktab),JPTS), INTENT(  out) ::   pab    ! thermal/haline expansion ratio 
    658719      ! 
    659720      INTEGER  ::   ji, jj, jk                ! dummy loop indices 
     
    670731      CASE( np_teos10, np_eos80 )                !==  polynomial TEOS-10 / EOS-80 ==! 
    671732         ! 
    672          DO_2D( 1, 1, 1, 1 ) 
     733         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) 
    673734            ! 
    674735            zh  = pdep(ji,jj) * r1_Z0                                  ! depth 
     
    723784      CASE( np_seos )                  !==  simplified EOS  ==! 
    724785         ! 
    725          DO_2D( 1, 1, 1, 1 ) 
     786         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) 
    726787            ! 
    727788            zt    = pts  (ji,jj,jp_tem) - 10._wp   ! pot. temperature anomaly (t-T0) 
     
    748809      IF( ln_timing )   CALL timing_stop('rab_2d') 
    749810      ! 
    750    END SUBROUTINE rab_2d 
     811   END SUBROUTINE rab_2d_t 
    751812 
    752813 
     
    849910 
    850911   SUBROUTINE bn2( pts, pab, pn2, Kmm ) 
     912      !! 
     913      INTEGER                              , INTENT(in   ) ::  Kmm   ! time level index 
     914      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in   ) ::  pts   ! pot. temperature and salinity   [Celsius,psu] 
     915      REAL(wp), DIMENSION(:,:,:,:)         , INTENT(in   ) ::  pab   ! thermal/haline expansion coef.  [Celsius-1,psu-1] 
     916      REAL(wp), DIMENSION(:,:,:)           , INTENT(  out) ::  pn2   ! Brunt-Vaisala frequency squared [1/s^2] 
     917      !! 
     918      CALL bn2_t( pts, pab, is_tile(pab), pn2, is_tile(pn2), Kmm ) 
     919   END SUBROUTINE bn2 
     920 
     921 
     922   SUBROUTINE bn2_t( pts, pab, ktab, pn2, ktn2, Kmm ) 
    851923      !!---------------------------------------------------------------------- 
    852924      !!                  ***  ROUTINE bn2  *** 
    853925      !! 
    854       !! ** Purpose :   Compute the local Brunt-Vaisala frequency at the  
     926      !! ** Purpose :   Compute the local Brunt-Vaisala frequency at the 
    855927      !!                time-step of the input arguments 
    856928      !! 
     
    859931      !!      N.B. N^2 is set one for all to zero at jk=1 in istate module. 
    860932      !! 
    861       !! ** Action  :   pn2 : square of the brunt-vaisala frequency at w-point  
    862       !! 
    863       !!---------------------------------------------------------------------- 
    864       INTEGER                              , INTENT(in   ) ::   Kmm   ! time level index 
    865       REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in   ) ::  pts   ! pot. temperature and salinity   [Celsius,psu] 
    866       REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in   ) ::  pab   ! thermal/haline expansion coef.  [Celsius-1,psu-1] 
    867       REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(  out) ::  pn2   ! Brunt-Vaisala frequency squared [1/s^2] 
     933      !! ** Action  :   pn2 : square of the brunt-vaisala frequency at w-point 
     934      !! 
     935      !!---------------------------------------------------------------------- 
     936      INTEGER                                , INTENT(in   ) ::  Kmm   ! time level index 
     937      INTEGER                                , INTENT(in   ) ::  ktab, ktn2 
     938      REAL(wp), DIMENSION(jpi,jpj,  jpk,jpts), INTENT(in   ) ::  pts   ! pot. temperature and salinity   [Celsius,psu] 
     939      REAL(wp), DIMENSION(A2D_T(ktab),JPK,JPTS), INTENT(in   ) ::  pab   ! thermal/haline expansion coef.  [Celsius-1,psu-1] 
     940      REAL(wp), DIMENSION(A2D_T(ktn2),JPK     ), INTENT(  out) ::  pn2   ! Brunt-Vaisala frequency squared [1/s^2] 
    868941      ! 
    869942      INTEGER  ::   ji, jj, jk      ! dummy loop indices 
     
    873946      IF( ln_timing )   CALL timing_start('bn2') 
    874947      ! 
    875       DO_3D( 1, 1, 1, 1, 2, jpkm1 )      ! interior points only (2=< jk =< jpkm1 ); surface and bottom value set to zero one for all in istate.F90 
     948      DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 2, jpkm1 )      ! interior points only (2=< jk =< jpkm1 ); surface and bottom value set to zero one for all in istate.F90 
    876949         zrw =   ( gdepw(ji,jj,jk  ,Kmm) - gdept(ji,jj,jk,Kmm) )   & 
    877             &  / ( gdept(ji,jj,jk-1,Kmm) - gdept(ji,jj,jk,Kmm) )  
    878             ! 
    879          zaw = pab(ji,jj,jk,jp_tem) * (1. - zrw) + pab(ji,jj,jk-1,jp_tem) * zrw  
     950            &  / ( gdept(ji,jj,jk-1,Kmm) - gdept(ji,jj,jk,Kmm) ) 
     951            ! 
     952         zaw = pab(ji,jj,jk,jp_tem) * (1. - zrw) + pab(ji,jj,jk-1,jp_tem) * zrw 
    880953         zbw = pab(ji,jj,jk,jp_sal) * (1. - zrw) + pab(ji,jj,jk-1,jp_sal) * zrw 
    881954         ! 
     
    889962      IF( ln_timing )   CALL timing_stop('bn2') 
    890963      ! 
    891    END SUBROUTINE bn2 
     964   END SUBROUTINE bn2_t 
    892965 
    893966 
     
    9491022 
    9501023 
    951    SUBROUTINE  eos_fzp_2d( psal, ptf, pdep ) 
     1024   SUBROUTINE eos_fzp_2d( psal, ptf, pdep ) 
     1025      !! 
     1026      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   )           ::   psal   ! salinity   [psu] 
     1027      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ), OPTIONAL ::   pdep   ! depth      [m] 
     1028      REAL(wp), DIMENSION(:,:)    , INTENT(out  )           ::   ptf    ! freezing temperature [Celsius] 
     1029      !! 
     1030      CALL eos_fzp_2d_t( psal, ptf, is_tile(ptf), pdep ) 
     1031   END SUBROUTINE eos_fzp_2d 
     1032 
     1033 
     1034   SUBROUTINE  eos_fzp_2d_t( psal, ptf, kttf, pdep ) 
    9521035      !!---------------------------------------------------------------------- 
    9531036      !!                 ***  ROUTINE eos_fzp  *** 
     
    9611044      !! Reference  :   UNESCO tech. papers in the marine science no. 28. 1978 
    9621045      !!---------------------------------------------------------------------- 
    963       REAL(wp), DIMENSION(jpi,jpj), INTENT(in   )           ::   psal   ! salinity   [psu] 
    964       REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ), OPTIONAL ::   pdep   ! depth      [m] 
    965       REAL(wp), DIMENSION(jpi,jpj), INTENT(out  )           ::   ptf    ! freezing temperature [Celsius] 
     1046      INTEGER                       , INTENT(in   )           ::   kttf 
     1047      REAL(wp), DIMENSION(jpi,jpj)  , INTENT(in   )           ::   psal   ! salinity   [psu] 
     1048      REAL(wp), DIMENSION(jpi,jpj)  , INTENT(in   ), OPTIONAL ::   pdep   ! depth      [m] 
     1049      REAL(wp), DIMENSION(A2D_T(kttf)), INTENT(out  )           ::   ptf    ! freezing temperature [Celsius] 
    9661050      ! 
    9671051      INTEGER  ::   ji, jj          ! dummy loop indices 
     
    9941078         CALL ctl_stop( 'eos_fzp_2d:', ctmp1 ) 
    9951079         ! 
    996       END SELECT       
    997       ! 
    998   END SUBROUTINE eos_fzp_2d 
     1080      END SELECT 
     1081      ! 
     1082  END SUBROUTINE eos_fzp_2d_t 
    9991083 
    10001084 
     
    10511135      !! ** Purpose :   Calculates nonlinear anomalies of alpha_PE, beta_PE and PE at T-points 
    10521136      !! 
    1053       !! ** Method  :   PE is defined analytically as the vertical  
     1137      !! ** Method  :   PE is defined analytically as the vertical 
    10541138      !!                   primitive of EOS times -g integrated between 0 and z>0. 
    10551139      !!                pen is the nonlinear bsq-PE anomaly: pen = ( PE - rho0 gz ) / rho0 gz - rd 
    1056       !!                                                      = 1/z * /int_0^z rd dz - rd  
     1140      !!                                                      = 1/z * /int_0^z rd dz - rd 
    10571141      !!                                where rd is the density anomaly (see eos_rhd function) 
    10581142      !!                ab_pe are partial derivatives of PE anomaly with respect to T and S: 
     
    11181202               ! 
    11191203            zn  = ( zn2 * zh + zn1 ) * zh + zn0 
    1120             !                               
     1204            ! 
    11211205            pab_pe(ji,jj,jk,jp_tem) = zn * zh * r1_rho0 * ztm 
    11221206            ! 
     
    11331217               ! 
    11341218            zn  = ( zn2 * zh + zn1 ) * zh + zn0 
    1135             !                               
     1219            ! 
    11361220            pab_pe(ji,jj,jk,jp_sal) = zn / zs * zh * r1_rho0 * ztm 
    11371221            ! 
     
    12131297         IF(lwp) WRITE(numout,*) '   ==>>>   use of TEOS-10 equation of state (cons. temp. and abs. salinity)' 
    12141298         ! 
    1215          l_useCT = .TRUE.                          ! model temperature is Conservative temperature  
     1299         l_useCT = .TRUE.                          ! model temperature is Conservative temperature 
    12161300         ! 
    12171301         rdeltaS = 32._wp 
     
    15941678 
    15951679         r1_S0  = 0.875_wp/35.16504_wp   ! Used to convert CT in potential temperature when using bulk formulae (eos_pt_from_ct) 
    1596           
     1680 
    15971681         IF(lwp) THEN 
    15981682            WRITE(numout,*) 
     
    16181702      END SELECT 
    16191703      ! 
    1620       rho0_rcp    = rho0 * rcp  
     1704      rho0_rcp    = rho0 * rcp 
    16211705      r1_rho0     = 1._wp / rho0 
    16221706      r1_rcp      = 1._wp / rcp 
    1623       r1_rho0_rcp = 1._wp / rho0_rcp  
     1707      r1_rho0_rcp = 1._wp / rho0_rcp 
    16241708      ! 
    16251709      IF(lwp) THEN 
Note: See TracChangeset for help on using the changeset viewer.