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 5870 for branches/2015/dev_r5803_NOC_WAD/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

Ignore:
Timestamp:
2015-11-09T18:33:54+01:00 (9 years ago)
Author:
acc
Message:

Branch 2015/dev_r5803_NOC_WAD. Merge in trunk changes from 5803 to 5869 in preparation for merge. Also tidied and reorganised some wetting and drying code. Renamed wadlmt.F90 to wetdry.F90. Wetting drying code changes restricted to domzgr.F90, domvvl.F90 nemogcm.F90 sshwzv.F90, dynspg_ts.F90, wetdry.F90 and dynhpg.F90. Code passes full SETTE tests with ln_wd=.false.. Still awaiting test case for checking with ln_wd=.false.

Location:
branches/2015/dev_r5803_NOC_WAD/NEMOGCM/NEMO/OPA_SRC/BDY
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • branches/2015/dev_r5803_NOC_WAD/NEMOGCM/NEMO/OPA_SRC/BDY/bdy_oce.F90

    r4699 r5870  
    88   !!            3.3  !  2010-09  (D. Storkey) add ice boundary conditions 
    99   !!            3.4  !  2011     (D. Storkey) rewrite in preparation for OBC-BDY merge 
    10    !!            3.6  !  2012-01  (C. Rousset) add ice boundary conditions for lim3 
     10   !!            3.6  !  2014-01  (C. Rousset) add ice boundary conditions for lim3 
    1111   !!---------------------------------------------------------------------- 
    1212#if defined key_bdy  
     
    2222 
    2323   TYPE, PUBLIC ::   OBC_INDEX    !: Indices and weights which define the open boundary 
    24       INTEGER,          DIMENSION(jpbgrd) ::  nblen 
    25       INTEGER,          DIMENSION(jpbgrd) ::  nblenrim 
    26       INTEGER, POINTER, DIMENSION(:,:)   ::  nbi 
    27       INTEGER, POINTER, DIMENSION(:,:)   ::  nbj 
    28       INTEGER, POINTER, DIMENSION(:,:)   ::  nbr 
    29       INTEGER, POINTER, DIMENSION(:,:)   ::  nbmap 
    30       REAL(wp)   , POINTER, DIMENSION(:,:)   ::  nbw 
    31       REAL(wp)   , POINTER, DIMENSION(:,:)   ::  nbd 
    32       REAL(wp)   , POINTER, DIMENSION(:,:)   ::  nbdout 
    33       REAL(wp)   , POINTER, DIMENSION(:,:)   ::  flagu 
    34       REAL(wp)   , POINTER, DIMENSION(:,:)   ::  flagv 
     24      INTEGER ,          DIMENSION(jpbgrd) ::  nblen 
     25      INTEGER ,          DIMENSION(jpbgrd) ::  nblenrim 
     26      INTEGER , POINTER, DIMENSION(:,:)    ::  nbi 
     27      INTEGER , POINTER, DIMENSION(:,:)    ::  nbj 
     28      INTEGER , POINTER, DIMENSION(:,:)    ::  nbr 
     29      INTEGER , POINTER, DIMENSION(:,:)    ::  nbmap 
     30      REAL(wp), POINTER, DIMENSION(:,:)    ::  nbw 
     31      REAL(wp), POINTER, DIMENSION(:,:)    ::  nbd 
     32      REAL(wp), POINTER, DIMENSION(:,:)    ::  nbdout 
     33      REAL(wp), POINTER, DIMENSION(:,:)    ::  flagu 
     34      REAL(wp), POINTER, DIMENSION(:,:)    ::  flagv 
    3535   END TYPE OBC_INDEX 
    3636 
     
    4141 
    4242   TYPE, PUBLIC ::   OBC_DATA     !: Storage for external data 
    43       INTEGER,       DIMENSION(2)     ::  nread 
    44       LOGICAL                         ::  ll_ssh 
    45       LOGICAL                         ::  ll_u2d 
    46       LOGICAL                         ::  ll_v2d 
    47       LOGICAL                         ::  ll_u3d 
    48       LOGICAL                         ::  ll_v3d 
    49       LOGICAL                         ::  ll_tem 
    50       LOGICAL                         ::  ll_sal 
    51       REAL(wp), POINTER, DIMENSION(:)     ::  ssh 
    52       REAL(wp), POINTER, DIMENSION(:)     ::  u2d 
    53       REAL(wp), POINTER, DIMENSION(:)     ::  v2d 
    54       REAL(wp), POINTER, DIMENSION(:,:)   ::  u3d 
    55       REAL(wp), POINTER, DIMENSION(:,:)   ::  v3d 
    56       REAL(wp), POINTER, DIMENSION(:,:)   ::  tem 
    57       REAL(wp), POINTER, DIMENSION(:,:)   ::  sal 
     43      INTEGER          , DIMENSION(2)   ::  nread 
     44      LOGICAL                           ::  ll_ssh 
     45      LOGICAL                           ::  ll_u2d 
     46      LOGICAL                           ::  ll_v2d 
     47      LOGICAL                           ::  ll_u3d 
     48      LOGICAL                           ::  ll_v3d 
     49      LOGICAL                           ::  ll_tem 
     50      LOGICAL                           ::  ll_sal 
     51      REAL(wp), POINTER, DIMENSION(:)   ::  ssh 
     52      REAL(wp), POINTER, DIMENSION(:)   ::  u2d 
     53      REAL(wp), POINTER, DIMENSION(:)   ::  v2d 
     54      REAL(wp), POINTER, DIMENSION(:,:) ::  u3d 
     55      REAL(wp), POINTER, DIMENSION(:,:) ::  v3d 
     56      REAL(wp), POINTER, DIMENSION(:,:) ::  tem 
     57      REAL(wp), POINTER, DIMENSION(:,:) ::  sal 
    5858#if defined key_lim2 
    59       LOGICAL                         ::  ll_frld 
    60       LOGICAL                         ::  ll_hicif 
    61       LOGICAL                         ::  ll_hsnif 
    62       REAL(wp), POINTER, DIMENSION(:)     ::  frld 
    63       REAL(wp), POINTER, DIMENSION(:)     ::  hicif 
    64       REAL(wp), POINTER, DIMENSION(:)     ::  hsnif 
     59      LOGICAL                           ::   ll_frld 
     60      LOGICAL                           ::   ll_hicif 
     61      LOGICAL                           ::   ll_hsnif 
     62      REAL(wp), POINTER, DIMENSION(:)   ::   frld 
     63      REAL(wp), POINTER, DIMENSION(:)   ::   hicif 
     64      REAL(wp), POINTER, DIMENSION(:)   ::   hsnif 
    6565#elif defined key_lim3 
    66       LOGICAL                         ::  ll_a_i 
    67       LOGICAL                         ::  ll_ht_i 
    68       LOGICAL                         ::  ll_ht_s 
    69       REAL, POINTER, DIMENSION(:,:)   ::  a_i   !: now ice leads fraction climatology 
    70       REAL, POINTER, DIMENSION(:,:)   ::  ht_i  !: Now ice  thickness climatology 
    71       REAL, POINTER, DIMENSION(:,:)   ::  ht_s  !: now snow thickness 
     66      LOGICAL                           ::   ll_a_i 
     67      LOGICAL                           ::   ll_ht_i 
     68      LOGICAL                           ::   ll_ht_s 
     69      REAL(wp), POINTER, DIMENSION(:,:) ::   a_i    !: now ice leads fraction climatology 
     70      REAL(wp), POINTER, DIMENSION(:,:) ::   ht_i   !: Now ice  thickness climatology 
     71      REAL(wp), POINTER, DIMENSION(:,:) ::   ht_s   !: now snow thickness 
    7272#endif 
    7373   END TYPE OBC_DATA 
     
    9999   INTEGER, DIMENSION(jp_bdy)           ::   nn_tra_dta     !: = 0 use the initial state as bdy dta ;  
    100100                                                            !: = 1 read it in a NetCDF file 
    101    LOGICAL, DIMENSION(jp_bdy) ::   ln_tra_dmp               !: =T Tracer damping 
    102    LOGICAL, DIMENSION(jp_bdy) ::   ln_dyn3d_dmp             !: =T Baroclinic velocity damping 
    103    REAL(wp),    DIMENSION(jp_bdy) ::   rn_time_dmp              !: Damping time scale in days 
    104    REAL(wp),    DIMENSION(jp_bdy) ::   rn_time_dmp_out          !: Damping time scale in days at radiation outflow points 
     101   LOGICAL , DIMENSION(jp_bdy) ::   ln_tra_dmp              !: =T Tracer damping 
     102   LOGICAL , DIMENSION(jp_bdy) ::   ln_dyn3d_dmp            !: =T Baroclinic velocity damping 
     103   REAL(wp), DIMENSION(jp_bdy) ::   rn_time_dmp             !: Damping time scale in days 
     104   REAL(wp), DIMENSION(jp_bdy) ::   rn_time_dmp_out         !: Damping time scale in days at radiation outflow points 
    105105 
    106106   CHARACTER(len=20), DIMENSION(jp_bdy) ::   cn_ice_lim       ! Choice of boundary condition for sea ice variables  
    107    INTEGER, DIMENSION(jp_bdy)           ::   nn_ice_lim_dta   !: = 0 use the initial state as bdy dta ;  
     107   INTEGER , DIMENSION(jp_bdy)          ::   nn_ice_lim_dta   !: = 0 use the initial state as bdy dta ;  
    108108                                                              !: = 1 read it in a NetCDF file 
    109    REAL(wp),    DIMENSION(jp_bdy) ::   rn_ice_tem             !: choice of the temperature of incoming sea ice 
    110    REAL(wp),    DIMENSION(jp_bdy) ::   rn_ice_sal             !: choice of the salinity    of incoming sea ice 
    111    REAL(wp),    DIMENSION(jp_bdy) ::   rn_ice_age             !: choice of the age         of incoming sea ice 
     109   REAL(wp), DIMENSION(jp_bdy) ::   rn_ice_tem              !: choice of the temperature of incoming sea ice 
     110   REAL(wp), DIMENSION(jp_bdy) ::   rn_ice_sal              !: choice of the salinity    of incoming sea ice 
     111   REAL(wp), DIMENSION(jp_bdy) ::   rn_ice_age              !: choice of the age         of incoming sea ice 
    112112   ! 
    113113    
  • branches/2015/dev_r5803_NOC_WAD/NEMOGCM/NEMO/OPA_SRC/BDY/bdyice_lim.F90

    r5656 r5870  
    5959      !! 
    6060      !!---------------------------------------------------------------------- 
    61       INTEGER, INTENT( in ) :: kt     ! Main time step counter 
    62       INTEGER               :: ib_bdy ! Loop index 
    63  
     61      INTEGER, INTENT( in ) ::   kt   ! Main time step counter 
     62      ! 
     63      INTEGER ::   ib_bdy   ! Loop index 
     64      !!---------------------------------------------------------------------- 
     65      ! 
    6466#if defined key_lim3 
    6567      CALL lim_var_glo2eqv 
    6668#endif 
    67  
     69      ! 
    6870      DO ib_bdy=1, nb_bdy 
    69  
     71         ! 
    7072         SELECT CASE( cn_ice_lim(ib_bdy) ) 
    7173         CASE('none') 
     
    7678            CALL ctl_stop( 'bdy_ice_lim : unrecognised option for open boundaries for ice fields' ) 
    7779         END SELECT 
    78  
     80         ! 
    7981      END DO 
    80  
     82      ! 
    8183#if defined key_lim3 
    8284      CALL lim_var_zapsmall 
    8385      CALL lim_var_agg(1) 
    8486#endif 
    85  
     87      ! 
    8688   END SUBROUTINE bdy_ice_lim 
     89 
    8790 
    8891   SUBROUTINE bdy_ice_frs( idx, dta, kt, ib_bdy ) 
     
    9699      !!             dimensional baroclinic ocean model with realistic topography. Tellus, 365-382. 
    97100      !!------------------------------------------------------------------------------ 
    98       TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices 
    99       TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data 
    100       INTEGER,         INTENT(in) ::   kt   ! main time-step counter 
     101      TYPE(OBC_INDEX), INTENT(in) ::   idx     ! OBC indices 
     102      TYPE(OBC_DATA),  INTENT(in) ::   dta     ! OBC external data 
     103      INTEGER,         INTENT(in) ::   kt      ! main time-step counter 
    101104      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index 
    102  
     105      ! 
    103106      INTEGER  ::   jpbound            ! 0 = incoming ice 
    104                                        ! 1 = outgoing ice 
     107      !                                ! 1 = outgoing ice 
    105108      INTEGER  ::   jb, jk, jgrd, jl   ! dummy loop indices 
    106109      INTEGER  ::   ji, jj, ii, ij     ! local scalar 
     
    111114     USE ice_2, vt_i => hicm 
    112115#endif 
    113  
    114       !!------------------------------------------------------------------------------ 
    115       ! 
    116       IF( nn_timing == 1 ) CALL timing_start('bdy_ice_frs') 
     116      !!------------------------------------------------------------------------------ 
     117      ! 
     118      IF( nn_timing == 1 )   CALL timing_start('bdy_ice_frs') 
    117119      ! 
    118120      jgrd = 1      ! Everything is at T-points here 
     
    181183            ! condition on ice thickness depends on the ice velocity 
    182184            ! if velocity is outward (strictly), then ice thickness, volume... must be equal to adjacent values 
    183             jpbound = 0; ii = ji; ij = jj; 
    184  
     185            jpbound = 0   ;   ii = ji   ;   ij = jj 
     186            ! 
    185187            IF( u_ice(ji+1,jj  ) < 0. .AND. umask(ji-1,jj  ,1) == 0. ) jpbound = 1; ii = ji+1; ij = jj 
    186188            IF( u_ice(ji-1,jj  ) > 0. .AND. umask(ji+1,jj  ,1) == 0. ) jpbound = 1; ii = ji-1; ij = jj 
    187189            IF( v_ice(ji  ,jj+1) < 0. .AND. vmask(ji  ,jj-1,1) == 0. ) jpbound = 1; ii = ji  ; ij = jj+1 
    188190            IF( v_ice(ji  ,jj-1) > 0. .AND. vmask(ji  ,jj+1,1) == 0. ) jpbound = 1; ii = ji  ; ij = jj-1 
    189  
     191            ! 
    190192            IF( nn_ice_lim_dta(ib_bdy) == 0 ) jpbound = 0; ii = ji; ij = jj   ! case ice boundaries = initial conditions 
    191                                                                               !      do not make state variables dependent on velocity 
    192                 
    193  
     193            !                                                                 !      do not make state variables dependent on velocity 
     194            ! 
    194195            rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ii,ij) - 0.01 ) ) ! 0 if no ice 
    195  
     196            ! 
    196197            ! concentration and thickness 
    197198            a_i (ji,jj,jl) = a_i (ii,ij,jl) * rswitch 
    198199            ht_i(ji,jj,jl) = ht_i(ii,ij,jl) * rswitch 
    199200            ht_s(ji,jj,jl) = ht_s(ii,ij,jl) * rswitch 
    200  
     201            ! 
    201202            ! Ice and snow volumes 
    202203            v_i(ji,jj,jl) = ht_i(ji,jj,jl) * a_i(ji,jj,jl) 
    203204            v_s(ji,jj,jl) = ht_s(ji,jj,jl) * a_i(ji,jj,jl) 
    204  
     205            ! 
    205206            SELECT CASE( jpbound ) 
    206  
    207             CASE( 0 ) ! velocity is inward 
    208  
     207            ! 
     208            CASE( 0 )   ! velocity is inward 
     209               ! 
    209210               ! Ice salinity, age, temperature 
    210211               sm_i(ji,jj,jl)   = rswitch * rn_ice_sal(ib_bdy)  + ( 1.0 - rswitch ) * rn_simin 
     
    218219                  s_i(ji,jj,jk,jl) = rswitch * rn_ice_sal(ib_bdy) + ( 1.0 - rswitch ) * rn_simin 
    219220               END DO 
    220                 
    221             CASE( 1 ) ! velocity is outward 
    222   
     221               ! 
     222            CASE( 1 )   ! velocity is outward 
     223               ! 
    223224               ! Ice salinity, age, temperature 
    224225               sm_i(ji,jj,jl)   = rswitch * sm_i(ii,ij,jl)  + ( 1.0 - rswitch ) * rn_simin 
     
    232233                  s_i(ji,jj,jk,jl) = rswitch * s_i(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rn_simin 
    233234               END DO 
    234  
     235               ! 
    235236            END SELECT 
    236  
    237             ! if salinity is constant, then overwrite rn_ice_sal 
    238             IF( nn_icesal == 1 ) THEN 
    239                sm_i(ji,jj,jl)   = rn_icesal 
     237            ! 
     238            IF( nn_icesal == 1 ) THEN     ! constant salinity : overwrite rn_ice_sal 
     239               sm_i(ji,jj  ,jl) = rn_icesal 
    240240               s_i (ji,jj,:,jl) = rn_icesal 
    241241            ENDIF 
    242  
     242            ! 
    243243            ! contents 
    244244            smv_i(ji,jj,jl)  = MIN( sm_i(ji,jj,jl) , sss_m(ji,jj) ) * v_i(ji,jj,jl) 
     
    259259               e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) * a_i(ji,jj,jl) * ht_i(ji,jj,jl) * r1_nlay_i 
    260260            END DO 
    261  
     261            ! 
    262262         END DO 
    263   
     263         ! 
    264264         CALL lbc_bdy_lnk(  a_i(:,:,jl), 'T', 1., ib_bdy ) 
    265265         CALL lbc_bdy_lnk( ht_i(:,:,jl), 'T', 1., ib_bdy ) 
     
    267267         CALL lbc_bdy_lnk(  v_i(:,:,jl), 'T', 1., ib_bdy ) 
    268268         CALL lbc_bdy_lnk(  v_s(:,:,jl), 'T', 1., ib_bdy ) 
    269   
     269         ! 
    270270         CALL lbc_bdy_lnk( smv_i(:,:,jl), 'T', 1., ib_bdy ) 
    271271         CALL lbc_bdy_lnk(  sm_i(:,:,jl), 'T', 1., ib_bdy ) 
     
    280280            CALL lbc_bdy_lnk(e_i(:,:,jk,jl), 'T', 1., ib_bdy ) 
    281281         END DO 
    282  
     282         ! 
    283283      END DO !jl 
    284      
     284      ! 
    285285#endif 
    286286      !       
    287       IF( nn_timing == 1 ) CALL timing_stop('bdy_ice_frs') 
     287      IF( nn_timing == 1 )   CALL timing_stop('bdy_ice_frs') 
    288288      ! 
    289289   END SUBROUTINE bdy_ice_frs 
     
    300300      !! 2013-06 : C. Rousset 
    301301      !!------------------------------------------------------------------------------ 
    302       !! 
    303302      CHARACTER(len=1), INTENT(in)  ::   cd_type   ! nature of velocity grid-points 
     303      ! 
    304304      INTEGER  ::   jb, jgrd           ! dummy loop indices 
    305305      INTEGER  ::   ji, jj             ! local scalar 
    306306      INTEGER  ::   ib_bdy             ! Loop index 
    307307      REAL(wp) ::   zmsk1, zmsk2, zflag 
    308      !!------------------------------------------------------------------------------ 
     308      !!------------------------------------------------------------------------------ 
    309309      ! 
    310310      IF( nn_timing == 1 ) CALL timing_start('bdy_ice_lim_dyn') 
     
    313313         ! 
    314314         SELECT CASE( cn_ice_lim(ib_bdy) ) 
    315  
     315         ! 
    316316         CASE('none') 
    317  
    318317            CYCLE 
    319              
     318            ! 
    320319         CASE('frs') 
    321              
     320            ! 
    322321            IF( nn_ice_lim_dta(ib_bdy) == 0 ) CYCLE            ! case ice boundaries = initial conditions  
    323                                                                !      do not change ice velocity (it is only computed by rheology) 
    324   
     322            !                                                  !      do not change ice velocity (it is only computed by rheology) 
    325323            SELECT CASE ( cd_type ) 
    326                 
    327             CASE ( 'U' ) 
    328                 
     324            !      
     325            CASE ( 'U' )   
    329326               jgrd = 2      ! u velocity 
    330327               DO jb = 1, idx_bdy(ib_bdy)%nblen(jgrd) 
     
    332329                  jj    = idx_bdy(ib_bdy)%nbj(jb,jgrd) 
    333330                  zflag = idx_bdy(ib_bdy)%flagu(jb,jgrd) 
    334                    
     331                  ! 
    335332                  IF ( ABS( zflag ) == 1. ) THEN  ! eastern and western boundaries 
    336333                     ! one of the two zmsk is always 0 (because of zflag) 
    337334                     zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji+1,jj) ) ) ! 0 if no ice 
    338335                     zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji-1,jj) ) ) ! 0 if no ice 
    339                       
     336                      
    340337                     ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then u_ice = u_oce) 
    341338                     u_ice (ji,jj) = u_ice(ji+1,jj) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + & 
     
    349346                  rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ji,jj) - 0.01_wp ) ) ! 0 if no ice 
    350347                  u_ice(ji,jj) = rswitch * u_ice(ji,jj) 
    351                    
    352                ENDDO 
    353                 
     348                  ! 
     349               END DO 
    354350               CALL lbc_bdy_lnk( u_ice(:,:), 'U', -1., ib_bdy ) 
    355                 
     351               ! 
    356352            CASE ( 'V' ) 
    357                 
    358353               jgrd = 3      ! v velocity 
    359354               DO jb = 1, idx_bdy(ib_bdy)%nblen(jgrd) 
     
    361356                  jj    = idx_bdy(ib_bdy)%nbj(jb,jgrd) 
    362357                  zflag = idx_bdy(ib_bdy)%flagv(jb,jgrd) 
    363                    
     358                  ! 
    364359                  IF ( ABS( zflag ) == 1. ) THEN  ! northern and southern boundaries 
    365360                     ! one of the two zmsk is always 0 (because of zflag) 
    366361                     zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj+1) ) ) ! 0 if no ice 
    367362                     zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj-1) ) ) ! 0 if no ice 
    368                       
     363                      
    369364                     ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then u_ice = u_oce) 
    370365                     v_ice (ji,jj) = v_ice(ji,jj+1) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + & 
     
    378373                  rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ji,jj) - 0.01 ) ) ! 0 if no ice 
    379374                  v_ice(ji,jj) = rswitch * v_ice(ji,jj) 
    380                    
    381                ENDDO 
    382                 
     375                  ! 
     376               END DO 
    383377               CALL lbc_bdy_lnk( v_ice(:,:), 'V', -1., ib_bdy ) 
    384                    
     378               ! 
    385379            END SELECT 
    386              
     380            ! 
    387381         CASE DEFAULT 
    388382            CALL ctl_stop( 'bdy_ice_lim_dyn : unrecognised option for open boundaries for ice fields' ) 
    389383         END SELECT 
    390           
    391       ENDDO 
    392  
     384         ! 
     385      END DO 
     386      ! 
    393387      IF( nn_timing == 1 ) CALL timing_stop('bdy_ice_lim_dyn') 
    394        
     388      ! 
    395389    END SUBROUTINE bdy_ice_lim_dyn 
    396390 
  • branches/2015/dev_r5803_NOC_WAD/NEMOGCM/NEMO/OPA_SRC/BDY/bdyini.F90

    r5656 r5870  
    7676      INTEGER  ::   ib_bdy, ii, ij, ik, igrd, ib, ir, iseg ! dummy loop indices 
    7777      INTEGER  ::   icount, icountr, ibr_max, ilen1, ibm1  ! local integers 
    78       INTEGER  ::   iwe, ies, iso, ino, inum, id_dummy         !   -       - 
     78      INTEGER  ::   iwe, ies, iso, ino, inum, id_dummy     !   -       - 
    7979      INTEGER  ::   igrd_start, igrd_end, jpbdta           !   -       - 
    8080      INTEGER  ::   jpbdtau, jpbdtas                       !   -       - 
  • branches/2015/dev_r5803_NOC_WAD/NEMOGCM/NEMO/OPA_SRC/BDY/bdyvol.F90

    r5643 r5870  
    1515   !!   'key_dynspg_flt'                              filtered free surface 
    1616   !!---------------------------------------------------------------------- 
    17    USE timing          ! Timing 
    1817   USE oce             ! ocean dynamics and tracers  
    19    USE sbcisf          ! ice shelf 
     18   USE bdy_oce         ! ocean open boundary conditions 
     19   USE sbc_oce         ! ocean surface boundary conditions 
    2020   USE dom_oce         ! ocean space and time domain  
    2121   USE phycst          ! physical constants 
    22    USE bdy_oce         ! ocean open boundary conditions 
     22   USE sbcisf          ! ice shelf 
     23   ! 
     24   USE in_out_manager  ! I/O manager 
    2325   USE lib_mpp         ! for mppsum 
    24    USE in_out_manager  ! I/O manager 
    25    USE sbc_oce         ! ocean surface boundary conditions 
     26   USE timing          ! Timing 
     27   USE lib_fortran     ! Fortran routines library 
    2628 
    2729   IMPLICIT NONE 
     
    3335#  include "domzgr_substitute.h90" 
    3436   !!---------------------------------------------------------------------- 
    35    !! NEMO/OPA 3.3 , NEMO Consortium (2010) 
     37   !! NEMO/OPA 3.6 , NEMO Consortium (2014) 
    3638   !! $Id$  
    3739   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) 
     
    7880      TYPE(OBC_INDEX), POINTER :: idx 
    7981      !!----------------------------------------------------------------------------- 
    80  
    81       IF( nn_timing == 1 ) CALL timing_start('bdy_vol') 
    82  
     82      ! 
     83      IF( nn_timing == 1 )   CALL timing_start('bdy_vol') 
     84      ! 
    8385      IF( ln_vol ) THEN 
    84  
     86      ! 
    8587      IF( kt == nit000 ) THEN  
    8688         IF(lwp) WRITE(numout,*) 
     
    9193      ! Calculate the cumulate surface Flux z_cflxemp (m3/s) over all the domain 
    9294      ! ----------------------------------------------------------------------- 
    93       z_cflxemp = SUM ( ( emp(:,:)-rnf(:,:)+fwfisf(:,:) ) * bdytmask(:,:) * e1t(:,:) * e2t(:,:) ) / rau0 
     95!!gm replace these lines : 
     96      z_cflxemp = SUM ( ( emp(:,:)-rnf(:,:)+fwfisf(:,:) ) * bdytmask(:,:) * e1e2t(:,:) ) / rau0 
    9497      IF( lk_mpp )   CALL mpp_sum( z_cflxemp )     ! sum over the global domain 
     98!!gm   by : 
     99!!gm      z_cflxemp = glob_sum(  ( emp(:,:)-rnf(:,:)+fwfisf(:,:) ) * bdytmask(:,:) * e1e2t(:,:)  ) / rau0 
     100!!gm 
    95101 
    96102      ! Transport through the unstructured open boundary 
    97103      ! ------------------------------------------------ 
    98       zubtpecor = 0.e0 
     104      zubtpecor = 0._wp 
    99105      DO ib_bdy = 1, nb_bdy 
    100106         idx => idx_bdy(ib_bdy) 
    101  
     107         ! 
    102108         jgrd = 2                               ! cumulate u component contribution first  
    103109         DO jb = 1, idx%nblenrim(jgrd) 
     
    116122            END DO 
    117123         END DO 
    118  
     124         ! 
    119125      END DO 
    120126      IF( lk_mpp )   CALL mpp_sum( zubtpecor )   ! sum over the global domain 
     
    123129      ! ------------------------------ 
    124130      IF( nn_volctl==1 ) THEN   ;   zubtpecor = ( zubtpecor - z_cflxemp) / bdysurftot  
    125       ELSE                   ;   zubtpecor =   zubtpecor             / bdysurftot 
     131      ELSE                      ;   zubtpecor =   zubtpecor             / bdysurftot 
    126132      END IF 
    127133 
    128134      ! Correction of the total velocity on the unstructured boundary to respect the mass flux conservation 
    129135      ! ------------------------------------------------------------- 
    130       ztranst = 0.e0 
     136      ztranst = 0._wp 
    131137      DO ib_bdy = 1, nb_bdy 
    132138         idx => idx_bdy(ib_bdy) 
    133  
     139         ! 
    134140         jgrd = 2                               ! correct u component 
    135141         DO jb = 1, idx%nblenrim(jgrd) 
     
    150156            END DO 
    151157         END DO 
    152  
     158         ! 
    153159      END DO 
    154160      IF( lk_mpp )   CALL mpp_sum( ztranst )   ! sum over the global domain 
     
    169175      ! 
    170176      END IF ! ln_vol 
    171  
     177      ! 
    172178   END SUBROUTINE bdy_vol 
    173179 
Note: See TracChangeset for help on using the changeset viewer.