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 12377 for NEMO/trunk/src/OCE/FLO – NEMO

Ignore:
Timestamp:
2020-02-12T15:39:06+01:00 (4 years ago)
Author:
acc
Message:

The big one. Merging all 2019 developments from the option 1 branch back onto the trunk.

This changeset reproduces 2019/dev_r11943_MERGE_2019 on the trunk using a 2-URL merge
onto a working copy of the trunk. I.e.:

svn merge --ignore-ancestry \

svn+ssh://acc@forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/NEMO/trunk \
svn+ssh://acc@forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/NEMO/branches/2019/dev_r11943_MERGE_2019 ./

The --ignore-ancestry flag avoids problems that may otherwise arise from the fact that
the merge history been trunk and branch may have been applied in a different order but
care has been taken before this step to ensure that all applicable fixes and updates
are present in the merge branch.

The trunk state just before this step has been branched to releases/release-4.0-HEAD
and that branch has been immediately tagged as releases/release-4.0.2. Any fixes
or additions in response to tickets on 4.0, 4.0.1 or 4.0.2 should be done on
releases/release-4.0-HEAD. From now on future 'point' releases (e.g. 4.0.2) will
remain unchanged with periodic releases as needs demand. Note release-4.0-HEAD is a
transitional naming convention. Future full releases, say 4.2, will have a release-4.2
branch which fulfills this role and the first point release (e.g. 4.2.0) will be made
immediately following the release branch creation.

2020 developments can be started from any trunk revision later than this one.

Location:
NEMO/trunk
Files:
6 edited

Legend:

