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 10970 for NEMO/branches – NEMO

Changeset 10970 for NEMO/branches


Ignore:
Timestamp:
2019-05-13T14:02:19+02:00 (5 years ago)
Author:
davestorkey
Message:

2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps : CRS and FLO. Only tested compilation. Note that base code doesn't compile with key_floats (#2279), so changes to FLO not really tested at all.

Location:
NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE
Files:
9 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/CRS/crsfld.F90

    r10425 r10970  
    4040CONTAINS 
    4141 
    42    SUBROUTINE crs_fld( kt ) 
     42   SUBROUTINE crs_fld( kt, Kmm ) 
    4343      !!--------------------------------------------------------------------- 
    4444      !!                  ***  ROUTINE crs_fld  *** 
     
    5454      !!---------------------------------------------------------------------- 
    5555      INTEGER, INTENT(in) ::   kt   ! ocean time-step index 
     56      INTEGER, INTENT(in) ::   Kmm  ! time level index 
    5657      ! 
    5758      INTEGER  ::   ji, jj, jk        ! dummy loop indices 
     
    6768 
    6869      ! Depth work arrrays 
    69       ze3t(:,:,:) = e3t_n(:,:,:) 
    70       ze3u(:,:,:) = e3u_n(:,:,:) 
    71       ze3v(:,:,:) = e3v_n(:,:,:) 
    72       ze3w(:,:,:) = e3w_n(:,:,:) 
     70      ze3t(:,:,:) = e3t(:,:,:,Kmm) 
     71      ze3u(:,:,:) = e3u(:,:,:,Kmm) 
     72      ze3v(:,:,:) = e3v(:,:,:,Kmm) 
     73      ze3w(:,:,:) = e3w(:,:,:,Kmm) 
    7374 
    7475      IF( kt == nit000  ) THEN 
     
    9697 
    9798      !  Temperature 
    98       zt(:,:,:) = tsn(:,:,:,jp_tem)  ;      zt_crs(:,:,:) = 0._wp 
     99      zt(:,:,:) = ts(:,:,:,jp_tem,Kmm)  ;      zt_crs(:,:,:) = 0._wp 
    99100      CALL crs_dom_ope( zt, 'VOL', 'T', tmask, zt_crs, p_e12=e1e2t, p_e3=ze3t, psgn=1.0 ) 
    100101      tsn_crs(:,:,:,jp_tem) = zt_crs(:,:,:) 
     
    105106       
    106107      !  Salinity 
    107       zs(:,:,:) = tsn(:,:,:,jp_sal)  ;      zs_crs(:,:,:) = 0._wp 
     108      zs(:,:,:) = ts(:,:,:,jp_sal,Kmm)  ;      zs_crs(:,:,:) = 0._wp 
    108109      CALL crs_dom_ope( zs, 'VOL', 'T', tmask, zs_crs, p_e12=e1e2t, p_e3=ze3t, psgn=1.0 ) 
    109110      tsn_crs(:,:,:,jp_sal) = zt_crs(:,:,:) 
     
    113114 
    114115      !  U-velocity 
    115       CALL crs_dom_ope( un, 'SUM', 'U', umask, un_crs, p_e12=e2u, p_e3=ze3u, p_surf_crs=e2e3u_msk, psgn=-1.0 ) 
     116      CALL crs_dom_ope( uu(:,:,:,Kmm), 'SUM', 'U', umask, un_crs, p_e12=e2u, p_e3=ze3u, p_surf_crs=e2e3u_msk, psgn=-1.0 ) 
    116117      ! 
    117118      zt(:,:,:) = 0._wp     ;    zs(:,:,:) = 0._wp  ;   zt_crs(:,:,:) = 0._wp   ;    zs_crs(:,:,:) = 0._wp 
     
    119120         DO jj = 2, jpjm1 
    120121            DO ji = 2, jpim1    
    121                zt(ji,jj,jk)  = un(ji,jj,jk) * 0.5 * ( tsn(ji,jj,jk,jp_tem) + tsn(ji+1,jj,jk,jp_tem) )  
    122                zs(ji,jj,jk)  = un(ji,jj,jk) * 0.5 * ( tsn(ji,jj,jk,jp_sal) + tsn(ji+1,jj,jk,jp_sal) )  
     122               zt(ji,jj,jk)  = uu(ji,jj,jk,Kmm) * 0.5 * ( ts(ji,jj,jk,jp_tem,Kmm) + ts(ji+1,jj,jk,jp_tem,Kmm) )  
     123               zs(ji,jj,jk)  = uu(ji,jj,jk,Kmm) * 0.5 * ( ts(ji,jj,jk,jp_sal,Kmm) + ts(ji+1,jj,jk,jp_sal,Kmm) )  
    123124            END DO 
    124125         END DO 
     
    132133 
    133134      !  V-velocity 
    134       CALL crs_dom_ope( vn, 'SUM', 'V', vmask, vn_crs, p_e12=e1v, p_e3=ze3v, p_surf_crs=e1e3v_msk, psgn=-1.0 ) 
     135      CALL crs_dom_ope( vv(:,:,:,Kmm), 'SUM', 'V', vmask, vn_crs, p_e12=e1v, p_e3=ze3v, p_surf_crs=e1e3v_msk, psgn=-1.0 ) 
    135136      !                                                                                  
    136137      zt(:,:,:) = 0._wp     ;    zs(:,:,:) = 0._wp  ;   zt_crs(:,:,:) = 0._wp   ;    zs_crs(:,:,:) = 0._wp 
     
    138139         DO jj = 2, jpjm1 
    139140            DO ji = 2, jpim1    
    140                zt(ji,jj,jk)  = vn(ji,jj,jk) * 0.5 * ( tsn(ji,jj,jk,jp_tem) + tsn(ji,jj+1,jk,jp_tem) )  
    141                zs(ji,jj,jk)  = vn(ji,jj,jk) * 0.5 * ( tsn(ji,jj,jk,jp_sal) + tsn(ji,jj+1,jk,jp_sal) )  
     141               zt(ji,jj,jk)  = vv(ji,jj,jk,Kmm) * 0.5 * ( ts(ji,jj,jk,jp_tem,Kmm) + ts(ji,jj+1,jk,jp_tem,Kmm) )  
     142               zs(ji,jj,jk)  = vv(ji,jj,jk,Kmm) * 0.5 * ( ts(ji,jj,jk,jp_sal,Kmm) + ts(ji,jj+1,jk,jp_sal,Kmm) )  
    142143            END DO 
    143144         END DO 
     
    155156            DO jj = 2, jpjm1 
    156157               DO ji = fs_2, fs_jpim1   ! vector opt. 
    157                   zztmp  = r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk) 
     158                  zztmp  = r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm) 
    158159                  z3d(ji,jj,jk) = 0.25_wp * zztmp * (                                    & 
    159                      &            un(ji-1,jj,jk)**2 * e2u(ji-1,jj) * e3u_n(ji-1,jj,jk)   & 
    160                      &          + un(ji  ,jj,jk)**2 * e2u(ji  ,jj) * e3u_n(ji  ,jj,jk)   & 
    161                      &          + vn(ji,jj-1,jk)**2 * e1v(ji,jj-1) * e3v_n(ji,jj-1,jk)   & 
    162                      &          + vn(ji,jj  ,jk)**2 * e1v(ji,jj  ) * e3v_n(ji,jj  ,jk)   ) 
     160                     &            uu(ji-1,jj,jk,Kmm)**2 * e2u(ji-1,jj) * e3u(ji-1,jj,jk,Kmm)   & 
     161                     &          + uu(ji  ,jj,jk,Kmm)**2 * e2u(ji  ,jj) * e3u(ji  ,jj,jk,Kmm)   & 
     162                     &          + vv(ji,jj-1,jk,Kmm)**2 * e1v(ji,jj-1) * e3v(ji,jj-1,jk,Kmm)   & 
     163                     &          + vv(ji,jj  ,jk,Kmm)**2 * e1v(ji,jj  ) * e3v(ji,jj  ,jk,Kmm)   ) 
    163164               END DO 
    164165            END DO 
     
    191192      !  W-velocity 
    192193      IF( ln_crs_wn ) THEN 
    193          CALL crs_dom_ope( wn, 'SUM', 'W', tmask, wn_crs, p_e12=e1e2t, p_surf_crs=e1e2w_msk, psgn=1.0 ) 
    194        !  CALL crs_dom_ope( wn, 'VOL', 'W', tmask, wn_crs, p_e12=e1e2t, p_e3=ze3w ) 
     194         CALL crs_dom_ope( ww, 'SUM', 'W', tmask, wn_crs, p_e12=e1e2t, p_surf_crs=e1e2w_msk, psgn=1.0 ) 
     195       !  CALL crs_dom_ope( ww, 'VOL', 'W', tmask, wn_crs, p_e12=e1e2t, p_e3=ze3w ) 
    195196      ELSE 
    196197        wn_crs(:,:,jpk) = 0._wp 
     
    219220       
    220221      !  sbc fields   
    221       CALL crs_dom_ope( sshn , 'VOL', 'T', tmask, sshn_crs , p_e12=e1e2t, p_e3=ze3t           , psgn=1.0 )   
     222      CALL crs_dom_ope( ssh(:,:,Kmm) , 'VOL', 'T', tmask, sshn_crs , p_e12=e1e2t, p_e3=ze3t           , psgn=1.0 )   
    222223      CALL crs_dom_ope( utau , 'SUM', 'U', umask, utau_crs , p_e12=e2u  , p_surf_crs=e2u_crs  , psgn=1.0 ) 
    223224      CALL crs_dom_ope( vtau , 'SUM', 'V', vmask, vtau_crs , p_e12=e1v  , p_surf_crs=e1v_crs  , psgn=1.0 ) 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/CRS/crsini.F90

    r10068 r10970  
    3535CONTAINS 
    3636    
    37    SUBROUTINE crs_init  
     37   SUBROUTINE crs_init( Kmm ) 
    3838      !!------------------------------------------------------------------- 
    3939      !!                     *** SUBROUTINE crs_init 
     
    6868      !!               - Read in pertinent data ? 
    6969      !!------------------------------------------------------------------- 
     70      INTEGER, INTENT(in) :: Kmm   ! time level index 
     71      ! 
    7072      INTEGER  :: ji,jj,jk      ! dummy indices 
    7173      INTEGER  :: ierr                                ! allocation error status 
     
    98100        WRITE(numout,*) '      create a mesh file (=T)               ln_msh_crs = ', ln_msh_crs 
    99101        WRITE(numout,*) '      type of Kz coarsening (0,1,2)         nn_crs_kz  = ', nn_crs_kz 
    100         WRITE(numout,*) '      wn coarsened or computed using hdivn  ln_crs_wn  = ', ln_crs_wn 
     102        WRITE(numout,*) '      ww coarsened or computed using hdiv  ln_crs_wn  = ', ln_crs_wn 
    101103     ENDIF 
    102104               
     
    174176      
    175177     ! 
    176      ze3t(:,:,:) = e3t_n(:,:,:) 
    177      ze3u(:,:,:) = e3u_n(:,:,:) 
    178      ze3v(:,:,:) = e3v_n(:,:,:) 
    179      ze3w(:,:,:) = e3w_n(:,:,:) 
     178     ze3t(:,:,:) = e3t(:,:,:,Kmm) 
     179     ze3u(:,:,:) = e3u(:,:,:,Kmm) 
     180     ze3v(:,:,:) = e3v(:,:,:,Kmm) 
     181     ze3w(:,:,:) = e3w(:,:,:,Kmm) 
    180182 
    181183     !    3.d.2   Surfaces  
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/FLO/flo4rk.F90

    r10068 r10970  
    3535CONTAINS 
    3636 
    37    SUBROUTINE flo_4rk( kt ) 
     37   SUBROUTINE flo_4rk( kt, Kbb, Kmm ) 
    3838      !!---------------------------------------------------------------------- 
    3939      !!                  ***  ROUTINE flo_4rk  *** 
     
    4747      !!       floats and the grid defined on the domain. 
    4848      !!---------------------------------------------------------------------- 
    49       INTEGER, INTENT(in) ::   kt   ! ocean time-step index 
     49      INTEGER, INTENT(in) ::   kt         ! ocean time-step index 
     50      INTEGER, INTENT(in) ::   Kbb, Kmm   ! ocean time level indices 
    5051      !! 
    5152      INTEGER ::  jfl, jind           ! dummy loop indices 
     
    127128       
    128129         ! for each step we compute the compute the velocity with Lagrange interpolation 
    129          CALL flo_interp( zgifl, zgjfl, zgkfl, zufl, zvfl, zwfl, jind ) 
     130         CALL flo_interp( Kbb, Kmm, zgifl, zgjfl, zgkfl, zufl, zvfl, zwfl, jind ) 
    130131          
    131132         ! computation of Runge-Kutta factor 
     
    155156 
    156157 
    157    SUBROUTINE flo_interp( pxt , pyt , pzt ,      & 
     158   SUBROUTINE flo_interp( Kbb, Kmm,              & 
     159      &                   pxt , pyt , pzt ,      & 
    158160      &                   pufl, pvfl, pwfl, ki ) 
    159161      !!---------------------------------------------------------------------- 
     
    167169      !!      integrated with RK method. 
    168170      !!---------------------------------------------------------------------- 
     171      INTEGER                    , INTENT(in   ) ::   Kbb, Kmm           ! ocean time level indices 
    169172      REAL(wp) , DIMENSION(jpnfl), INTENT(in   ) ::   pxt , pyt , pzt    ! position of the float 
    170173      REAL(wp) , DIMENSION(jpnfl), INTENT(  out) ::   pufl, pvfl, pwfl   ! velocity at this position 
     
    248251               DO jind3 = 1, 4 
    249252                  ztufl(jfl,jind1,jind2,jind3) =   & 
    250                      &   (  tcoef1(ki) * ub(iidu(jfl,jind1),ijdu(jfl,jind2),ikdu(jfl,jind3)) +   & 
    251                      &      tcoef2(ki) * un(iidu(jfl,jind1),ijdu(jfl,jind2),ikdu(jfl,jind3)) )   & 
     253                     &   (  tcoef1(ki) * uu(iidu(jfl,jind1),ijdu(jfl,jind2),ikdu(jfl,jind3),Kbb) +   & 
     254                     &      tcoef2(ki) * uu(iidu(jfl,jind1),ijdu(jfl,jind2),ikdu(jfl,jind3),Kmm) )   & 
    252255                     &      / e1u(iidu(jfl,jind1),ijdu(jfl,jind2))  
    253256               END DO 
     
    332335               DO jind3 = 1 ,4 
    333336                  ztvfl(jfl,jind1,jind2,jind3)=   & 
    334                      &   ( tcoef1(ki) * vb(iidv(jfl,jind1),ijdv(jfl,jind2),ikdv(jfl,jind3))  +   & 
    335                      &     tcoef2(ki) * vn(iidv(jfl,jind1),ijdv(jfl,jind2),ikdv(jfl,jind3)) )    &  
     337                     &   ( tcoef1(ki) * vv(iidv(jfl,jind1),ijdv(jfl,jind2),ikdv(jfl,jind3),Kbb)  +   & 
     338                     &     tcoef2(ki) * vv(iidv(jfl,jind1),ijdv(jfl,jind2),ikdv(jfl,jind3),Kmm) )    &  
    336339                     &     / e2v(iidv(jfl,jind1),ijdv(jfl,jind2)) 
    337340               END DO 
     
    424427                  ztwfl(jfl,jind1,jind2,jind3)=   & 
    425428                     &   ( tcoef1(ki) * wb(iidw(jfl,jind1),ijdw(jfl,jind2),ikdw(jfl,jind3))+   & 
    426                      &     tcoef2(ki) * wn(iidw(jfl,jind1),ijdw(jfl,jind2),ikdw(jfl,jind3)) )  & 
    427                      &   / e3w_n(iidw(jfl,jind1),ijdw(jfl,jind2),ikdw(jfl,jind3)) 
     429                     &     tcoef2(ki) * ww(iidw(jfl,jind1),ijdw(jfl,jind2),ikdw(jfl,jind3)) )  & 
     430                     &   / e3w(iidw(jfl,jind1),ijdw(jfl,jind2),ikdw(jfl,jind3),Kmm) 
    428431               END DO 
    429432            END DO 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/FLO/floats.F90

    r10068 r10970  
    3939CONTAINS 
    4040 
    41    SUBROUTINE flo_stp( kt ) 
     41   SUBROUTINE flo_stp( kt, Kbb, Kmm ) 
    4242      !!---------------------------------------------------------------------- 
    4343      !!                   ***  ROUTINE flo_stp  *** 
     
    5050      !!        if ln_flork4 =T 
    5151      !!---------------------------------------------------------------------- 
    52       INTEGER, INTENT( in  ) ::   kt   ! ocean time step 
     52      INTEGER, INTENT( in  ) ::   kt        ! ocean time step 
     53      INTEGER, INTENT( in  ) ::   Kbb, Kmm  ! ocean time level indices  
    5354      !!---------------------------------------------------------------------- 
    5455      ! 
    5556      IF( ln_timing )   CALL timing_start('flo_stp') 
    5657      ! 
    57       IF( ln_flork4 ) THEN   ;   CALL flo_4rk( kt )        ! Trajectories using a 4th order Runge Kutta scheme 
    58       ELSE                   ;   CALL flo_blk( kt )        ! Trajectories using Blanke' algorithme 
     58      IF( ln_flork4 ) THEN   ;   CALL flo_4rk( kt, Kbb, Kmm )  ! Trajectories using a 4th order Runge Kutta scheme 
     59      ELSE                   ;   CALL flo_blk( kt, Kbb, Kmm )  ! Trajectories using Blanke' algorithme 
    5960      ENDIF 
    6061      ! 
    6162      IF( lk_mpp )   CALL mppsync   ! synchronization of all the processor 
    6263      ! 
    63       CALL flo_wri( kt )      ! trajectories ouput  
     64      CALL flo_wri( kt, Kmm ) ! trajectories ouput  
    6465      ! 
    6566      CALL flo_rst( kt )      ! trajectories restart 
    6667      ! 
    67       wb(:,:,:) = wn(:,:,:)         ! Save the old vertical velocity field 
     68      wb(:,:,:) = ww(:,:,:)         ! Save the old vertical velocity field 
    6869      ! 
    6970      IF( ln_timing )   CALL timing_stop('flo_stp') 
     
    7273 
    7374 
    74    SUBROUTINE flo_init 
     75   SUBROUTINE flo_init( Kmm ) 
    7576      !!---------------------------------------------------------------- 
    7677      !!                 ***  ROUTINE flo_init  *** 
     
    7879      !! ** Purpose :   Read the namelist of floats 
    7980      !!---------------------------------------------------------------------- 
     81      INTEGER, INTENT(in) :: Kmm       ! ocean time level index 
     82      ! 
    8083      INTEGER ::   jfl 
    8184      INTEGER ::   ios                 ! Local integer output status for namelist read 
     
    130133      END DO 
    131134      ! 
    132       CALL flo_dom                  ! compute/read initial position of floats 
     135      CALL flo_dom( Kmm )           ! compute/read initial position of floats 
    133136      ! 
    134       wb(:,:,:) = wn(:,:,:)         ! set wb for computation of floats trajectories at the first time step 
     137      wb(:,:,:) = ww(:,:,:)         ! set wb for computation of floats trajectories at the first time step 
    135138      ! 
    136139   END SUBROUTINE flo_init 
     
    141144   !!---------------------------------------------------------------------- 
    142145CONTAINS 
    143    SUBROUTINE flo_stp( kt )          ! Empty routine 
     146   SUBROUTINE flo_stp( kt, Kbb, Kmm )  ! Empty routine 
    144147      IMPLICIT NONE 
    145148      INTEGER, INTENT( in ) :: kt 
     149      INTEGER, INTENT( in ) :: Kbb, Kmm 
    146150      WRITE(*,*) 'flo_stp: You should not have seen this print! error?', kt 
    147151   END SUBROUTINE flo_stp 
    148    SUBROUTINE flo_init          ! Empty routine 
     152   SUBROUTINE flo_init( Kmm )          ! Empty routine 
    149153      IMPLICIT NONE 
     154      INTEGER, INTENT( in ) :: Kmm 
    150155   END SUBROUTINE flo_init 
    151156#endif 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/FLO/floblk.F90

    r10425 r10970  
    2929CONTAINS 
    3030 
    31    SUBROUTINE flo_blk( kt ) 
     31   SUBROUTINE flo_blk( kt, Kbb, Kmm ) 
    3232      !!--------------------------------------------------------------------- 
    3333      !!                  ***  ROUTINE flo_blk  *** 
     
    4040      !!      of the floats and the grid defined on the domain. 
    4141      !!---------------------------------------------------------------------- 
    42       INTEGER, INTENT( in  ) ::   kt ! ocean time step 
     42      INTEGER, INTENT( in  ) ::   kt       ! ocean time step 
     43      INTEGER, INTENT( in  ) ::   Kbb, Kmm ! ocean time level indices 
    4344      !! 
    4445      INTEGER :: jfl              ! dummy loop arguments 
     
    112113            ! compute the transport across the mesh where the float is.             
    113114!!bug (gm) change e3t into e3. but never checked  
    114             zsurfx(1) = e2u(iiloc(jfl)-1,ijloc(jfl)  ) * e3u_n(iiloc(jfl)-1,ijloc(jfl)  ,-ikl(jfl)) 
    115             zsurfx(2) = e2u(iiloc(jfl)  ,ijloc(jfl)  ) * e3u_n(iiloc(jfl)  ,ijloc(jfl)  ,-ikl(jfl)) 
    116             zsurfy(1) = e1v(iiloc(jfl)  ,ijloc(jfl)-1) * e3v_n(iiloc(jfl)  ,ijloc(jfl)-1,-ikl(jfl)) 
    117             zsurfy(2) = e1v(iiloc(jfl)  ,ijloc(jfl)  ) * e3v_n(iiloc(jfl)  ,ijloc(jfl)  ,-ikl(jfl)) 
     115            zsurfx(1) = e2u(iiloc(jfl)-1,ijloc(jfl)  ) * e3u(iiloc(jfl)-1,ijloc(jfl)  ,-ikl(jfl),Kmm) 
     116            zsurfx(2) = e2u(iiloc(jfl)  ,ijloc(jfl)  ) * e3u(iiloc(jfl)  ,ijloc(jfl)  ,-ikl(jfl),Kmm) 
     117            zsurfy(1) = e1v(iiloc(jfl)  ,ijloc(jfl)-1) * e3v(iiloc(jfl)  ,ijloc(jfl)-1,-ikl(jfl),Kmm) 
     118            zsurfy(2) = e1v(iiloc(jfl)  ,ijloc(jfl)  ) * e3v(iiloc(jfl)  ,ijloc(jfl)  ,-ikl(jfl),Kmm) 
    118119 
    119120            ! for a isobar float zsurfz is put to zero. The vertical velocity will be zero too. 
    120121            zsurfz =          e1e2t(iiloc(jfl),ijloc(jfl)) 
    121             zvol   = zsurfz * e3t_n(iiloc(jfl),ijloc(jfl),-ikl(jfl)) 
     122            zvol   = zsurfz * e3t(iiloc(jfl),ijloc(jfl),-ikl(jfl),Kmm) 
    122123 
    123124            ! 
    124             zuinfl =( ub(iiloc(jfl)-1,ijloc(jfl),-ikl(jfl)) + un(iiloc(jfl)-1,ijloc(jfl),-ikl(jfl)) )/2.*zsurfx(1) 
    125             zuoutfl=( ub(iiloc(jfl)  ,ijloc(jfl),-ikl(jfl)) + un(iiloc(jfl)  ,ijloc(jfl),-ikl(jfl)) )/2.*zsurfx(2) 
    126             zvinfl =( vb(iiloc(jfl),ijloc(jfl)-1,-ikl(jfl)) + vn(iiloc(jfl),ijloc(jfl)-1,-ikl(jfl)) )/2.*zsurfy(1) 
    127             zvoutfl=( vb(iiloc(jfl),ijloc(jfl)  ,-ikl(jfl)) + vn(iiloc(jfl),ijloc(jfl)  ,-ikl(jfl)) )/2.*zsurfy(2) 
     125            zuinfl =( uu(iiloc(jfl)-1,ijloc(jfl),-ikl(jfl),Kbb) + uu(iiloc(jfl)-1,ijloc(jfl),-ikl(jfl),Kmm) )/2.*zsurfx(1) 
     126            zuoutfl=( uu(iiloc(jfl)  ,ijloc(jfl),-ikl(jfl),Kbb) + uu(iiloc(jfl)  ,ijloc(jfl),-ikl(jfl),Kmm) )/2.*zsurfx(2) 
     127            zvinfl =( vv(iiloc(jfl),ijloc(jfl)-1,-ikl(jfl),Kbb) + vv(iiloc(jfl),ijloc(jfl)-1,-ikl(jfl),Kmm) )/2.*zsurfy(1) 
     128            zvoutfl=( vv(iiloc(jfl),ijloc(jfl)  ,-ikl(jfl),Kbb) + vv(iiloc(jfl),ijloc(jfl)  ,-ikl(jfl),Kmm) )/2.*zsurfy(2) 
    128129            zwinfl =-(wb(iiloc(jfl),ijloc(jfl),-(ikl(jfl)-1))    & 
    129                &   +  wn(iiloc(jfl),ijloc(jfl),-(ikl(jfl)-1)) )/2. *  zsurfz*nisobfl(jfl) 
     130               &   +  ww(iiloc(jfl),ijloc(jfl),-(ikl(jfl)-1)) )/2. *  zsurfz*nisobfl(jfl) 
    130131            zwoutfl=-(wb(iiloc(jfl),ijloc(jfl),- ikl(jfl)   )   & 
    131                &   +  wn(iiloc(jfl),ijloc(jfl),- ikl(jfl)   ) )/2. *  zsurfz*nisobfl(jfl) 
     132               &   +  ww(iiloc(jfl),ijloc(jfl),- ikl(jfl)   ) )/2. *  zsurfz*nisobfl(jfl) 
    132133             
    133134            ! interpolation of velocity field on the float initial position             
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/FLO/flodom.F90

    r10425 r10970  
    4444CONTAINS 
    4545 
    46    SUBROUTINE flo_dom 
     46   SUBROUTINE flo_dom( Kmm ) 
    4747      !! --------------------------------------------------------------------- 
    4848      !!                  ***  ROUTINE flo_dom  *** 
     
    5353      !!               the longitude (degree) and the depth (m). 
    5454      !!----------------------------------------------------------------------       
     55      INTEGER, INTENT(in) ::  Kmm    ! ocean time level index 
     56      ! 
    5557      INTEGER            ::   jfl    ! dummy loop   
    5658      INTEGER            ::   inum   ! logical unit for file read 
     
    9496                CALL flo_add_new_ariane_floats(jpnrstflo+1,jpnfl)  
    9597            ELSE                 !Add new floats with long/lat convention 
    96                 CALL flo_add_new_floats(jpnrstflo+1,jpnfl) 
     98                CALL flo_add_new_floats(Kmm,jpnrstflo+1,jpnfl) 
    9799            ENDIF 
    98100         ENDIF 
     
    106108            CALL flo_add_new_ariane_floats(1,jpnfl) 
    107109         ELSE                      !Add new floats with long/lat convention 
    108             CALL flo_add_new_floats(1,jpnfl) 
     110            CALL flo_add_new_floats(Kmm,1,jpnfl) 
    109111         ENDIF 
    110112 
     
    113115   END SUBROUTINE flo_dom 
    114116 
    115    SUBROUTINE flo_add_new_floats(kfl_start, kfl_end) 
     117   SUBROUTINE flo_add_new_floats(Kmm, kfl_start, kfl_end) 
    116118      !! ------------------------------------------------------------- 
    117119      !!                 ***  SUBROUTINE add_new_arianefloats  *** 
     
    128130      !! ** Method  :  
    129131      !!---------------------------------------------------------------------- 
     132      INTEGER, INTENT(in) :: Kmm 
    130133      INTEGER, INTENT(in) :: kfl_start, kfl_end 
    131134      !! 
     
    174177                  ihtest(jfl) = ihtest(jfl)+1 
    175178                  DO jk = 1, jpk-1 
    176                      IF( (gdepw_n(ji,jj,jk) <= flzz(jfl)) .AND. (gdepw_n(ji,jj,jk+1) > flzz(jfl)) ) THEN 
     179                     IF( (gdepw(ji,jj,jk,Kmm) <= flzz(jfl)) .AND. (gdepw(ji,jj,jk+1,Kmm) > flzz(jfl)) ) THEN 
    177180                        ikmfl(jfl) = jk 
    178181                        ivtest(jfl) = ivtest(jfl) + 1 
     
    236239            zgifl(jfl)= (iimfl(jfl)-0.5) + zdxab/e1u(iimfl(jfl)-1,ijmfl(jfl)) + (mig(1)-1) 
    237240            zgjfl(jfl)= (ijmfl(jfl)-0.5) + zdyad/e2v(iimfl(jfl),ijmfl(jfl)-1) + (mjg(1)-1) 
    238             zgkfl(jfl) = (( gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1) - flzz(jfl) )* ikmfl(jfl))   & 
    239                &                 / (  gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1)                              & 
    240                &                    - gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl) ) )                             & 
    241                &                 + (( flzz(jfl)-gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)) ) *(ikmfl(jfl)+1))   & 
    242                &                 / (  gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1)                              & 
    243                &                    - gdepw_n(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)) ) 
     241            zgkfl(jfl) = (( gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1,Kmm) - flzz(jfl) )* ikmfl(jfl))   & 
     242               &                 / (  gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1,Kmm)                              & 
     243               &                    - gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl) ,Kmm) )                             & 
     244               &                 + (( flzz(jfl)-gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl),Kmm) ) *(ikmfl(jfl)+1))   & 
     245               &                 / (  gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl)+1,Kmm)                              & 
     246               &                    - gdepw(iimfl(jfl),ijmfl(jfl),ikmfl(jfl),Kmm) ) 
    244247         ELSE 
    245248            zgifl(jfl) = 0.e0 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/FLO/flowri.F90

    r10425 r10970  
    5555   END FUNCTION flo_wri_alloc 
    5656 
    57    SUBROUTINE flo_wri( kt ) 
     57   SUBROUTINE flo_wri( kt, Kmm ) 
    5858      !!--------------------------------------------------------------------- 
    5959      !!                  ***  ROUTINE flo_wri *** 
     
    6868      !!---------------------------------------------------------------------- 
    6969      !! * Arguments 
    70       INTEGER  :: kt                               ! time step 
     70      INTEGER, INTENT(in)  :: kt                 ! time step 
     71      INTEGER, INTENT(in)  :: Kmm                ! time level index 
    7172 
    7273      !! * Local declarations 
     
    120121               zlon(jfl) = (1.-zafl)*(1.-zbfl)*glamt(iafloc ,ibfloc ) + (1.-zafl) * zbfl * glamt(iafloc ,ib1floc)   & 
    121122                     +     zafl *(1.-zbfl)*glamt(ia1floc,ibfloc ) +     zafl  * zbfl * glamt(ia1floc,ib1floc) 
    122                zdep(jfl) = (1.-zcfl)*gdepw_n(iafloc,ibfloc,icfl ) + zcfl * gdepw_n(iafloc,ibfloc,ic1fl)      
     123               zdep(jfl) = (1.-zcfl)*gdepw(iafloc,ibfloc,icfl ,Kmm) + zcfl * gdepw(iafloc,ibfloc,ic1fl,Kmm)      
    123124 
    124125               !save temperature, salinity and density at this position 
    125                ztem(jfl) = tsn(iafloc,ibfloc,icfl,jp_tem) 
    126                zsal (jfl) = tsn(iafloc,ibfloc,icfl,jp_sal) 
     126               ztem(jfl) = ts(iafloc,ibfloc,icfl,jp_tem,Kmm) 
     127               zsal (jfl) = ts(iafloc,ibfloc,icfl,jp_sal,Kmm) 
    127128               zrho (jfl) = (rhd(iafloc,ibfloc,icfl)+1)*rau0 
    128129 
     
    141142            zlon(jfl) = (1.-zafl)*(1.-zbfl)*glamt(iafloc ,ibfloc ) + (1.-zafl) * zbfl * glamt(iafloc ,ib1floc)   & 
    142143                      +     zafl *(1.-zbfl)*glamt(ia1floc,ibfloc ) +     zafl  * zbfl * glamt(ia1floc,ib1floc) 
    143             zdep(jfl) = (1.-zcfl)*gdepw_n(iafloc,ibfloc,icfl ) + zcfl * gdepw_n(iafloc,ibfloc,ic1fl) 
    144  
    145             ztem(jfl) = tsn(iafloc,ibfloc,icfl,jp_tem) 
    146             zsal(jfl) = tsn(iafloc,ibfloc,icfl,jp_sal) 
     144            zdep(jfl) = (1.-zcfl)*gdepw(iafloc,ibfloc,icfl ,Kmm) + zcfl * gdepw(iafloc,ibfloc,ic1fl,Kmm) 
     145 
     146            ztem(jfl) = ts(iafloc,ibfloc,icfl,jp_tem,Kmm) 
     147            zsal(jfl) = ts(iafloc,ibfloc,icfl,jp_sal,Kmm) 
    147148            zrho(jfl) = (rhd(iafloc,ibfloc,icfl)+1)*rau0 
    148149           
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/nemogcm.F90

    r10969 r10970  
    421421                           CALL     wad_init        ! Wetting and drying options 
    422422                           CALL     dom_init("OPA") ! Domain 
    423       IF( ln_crs       )   CALL     crs_init        ! coarsened grid: domain initialization  
     423      IF( ln_crs       )   CALL     crs_init( Nnn ) ! coarsened grid: domain initialization  
    424424      IF( ln_ctl       )   CALL prt_ctl_init        ! Print control 
    425425       
     
    485485      
    486486      !                                      ! Diagnostics 
    487       IF( lk_floats    )   CALL     flo_init    ! drifting Floats 
     487      IF( lk_floats    )   CALL     flo_init( Nnn )    ! drifting Floats 
    488488      IF( ln_diacfl    )   CALL dia_cfl_init    ! Initialise CFL diagnostics 
    489489                           CALL dia_ptr_init    ! Poleward TRansports initialization 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/step.F90

    r10969 r10970  
    218218      ! diagnostics and outputs 
    219219      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    220       IF( lk_floats  )   CALL flo_stp ( kstp )        ! drifting Floats 
    221       IF( ln_diacfl  )   CALL dia_cfl ( kstp, Nnn )   ! Courant number diagnostics 
    222       IF( lk_diahth  )   CALL dia_hth ( kstp, Nnn )   ! Thermocline depth (20 degres isotherm depth) 
    223       IF( lk_diadct  )   CALL dia_dct ( kstp, Nnn )   ! Transports 
    224                          CALL dia_ar5 ( kstp, Nnn )   ! ar5 diag 
    225       IF( lk_diaharm )   CALL dia_harm( kstp, Nnn )   ! Tidal harmonic analysis 
    226                          CALL dia_wri ( kstp, Nnn )   ! ocean model: outputs 
    227       ! 
    228       IF( ln_crs     )   CALL crs_fld       ( kstp )  ! ocean model: online field coarsening & output 
     220      IF( lk_floats  )   CALL flo_stp ( kstp, Nbb, Nnn )   ! drifting Floats 
     221      IF( ln_diacfl  )   CALL dia_cfl ( kstp,      Nnn )   ! Courant number diagnostics 
     222      IF( lk_diahth  )   CALL dia_hth ( kstp,      Nnn )   ! Thermocline depth (20 degres isotherm depth) 
     223      IF( lk_diadct  )   CALL dia_dct ( kstp,      Nnn )   ! Transports 
     224                         CALL dia_ar5 ( kstp,      Nnn )   ! ar5 diag 
     225      IF( lk_diaharm )   CALL dia_harm( kstp,      Nnn )   ! Tidal harmonic analysis 
     226                         CALL dia_wri ( kstp,      Nnn )   ! ocean model: outputs 
     227      ! 
     228      IF( ln_crs     )   CALL crs_fld ( kstp,      Nnn )   ! ocean model: online field coarsening & output 
    229229       
    230230#if defined key_top 
Note: See TracChangeset for help on using the changeset viewer.