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 4596 for branches/2014/dev_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/DOM – NEMO

Ignore:
Timestamp:
2014-03-26T12:02:30+01:00 (10 years ago)
Author:
gm
Message:

#1260: LDF simplification + bilap iso-neutral for TRA and GYRE

Location:
branches/2014/dev_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/DOM
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • branches/2014/dev_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/DOM/domzgr.F90

    r4292 r4596  
    622622      !!              - update bathy : meter bathymetry (in meters) 
    623623      !!---------------------------------------------------------------------- 
    624       !! 
    625624      INTEGER ::   ji, jj, jl                    ! dummy loop indices 
    626625      INTEGER ::   icompt, ibtest, ikmax         ! temporary integers 
    627626      REAL(wp), POINTER, DIMENSION(:,:) ::  zbathy 
    628  
    629627      !!---------------------------------------------------------------------- 
    630628      ! 
     
    11151113      !! 
    11161114      !!---------------------------------------------------------------------- 
    1117       ! 
    11181115      INTEGER  ::   ji, jj, jk, jl           ! dummy loop argument 
    1119       INTEGER  ::   iip1, ijp1, iim1, ijm1   ! temporary integers 
     1116      INTEGER  ::   iip1, ijp1, iim1, ijm1   ! local integers 
    11201117      INTEGER  ::   ios                      ! Local integer output status for namelist read 
    1121       REAL(wp) ::   zrmax, ztaper   ! temporary scalars 
    1122       REAL(wp) ::   zrfact 
     1118      REAL(wp) ::   zrmax, ztaper, zrfact    ! local scalars 
    11231119      ! 
    11241120      REAL(wp), POINTER, DIMENSION(:,:  ) :: ztmpi1, ztmpi2, ztmpj1, ztmpj2 
     
    12831279         DO jj = 1, jpj 
    12841280            DO ji = 1, jpi 
    1285                ztaper = EXP( -(gphit(ji,jj)/8._wp)**2._wp ) 
     1281               ztaper = EXP( -(gphit(ji,jj)/8._wp)**2 ) 
    12861282               hbatt(ji,jj) = rn_sbot_max * ztaper + hbatt(ji,jj) * ( 1._wp - ztaper ) 
    12871283            END DO 
     
    15521548   END SUBROUTINE zgr_sco 
    15531549 
    1554 !!====================================================================== 
     1550 
    15551551   SUBROUTINE s_sh94() 
    1556  
    15571552      !!---------------------------------------------------------------------- 
    15581553      !!                  ***  ROUTINE s_sh94  *** 
     
    15651560      !! Reference : Song and Haidvogel 1994.  
    15661561      !!---------------------------------------------------------------------- 
    1567       ! 
    15681562      INTEGER  ::   ji, jj, jk           ! dummy loop argument 
    15691563      REAL(wp) ::   zcoeft, zcoefw   ! temporary scalars 
     
    16511645   END SUBROUTINE s_sh94 
    16521646 
     1647 
    16531648   SUBROUTINE s_sf12 
    1654  
    16551649      !!---------------------------------------------------------------------- 
    16561650      !!                  ***  ROUTINE s_sf12 ***  
     
    16661660      !! 
    16671661      !! 
    1668       !! Reference : Siddorn and Furner 2012 (submitted Ocean modelling). 
    1669       !!---------------------------------------------------------------------- 
    1670       ! 
     1662      !! Reference : Siddorn and Furner 2013 (Ocean modelling). 
     1663      !!---------------------------------------------------------------------- 
    16711664      INTEGER  ::   ji, jj, jk           ! dummy loop argument 
    16721665      REAL(wp) ::   zsmth               ! smoothing around critical depth 
     
    16741667      ! 
    16751668      REAL(wp), POINTER, DIMENSION(:,:,:) :: z_gsigw3, z_gsigt3, z_gsi3w3 
    1676       REAL(wp), POINTER, DIMENSION(:,:,:) :: z_esigt3, z_esigw3, z_esigtu3, z_esigtv3, z_esigtf3, z_esigwu3, z_esigwv3            
    1677  
     1669      REAL(wp), POINTER, DIMENSION(:,:,:) :: z_esigt3, z_esigw3, z_esigtu3, z_esigtv3, z_esigtf3, z_esigwu3, z_esigwv3 
     1670      !!---------------------------------------------------------------------- 
    16781671      ! 
    16791672      CALL wrk_alloc( jpi, jpj, jpk, z_gsigw3, z_gsigt3, z_gsi3w3                                      ) 
     
    17441737          END DO 
    17451738 
    1746         ENDDO   ! for all jj's 
    1747       ENDDO    ! for all ji's 
     1739        END DO   ! for all jj's 
     1740      END DO    ! for all ji's 
    17481741 
    17491742      DO ji=1,jpi-1 
     
    17731766          END DO 
    17741767 
    1775         ENDDO 
    1776       ENDDO 
     1768        END DO 
     1769      END DO 
    17771770      ! 
    17781771      CALL lbc_lnk(e3t_0 ,'T',1.) ; CALL lbc_lnk(e3u_0 ,'T',1.) 
     
    17881781   END SUBROUTINE s_sf12 
    17891782 
     1783 
    17901784   SUBROUTINE s_tanh() 
    1791  
    17921785      !!---------------------------------------------------------------------- 
    17931786      !!                  ***  ROUTINE s_tanh***  
     
    17991792      !! Reference : Madec, Lott, Delecluse and Crepon, 1996. JPO, 26, 1393-1408. 
    18001793      !!---------------------------------------------------------------------- 
    1801  
    18021794      INTEGER  ::   ji, jj, jk           ! dummy loop argument 
    18031795      REAL(wp) ::   zcoeft, zcoefw   ! temporary scalars 
    1804  
     1796      ! 
    18051797      REAL(wp), POINTER, DIMENSION(:) :: z_gsigw, z_gsigt, z_gsi3w 
    18061798      REAL(wp), POINTER, DIMENSION(:) :: z_esigt, z_esigw 
     1799      !!---------------------------------------------------------------------- 
    18071800 
    18081801      CALL wrk_alloc( jpk, z_gsigw, z_gsigt, z_gsi3w                                      ) 
     
    18621855   END SUBROUTINE s_tanh 
    18631856 
     1857 
    18641858   FUNCTION fssig( pk ) RESULT( pf ) 
    18651859      !!---------------------------------------------------------------------- 
     
    19321926      REAL(wp), INTENT(in   ) ::   pk1(jpk)       ! continuous "k" coordinate 
    19331927      REAL(wp)                ::   p_gamma(jpk)   ! stretched coordinate 
    1934       REAL(wp), INTENT(in   ) ::   pzb           ! Bottom box depth 
    1935       REAL(wp), INTENT(in   ) ::   pzs           ! surface box depth 
    1936       REAL(wp), INTENT(in   ) ::   psmth       ! Smoothing parameter 
     1928      REAL(wp), INTENT(in   ) ::   pzb            ! Bottom box depth 
     1929      REAL(wp), INTENT(in   ) ::   pzs            ! surface box depth 
     1930      REAL(wp), INTENT(in   ) ::   psmth          ! Smoothing parameter 
     1931      ! 
     1932      INTEGER                 ::   jk 
    19371933      REAL(wp)                ::   za1,za2,za3    ! local variables 
    19381934      REAL(wp)                ::   zn1,zn2        ! local variables 
    19391935      REAL(wp)                ::   za,zb,zx       ! local variables 
    1940       integer                 ::   jk 
    1941       !!---------------------------------------------------------------------- 
    1942       ! 
    1943  
     1936      !!---------------------------------------------------------------------- 
     1937      ! 
    19441938      zn1  =  1./(jpk-1.) 
    19451939      zn2  =  1. -  zn1 
    1946  
     1940      ! 
    19471941      za1 = (rn_alpha+2.0_wp)*zn1**(rn_alpha+1.0_wp)-(rn_alpha+1.0_wp)*zn1**(rn_alpha+2.0_wp)  
    19481942      za2 = (rn_alpha+2.0_wp)*zn2**(rn_alpha+1.0_wp)-(rn_alpha+1.0_wp)*zn2**(rn_alpha+2.0_wp) 
    19491943      za3 = (zn2**3.0_wp - za2)/( zn1**3.0_wp - za1) 
    1950       
     1944      ! 
    19511945      za = pzb - za3*(pzs-za1)-za2 
    19521946      za = za/( zn2-0.5_wp*(za2+zn2**2.0_wp) - za3*(zn1-0.5_wp*(za1+zn1**2.0_wp) ) ) 
    19531947      zb = (pzs - za1 - za*( zn1-0.5_wp*(za1+zn1**2.0_wp ) ) ) / (zn1**3.0_wp - za1) 
    19541948      zx = 1.0_wp-za/2.0_wp-zb 
    1955   
     1949      ! 
    19561950      DO jk = 1, jpk 
    19571951        p_gamma(jk) = za*(pk1(jk)*(1.0_wp-pk1(jk)/2.0_wp))+zb*pk1(jk)**3.0_wp +  & 
     
    19591953                    &      (rn_alpha+1.0_wp)*pk1(jk)**(rn_alpha+2.0_wp) ) 
    19601954        p_gamma(jk) = p_gamma(jk)*psmth+pk1(jk)*(1.0_wp-psmth) 
    1961       ENDDO  
    1962  
     1955      END DO  
    19631956      ! 
    19641957   END FUNCTION fgamma 
  • branches/2014/dev_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90

    r4370 r4596  
    2929   USE daymod          ! calendar 
    3030   USE eosbn2          ! eq. of state, Brunt Vaisala frequency (eos     routine) 
    31    USE ldftra_oce      ! ocean active tracers: lateral physics 
     31   USE ldftra          ! lateral physics: ocean active tracers 
    3232   USE zdf_oce         ! ocean vertical physics 
    3333   USE phycst          ! physical constants 
    3434   USE dtatsd          ! data temperature and salinity   (dta_tsd routine) 
    3535   USE dtauvd          ! data: U & V current             (dta_uvd routine) 
    36    USE in_out_manager  ! I/O manager 
    37    USE iom             ! I/O library 
    3836   USE zpshde          ! partial step: hor. derivative (zps_hde routine) 
    3937   USE eosbn2          ! equation of state            (eos bn2 routine) 
     
    4240   USE dynspg_flt      ! filtered free surface 
    4341   USE sol_oce         ! ocean solver variables 
     42   ! 
     43   USE in_out_manager  ! I/O manager 
     44   USE iom             ! I/O library 
    4445   USE lib_mpp         ! MPP library 
    4546   USE restart         ! restart 
     
    6869      !! ** Purpose :   Initialization of the dynamics and tracer fields. 
    6970      !!---------------------------------------------------------------------- 
    70       ! - ML - needed for initialization of e3t_b 
    71       INTEGER  ::  ji,jj,jk     ! dummy loop indices 
    72       REAL(wp), POINTER, DIMENSION(:,:,:,:)  ::  zuvd    ! U & V data workspace 
     71      INTEGER ::   ji, jj, jk     ! dummy loop indices 
     72      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   zuvd   ! U & V data workspace 
    7373      !!---------------------------------------------------------------------- 
    7474      ! 
    7575      IF( nn_timing == 1 )  CALL timing_start('istate_init') 
    7676      ! 
    77  
    7877      IF(lwp) WRITE(numout,*) 
    7978      IF(lwp) WRITE(numout,*) 'istate_ini : Initialization of the dynamics and tracers' 
    8079      IF(lwp) WRITE(numout,*) '~~~~~~~~~~' 
    8180 
    82       CALL dta_tsd_init                       ! Initialisation of T & S input data 
    83       IF( lk_c1d ) CALL dta_uvd_init          ! Initialization of U & V input data 
    84  
    85       rhd  (:,:,:  ) = 0.e0 
    86       rhop (:,:,:  ) = 0.e0 
    87       rn2  (:,:,:  ) = 0.e0  
    88       tsa  (:,:,:,:) = 0.e0     
     81                     CALL dta_tsd_init        ! Initialisation of T & S input data 
     82      IF( lk_c1d )   CALL dta_uvd_init        ! Initialization of U & V input data 
     83 
     84      rhd  (:,:,:  ) = 0._wp 
     85      rhop (:,:,:  ) = 0._wp 
     86      rn2  (:,:,:  ) = 0._wp 
     87      tsa  (:,:,:,:) = 0._wp    
    8988 
    9089      IF( ln_rstart ) THEN                    ! Restart from a file 
     
    103102         ub   (:,:,:) = 0._wp   ;   un   (:,:,:) = 0._wp 
    104103         vb   (:,:,:) = 0._wp   ;   vn   (:,:,:) = 0._wp   
    105          rotb (:,:,:) = 0._wp   ;   rotn (:,:,:) = 0._wp 
    106          hdivb(:,:,:) = 0._wp   ;   hdivn(:,:,:) = 0._wp 
     104                                    hdivn(:,:,:) = 0._wp 
    107105         ! 
    108106         IF( cp_cfg == 'eel' ) THEN 
     
    158156      ! 
    159157      ! 
    160       un_b(:,:) = 0._wp ; vn_b(:,:) = 0._wp 
    161       ub_b(:,:) = 0._wp ; vb_b(:,:) = 0._wp 
     158      un_b(:,:) = 0._wp   ;  vn_b(:,:) = 0._wp 
     159      ub_b(:,:) = 0._wp   ;  vb_b(:,:) = 0._wp 
    162160      ! 
    163161      DO jk = 1, jpkm1 
    164 #if defined key_vectopt_loop 
    165          DO jj = 1, 1         !Vector opt. => forced unrolling 
    166             DO ji = 1, jpij 
    167 #else  
    168162         DO jj = 1, jpj 
    169163            DO ji = 1, jpi 
    170 #endif                   
    171164               un_b(ji,jj) = un_b(ji,jj) + fse3u_n(ji,jj,jk) * un(ji,jj,jk) * umask(ji,jj,jk) 
    172165               vn_b(ji,jj) = vn_b(ji,jj) + fse3v_n(ji,jj,jk) * vn(ji,jj,jk) * vmask(ji,jj,jk) 
     
    188181      ! 
    189182   END SUBROUTINE istate_init 
     183 
    190184 
    191185   SUBROUTINE istate_t_s 
     
    201195      !! References :  Philander ??? 
    202196      !!---------------------------------------------------------------------- 
    203       INTEGER  :: ji, jj, jk 
    204       REAL(wp) ::   zsal = 35.50 
     197      INTEGER  ::   ji, jj, jk 
     198      REAL(wp) ::   zsal = 35.50_wp 
    205199      !!---------------------------------------------------------------------- 
    206200      ! 
     
    218212      ! 
    219213   END SUBROUTINE istate_t_s 
     214 
    220215 
    221216   SUBROUTINE istate_eel 
     
    231226      !!                and relative vorticity fields 
    232227      !!---------------------------------------------------------------------- 
    233       USE divcur     ! hor. divergence & rel. vorticity      (div_cur routine) 
     228      USE divhor     ! hor. divergence      (div_hor routine) 
    234229      USE iom 
    235   
     230      ! 
    236231      INTEGER  ::   inum              ! temporary logical unit 
    237232      INTEGER  ::   ji, jj, jk        ! dummy loop indices 
     
    280275            tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) 
    281276            ! 
    282             ! set the dynamics: U,V, hdiv, rot (and ssh if necessary) 
     277            ! set the dynamics: U,V, hdiv (and ssh if necessary) 
    283278            ! ---------------- 
    284279            ! Start EEL5 configuration with barotropic geostrophic velocities  
     
    316311            ENDIF 
    317312            ! 
    318             CALL div_cur( nit000 )                  ! horizontal divergence and relative vorticity (curl) 
     313            CALL div_hor( nit000 )                  ! horizontal divergence and relative vorticity (curl) 
    319314            ! N.B. the vertical velocity will be computed from the horizontal divergence field 
    320315            ! in istate by a call to wzv routine 
     
    369364      !! 
    370365      !! ** Method  : - set temprature field 
    371       !!              - set salinity field 
     366      !!              - set salinity   field 
    372367      !!---------------------------------------------------------------------- 
    373368      INTEGER :: ji, jj, jk  ! dummy loop indices 
     
    443438   END SUBROUTINE istate_gyre 
    444439 
     440 
    445441   SUBROUTINE istate_uvg 
    446442      !!---------------------------------------------------------------------- 
     
    455451      !!---------------------------------------------------------------------- 
    456452      USE dynspg          ! surface pressure gradient             (dyn_spg routine) 
    457       USE divcur          ! hor. divergence & rel. vorticity      (div_cur routine) 
     453      USE divhor          ! hor. divergence                       (div_hor routine) 
    458454      USE lbclnk          ! ocean lateral boundary condition (or mpp link) 
    459  
     455      ! 
    460456      INTEGER ::   ji, jj, jk        ! dummy loop indices 
    461457      INTEGER ::   indic             ! ??? 
     
    553549      un(:,:,:) = ub(:,:,:) 
    554550      vn(:,:,:) = vb(:,:,:) 
    555         
    556       ! Compute the divergence and curl 
    557  
    558       CALL div_cur( nit000 )            ! now horizontal divergence and curl 
    559  
    560       hdivb(:,:,:) = hdivn(:,:,:)       ! set the before to the now value 
    561       rotb (:,:,:) = rotn (:,:,:)       ! set the before to the now value 
     551      ! 
     552!!gm  Check  here call to div_hor should not be necessary 
     553!!gm         div_hor call runoffs  not sure they are defined at that level 
     554      CALL div_hor( nit000 )            ! now horizontal divergence 
    562555      ! 
    563556      CALL wrk_dealloc( jpi, jpj, jpk, zprn) 
Note: See TracChangeset for help on using the changeset viewer.