Unmodified
Added
Removed
  • NEMO/trunk

    • Property svn:externals
      •  

        old new  
        33^/utils/build/mk@HEAD         mk 
        44^/utils/tools@HEAD            tools 
        5 ^/vendors/AGRIF/dev@HEAD      ext/AGRIF 
         5^/vendors/AGRIF/dev_r11615_ENHANCE-04_namelists_as_internalfiles_agrif@HEAD      ext/AGRIF 
        66^/vendors/FCM@HEAD            ext/FCM 
        77^/vendors/IOIPSL@HEAD         ext/IOIPSL 
  • NEMO/trunk/src/OCE/FLO/flo4rk.F90

    r11536 r12377  
    3333CONTAINS 
    3434 
    35    SUBROUTINE flo_4rk( kt ) 
     35   SUBROUTINE flo_4rk( kt, Kbb, Kmm ) 
    3636      !!---------------------------------------------------------------------- 
    3737      !!                  ***  ROUTINE flo_4rk  *** 
     
    4545      !!       floats and the grid defined on the domain. 
    4646      !!---------------------------------------------------------------------- 
    47       INTEGER, INTENT(in) ::   kt   ! ocean time-step index 
     47      INTEGER, INTENT(in) ::   kt         ! ocean time-step index 
     48      INTEGER, INTENT(in) ::   Kbb, Kmm   ! ocean time level indices 
    4849      !! 
    4950      INTEGER ::  jfl, jind           ! dummy loop indices 
     
    125126       
    126127         ! for each step we compute the compute the velocity with Lagrange interpolation 
    127          CALL flo_interp( zgifl, zgjfl, zgkfl, zufl, zvfl, zwfl, jind ) 
     128         CALL flo_interp( Kbb, Kmm, zgifl, zgjfl, zgkfl, zufl, zvfl, zwfl, jind ) 
    128129          
    129130         ! computation of Runge-Kutta factor 
     
    153154 
    154155 
    155    SUBROUTINE flo_interp( pxt , pyt , pzt ,      & 
     156   SUBROUTINE flo_interp( Kbb, Kmm,              & 
     157      &                   pxt , pyt , pzt ,      & 
    156158      &                   pufl, pvfl, pwfl, ki ) 
    157159      !!---------------------------------------------------------------------- 
     
    165167      !!      integrated with RK method. 
    166168      !!---------------------------------------------------------------------- 
     169      INTEGER                    , INTENT(in   ) ::   Kbb, Kmm           ! ocean time level indices 
    167170      REAL(wp) , DIMENSION(jpnfl), INTENT(in   ) ::   pxt , pyt , pzt    ! position of the float 
    168171      REAL(wp) , DIMENSION(jpnfl), INTENT(  out) ::   pufl, pvfl, pwfl   ! velocity at this position 
     
    246249               DO jind3 = 1, 4 
    247250                  ztufl(jfl,jind1,jind2,jind3) =   & 
    248                      &   (  tcoef1(ki) * ub(iidu(jfl,jind1),ijdu(jfl,jind2),ikdu(jfl,jind3)) +   & 
    249                      &      tcoef2(ki) * un(iidu(jfl,jind1),ijdu(jfl,jind2),ikdu(jfl,jind3)) )   & 
     251                     &   (  tcoef1(ki) * uu(iidu(jfl,jind1),ijdu(jfl,jind2),ikdu(jfl,jind3),Kbb) +   & 
     252                     &      tcoef2(ki) * uu(iidu(jfl,jind1),ijdu(jfl,jind2),ikdu(jfl,jind3),Kmm) )   & 
    250253                     &      / e1u(iidu(jfl,jind1),ijdu(jfl,jind2))  
    251254               END DO 
     
    330333               DO jind3 = 1 ,4 
    331334                  ztvfl(jfl,jind1,jind2,jind3)=   & 
    332                      &   ( tcoef1(ki) * vb(iidv(jfl,jind1),ijdv(jfl,jind2),ikdv(jfl,jind3))  +   & 
    333                      &     tcoef2(ki) * vn(iidv(jfl,jind1),ijdv(jfl,jind2),ikdv(jfl,jind3)) )    &  
     335                     &   ( tcoef1(ki) * vv(iidv(jfl,jind1),ijdv(jfl,jind2),ikdv(jfl,jind3),Kbb)  +   & 
     336                     &     tcoef2(ki) * vv(iidv(jfl,jind1),ijdv(jfl,jind2),ikdv(jfl,jind3),Kmm) )    &  
    334337                     &     / e2v(iidv(jfl,jind1),ijdv(jfl,jind2)) 
    335338               END DO 
     
    422425                  ztwfl(jfl,jind1,jind2,jind3)=   & 
    423426                     &   ( tcoef1(ki) * wb(iidw(jfl,jind1),ijdw(jfl,jind2),ikdw(jfl,jind3))+   & 
    424                      &     tcoef2(ki) * wn(iidw(jfl,jind1),ijdw(jfl,jind2),ikdw(jfl,jind3)) )  & 
    425                      &   / e3w_n(iidw(jfl,jind1),ijdw(jfl,jind2),ikdw(jfl,jind3)) 
     427                     &     tcoef2(ki) * ww(iidw(jfl,jind1),ijdw(jfl,jind2),ikdw(jfl,jind3)) )  & 
     428                     &   / e3w(iidw(jfl,jind1),ijdw(jfl,jind2),ikdw(jfl,jind3),Kmm) 
    426429               END DO 
    427430            END DO 
  • NEMO/trunk/src/OCE/FLO/floats.F90

    r11536 r12377  
    3737CONTAINS 
    3838 
    39    SUBROUTINE flo_stp( kt ) 
     39   SUBROUTINE flo_stp( kt, Kbb, Kmm ) 
    4040      !!---------------------------------------------------------------------- 
    4141      !!                   ***  ROUTINE flo_stp  *** 
     
    4848      !!        if ln_flork4 =T 
    4949      !!---------------------------------------------------------------------- 
    50       INTEGER, INTENT( in  ) ::   kt   ! ocean time step 
     50      INTEGER, INTENT( in  ) ::   kt        ! ocean time step 
     51      INTEGER, INTENT( in  ) ::   Kbb, Kmm  ! ocean time level indices  
    5152      !!---------------------------------------------------------------------- 
    5253      ! 
    5354      IF( ln_timing )   CALL timing_start('flo_stp') 
    5455      ! 
    55       IF( ln_flork4 ) THEN   ;   CALL flo_4rk( kt )        ! Trajectories using a 4th order Runge Kutta scheme 
    56       ELSE                   ;   CALL flo_blk( kt )        ! Trajectories using Blanke' algorithme 
     56      IF( ln_flork4 ) THEN   ;   CALL flo_4rk( kt, Kbb, Kmm )  ! Trajectories using a 4th order Runge Kutta scheme 
     57      ELSE                   ;   CALL flo_blk( kt, Kbb, Kmm )  ! Trajectories using Blanke' algorithme 
    5758      ENDIF 
    5859      ! 
    5960      IF( lk_mpp )   CALL mppsync   ! synchronization of all the processor 
    6061      ! 
    61       CALL flo_wri( kt )      ! trajectories ouput  
     62      CALL flo_wri( kt, Kmm ) ! trajectories ouput  
    6263      ! 
    6364      CALL flo_rst( kt )      ! trajectories restart 
    6465      ! 
    65       wb(:,:,:) = wn(:,:,:)         ! Save the old vertical velocity field 
     66      wb(:,:,:) = ww(:,:,:)         ! Save the old vertical velocity field 
    6667      ! 
    6768      IF( ln_timing )   CALL timing_stop('flo_stp') 
     
    7071 
    7172 
    72    SUBROUTINE flo_init 
     73   SUBROUTINE flo_init( Kmm ) 
    7374      !!---------------------------------------------------------------- 
    7475      !!                 ***  ROUTINE flo_init  *** 
     
    7677      !! ** Purpose :   Read the namelist of floats 
    7778      !!---------------------------------------------------------------------- 
     79      INTEGER, INTENT(in) :: Kmm       ! ocean time level index 
     80      ! 
    7881      INTEGER ::   jfl 
    7982      INTEGER ::   ios                 ! Local integer output status for namelist read 
     
    8689      IF(lwp) WRITE(numout,*) '~~~~~~~' 
    8790 
    88       REWIND( numnam_ref )              ! Namelist namflo in reference namelist : Floats 
    8991      READ  ( numnam_ref, namflo, IOSTAT = ios, ERR = 901) 
    9092901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namflo in reference namelist' ) 
    9193 
    92       REWIND( numnam_cfg )              ! Namelist namflo in configuration namelist : Floats 
    9394      READ  ( numnam_cfg, namflo, IOSTAT = ios, ERR = 902 ) 
    9495902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namflo in configuration namelist' ) 
     
    130131         END DO 
    131132         ! 
    132          CALL flo_dom                  ! compute/read initial position of floats 
     133         CALL flo_dom( Kmm )           ! compute/read initial position of floats 
    133134         ! 
    134          wb(:,:,:) = wn(:,:,:)         ! set wb for computation of floats trajectories at the first time step 
     135         wb(:,:,:) = ww(:,:,:)         ! set wb for computation of floats trajectories at the first time step 
    135136         ! 
    136137      ENDIF 
    137       ! 
    138138   END SUBROUTINE flo_init 
    139139 
  • NEMO/trunk/src/OCE/FLO/floblk.F90

    r11536 r12377  
    2727CONTAINS 
    2828 
    29    SUBROUTINE flo_blk( kt ) 
     29   SUBROUTINE flo_blk( kt, Kbb, Kmm ) 
    3030      !!--------------------------------------------------------------------- 
    3131      !!                  ***  ROUTINE flo_blk  *** 
     
    3838      !!      of the floats and the grid defined on the domain. 
    3939      !!---------------------------------------------------------------------- 
    40       INTEGER, INTENT( in  ) ::   kt ! ocean time step 
     40      INTEGER, INTENT( in  ) ::   kt       ! ocean time step 
     41      INTEGER, INTENT( in  ) ::   Kbb, Kmm ! ocean time level indices 
    4142      !! 
    4243      INTEGER :: jfl              ! dummy loop arguments 
     
    110111            ! compute the transport across the mesh where the float is.             
    111112!!bug (gm) change e3t into e3. but never checked  
    112             zsurfx(1) = e2u(iiloc(jfl)-1,ijloc(jfl)  ) * e3u_n(iiloc(jfl)-1,ijloc(jfl)  ,-ikl(jfl)) 
    113             zsurfx(2) = e2u(iiloc(jfl)  ,ijloc(jfl)  ) * e3u_n(iiloc(jfl)  ,ijloc(jfl)  ,-ikl(jfl)) 
    114             zsurfy(1) = e1v(iiloc(jfl)  ,ijloc(jfl)-1) * e3v_n(iiloc(jfl)  ,ijloc(jfl)-1,-ikl(jfl)) 
    115             zsurfy(2) = e1v(iiloc(jfl)  ,ijloc(jfl)  ) * e3v_n(iiloc(jfl)  ,ijloc(jfl)  ,-ikl(jfl)) 
     113            zsurfx(1) = e2u(iiloc(jfl)-1,ijloc(jfl)  ) * e3u(iiloc(jfl)-1,ijloc(jfl)  ,-ikl(jfl),Kmm) 
     114            zsurfx(2) = e2u(iiloc(jfl)  ,ijloc(jfl)  ) * e3u(iiloc(jfl)  ,ijloc(jfl)  ,-ikl(jfl),Kmm) 
     115            zsurfy(1) = e1v(iiloc(jfl)  ,ijloc(jfl)-1) * e3v(iiloc(jfl)  ,ijloc(jfl)-1,-ikl(jfl),Kmm) 
     116            zsurfy(2) = e1v(iiloc(jfl)  ,ijloc(jfl)  ) * e3v(iiloc(jfl)  ,ijloc(jfl)  ,-ikl(jfl),Kmm) 
    116117 
    117118            ! for a isobar float zsurfz is put to zero. The vertical velocity will be zero too. 
    118119            zsurfz =          e1e2t(iiloc(jfl),ijloc(jfl)) 
    119             zvol   = zsurfz * e3t_n(iiloc(jfl),ijloc(jfl),-ikl(jfl)) 
     120            zvol   = zsurfz * e3t(iiloc(jfl),ijloc(jfl),-ikl(jfl),Kmm) 
    120121 
    121122            ! 
    122             zuinfl =( ub(iiloc(jfl)-1,ijloc(jfl),-ikl(jfl)) + un(iiloc(jfl)-1,ijloc(jfl),-ikl(jfl)) )/2.*zsurfx(1) 
    123             zuoutfl=( ub(iiloc(jfl)  ,ijloc(jfl),-ikl(jfl)) + un(iiloc(jfl)  ,ijloc(jfl),-ikl(jfl)) )/2.*zsurfx(2) 
    124             zvinfl =( vb(iiloc(jfl),ijloc(jfl)-1,-ikl(jfl)) + vn(iiloc(jfl),ijloc(jfl)-1,-ikl(jfl)) )/2.*zsurfy(1) 
    125             zvoutfl=( vb(iiloc(jfl),ijloc(jfl)  ,-ikl(jfl)) + vn(iiloc(jfl),ijloc(jfl)  ,-ikl(jfl)) )/2.*zsurfy(2) 
     123            zuinfl =( uu(iiloc(jfl)-1,ijloc(jfl),-ikl(jfl),Kbb) + uu(iiloc(jfl)-1,ijloc(jfl),-ikl(jfl),Kmm) )/2.*zsurfx(1) 
     124            zuoutfl=( uu(iiloc(jfl)  ,ijloc(jfl),-ikl(jfl),Kbb) + uu(iiloc(jfl)  ,ijloc(jfl),-ikl(jfl),Kmm) )/2.*zsurfx(2) 
     125            zvinfl =( vv(iiloc(jfl),ijloc(jfl)-1,-ikl(jfl),Kbb) + vv(iiloc(jfl),ijloc(jfl)-1,-ikl(jfl),Kmm) )/2.*zsurfy(1) 
     126            zvoutfl=( vv(iiloc(jfl),ijloc(jfl)  ,-ikl(jfl),Kbb) + vv(iiloc(jfl),ijloc(jfl)  ,-ikl(jfl),Kmm) )/2.*zsurfy(2) 
    126127            zwinfl =-(wb(iiloc(jfl),ijloc(jfl),-(ikl(jfl)-1))    & 
    127                &   +  wn(iiloc(jfl),ijloc(jfl),-(ikl(jfl)-1)) )/2. *  zsurfz*nisobfl(jfl) 
     128               &   +  ww(iiloc(jfl),ijloc(jfl),-(ikl(jfl)-1)) )/2. *  zsurfz*nisobfl(jfl) 
    128129            zwoutfl=-(wb(iiloc(jfl),ijloc(jfl),- ikl(jfl)   )   & 
    129                &   +  wn(iiloc(jfl),ijloc(jfl),- ikl(jfl)   ) )/2. *  zsurfz*nisobfl(jfl) 
     130               &   +  ww(iiloc(jfl),ijloc(jfl),- ikl(jfl)   ) )/2. *  zsurfz*nisobfl(jfl) 
    130131             
    131132            ! interpolation of velocity field on the float initial position             
  • NEMO/trunk/src/OCE/FLO/flodom.F90

    r11818 r12377  
    4040CONTAINS 
    4141 
    42    SUBROUTINE flo_dom 
     42   SUBROUTINE flo_dom( Kmm ) 
    4343      !! --------------------------------------------------------------------- 
    4444      !!                  ***  ROUTINE flo_dom  *** 
     
    4949      !!               the longitude (degree) and the depth (m). 
    5050      !!----------------------------------------------------------------------       
     51      INTEGER, INTENT(in) ::  Kmm    ! ocean time level index 
     52      ! 
    5153      INTEGER            ::   jfl    ! dummy loop   
    5254      INTEGER            ::   inum   ! logical unit for file read 
     
    9092                CALL flo_add_new_ariane_floats(jpnrstflo+1,jpnfl)  
    9193            ELSE                 !Add new floats with long/lat convention 
    92                 CALL flo_add_new_floats(jpnrstflo+1,jpnfl) 
     94                CALL flo_add_new_floats(Kmm,jpnrstflo+1,jpnfl) 
    9395            ENDIF 
    9496         ENDIF 
     
    102104            CALL flo_add_new_ariane_floats(1,jpnfl) 
    103105         ELSE                      !Add new floats with long/lat convention 
    104             CALL flo_add_new_floats(1,jpnfl) 
     106            CALL flo_add_new_floats(Kmm,1,jpnfl) 
    105107         ENDIF 
    106108 
     
    109111   END SUBROUTINE flo_dom 
    110112 
    111    SUBROUTINE flo_add_new_floats(kfl_start, kfl_end) 
     113   SUBROUTINE flo_add_new_floats(Kmm, kfl_start, kfl_end) 
    112114      !! ------------------------------------------------------------- 
    113115      !!                 ***  SUBROUTINE add_new_arianefloats  *** 
     
    124126      !! ** Method  :  
    125127      !!---------------------------------------------------------------------- 
     128      INTEGER, INTENT(in) :: Kmm 
    126129      INTEGER, INTENT(in) :: kfl_start, kfl_end 
    127130      !! 
     
    170173                  ihtest(jfl) = ihtest(jfl)+1 
    171174                  DO jk = 1, jpk-1 
    172                      IF( (gdepw_n(ji,jj,jk) <= flzz(jfl)) .AND. (gdepw_n(ji,jj,jk+1) > flzz(jfl)) ) THEN 
     175                     IF( (gdepw(ji,jj,jk,Kmm) <= flzz(jfl)) .AND. (gdepw(ji,jj,jk+1,Kmm) > flzz(jfl)) ) THEN 
    173176                        ikmfl(jfl) = jk 
    174177                        ivtest(jfl) = ivtest(jfl) + 1 
     
    232235            zgifl(jfl)= (iimfl(jfl)-0.5) + zdxab/e1u(iimfl(jfl)-1,ijmfl(jfl)) + (mig(1)-1) 
    233236            zgjfl(jfl)= (ijmfl(jfl)-0.5) + zdyad/e2v(iimfl(jfl),ijmfl(jfl)-1) + (mjg(1)-1) 
    234             zgkfl(jfl) = (( gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1) - flzz(jfl) )* ikmfl(jfl))   & 
    235                &                 / (  gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1)                              & 
    236                &                    - gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl) ) )                             & 
    237                &                 + (( flzz(jfl)-gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)) ) *(ikmfl(jfl)+1))   & 
    238                &                 / (  gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1)                              & 
    239                &                    - gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)) ) 
     237            zgkfl(jfl) = (( gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1,Kmm) - flzz(jfl) )* ikmfl(jfl))   & 
     238               &                 / (  gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1,Kmm)                              & 
     239               &                    - gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl) ,Kmm) )                             & 
     240               &                 + (( flzz(jfl)-gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl),Kmm) ) *(ikmfl(jfl)+1))   & 
     241               &                 / (  gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1,Kmm)                              & 
     242               &                    - gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl),Kmm) ) 
    240243         ELSE 
    241244            zgifl(jfl) = 0.e0 
  • NEMO/trunk/src/OCE/FLO/flowri.F90

    r11818 r12377  
    5151   END FUNCTION flo_wri_alloc 
    5252 
    53    SUBROUTINE flo_wri( kt ) 
     53   SUBROUTINE flo_wri( kt, Kmm ) 
    5454      !!--------------------------------------------------------------------- 
    5555      !!                  ***  ROUTINE flo_wri *** 
     
    6464      !!---------------------------------------------------------------------- 
    6565      !! * Arguments 
    66       INTEGER  :: kt                               ! time step 
     66      INTEGER, INTENT(in)  :: kt                 ! time step 
     67      INTEGER, INTENT(in)  :: Kmm                ! time level index 
    6768 
    6869      !! * Local declarations 
     
    116117               zlon(jfl) = (1.-zafl)*(1.-zbfl)*glamt(iafloc ,ibfloc ) + (1.-zafl) * zbfl * glamt(iafloc ,ib1floc)   & 
    117118                     +     zafl *(1.-zbfl)*glamt(ia1floc,ibfloc ) +     zafl  * zbfl * glamt(ia1floc,ib1floc) 
    118                zdep(jfl) = (1.-zcfl)*gdepw_n(iafloc,ibfloc,icfl ) + zcfl * gdepw_n(iafloc,ibfloc,ic1fl)      
     119               zdep(jfl) = (1.-zcfl)*gdepw(iafloc,ibfloc,icfl ,Kmm) + zcfl * gdepw(iafloc,ibfloc,ic1fl,Kmm)      
    119120 
    120121               !save temperature, salinity and density at this position 
    121                ztem(jfl) = tsn(iafloc,ibfloc,icfl,jp_tem) 
    122                zsal (jfl) = tsn(iafloc,ibfloc,icfl,jp_sal) 
     122               ztem(jfl) = ts(iafloc,ibfloc,icfl,jp_tem,Kmm) 
     123               zsal (jfl) = ts(iafloc,ibfloc,icfl,jp_sal,Kmm) 
    123124               zrho (jfl) = (rhd(iafloc,ibfloc,icfl)+1)*rau0 
    124125 
     
    137138            zlon(jfl) = (1.-zafl)*(1.-zbfl)*glamt(iafloc ,ibfloc ) + (1.-zafl) * zbfl * glamt(iafloc ,ib1floc)   & 
    138139                      +     zafl *(1.-zbfl)*glamt(ia1floc,ibfloc ) +     zafl  * zbfl * glamt(ia1floc,ib1floc) 
    139             zdep(jfl) = (1.-zcfl)*gdepw_n(iafloc,ibfloc,icfl ) + zcfl * gdepw_n(iafloc,ibfloc,ic1fl) 
    140  
    141             ztem(jfl) = tsn(iafloc,ibfloc,icfl,jp_tem) 
    142             zsal(jfl) = tsn(iafloc,ibfloc,icfl,jp_sal) 
     140            zdep(jfl) = (1.-zcfl)*gdepw(iafloc,ibfloc,icfl ,Kmm) + zcfl * gdepw(iafloc,ibfloc,ic1fl,Kmm) 
     141 
     142            ztem(jfl) = ts(iafloc,ibfloc,icfl,jp_tem,Kmm) 
     143            zsal(jfl) = ts(iafloc,ibfloc,icfl,jp_sal,Kmm) 
    143144            zrho(jfl) = (rhd(iafloc,ibfloc,icfl)+1)*rau0 
    144145           
Note: See TracChangeset for help on using the changeset viewer.