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 13024 for utils/tools_dev_r12970_AGRIF_CMEMS/DOMAINcfg/src/agrif_user.F90 – NEMO

Ignore:
Timestamp:
2020-06-03T16:26:23+02:00 (4 years ago)
Author:
rblod
Message:

First version of new nesting tools merged with domaincfg, see ticket #2129

File:
1 edited

Legend:

Unmodified
Added
Removed
  • utils/tools_dev_r12970_AGRIF_CMEMS/DOMAINcfg/src/agrif_user.F90

    r12414 r13024  
    11#if defined key_agrif 
    2 subroutine agrif_initworkspace() 
     2   SUBROUTINE agrif_initworkspace() 
    33      !!---------------------------------------------------------------------- 
    44      !!                 *** ROUTINE Agrif_InitWorkspace *** 
    55      !!---------------------------------------------------------------------- 
    6 end subroutine agrif_initworkspace 
    7 SUBROUTINE Agrif_InitValues 
     6   END SUBROUTINE agrif_initworkspace 
     7 
     8   SUBROUTINE Agrif_InitValues 
    89      !!---------------------------------------------------------------------- 
    910      !!                 *** ROUTINE Agrif_InitValues *** 
     
    1112      !! ** Purpose :: Declaration of variables to be interpolated 
    1213      !!---------------------------------------------------------------------- 
    13    USE Agrif_Util 
    14    USE dom_oce 
    15    USE nemogcm 
    16    USE domain 
    17    !! 
    18    IMPLICIT NONE 
    19     
    20   
    21    CALL nemo_init       !* Initializations of each fine grid 
    22  
    23    CALL dom_nam 
    24    CALL cfg_write         ! create the configuration file 
    25     
    26 END SUBROUTINE Agrif_InitValues 
    27  
    28 SUBROUTINE Agrif_InitValues_cont 
    29  
    30 use dom_oce 
    31     integer :: irafx, irafy 
    32     logical :: ln_perio 
    33     integer nx,ny 
    34  
    35 irafx = agrif_irhox() 
    36 irafy = agrif_irhoy() 
    37  
    38 nx=nlci ; ny=nlcj 
     14      USE Agrif_Util 
     15      USE dom_oce 
     16      USE nemogcm 
     17      USE domain 
     18      !! 
     19      IMPLICIT NONE 
     20 
     21      ! No temporal refinement 
     22      CALL Agrif_Set_coeffreft(1) 
     23 
     24      CALL nemo_init       !* Initializations of each fine grid 
     25 
     26      CALL dom_nam 
     27 
     28   END SUBROUTINE Agrif_InitValues 
     29 
     30   SUBROUTINE Agrif_InitValues_cont 
     31      !!---------------------------------------------------------------------- 
     32      !!                 *** ROUTINE Agrif_InitValues_cont *** 
     33      !! 
     34      !! ** Purpose :: Initialisation of variables to be interpolated 
     35      !!---------------------------------------------------------------------- 
     36      USE dom_oce 
     37      USE lbclnk 
     38      !! 
     39      IMPLICIT NONE 
     40      ! 
     41      INTEGER :: nx, ny 
     42      INTEGER :: irafx, irafy 
     43      LOGICAL :: ln_perio 
     44      ! 
     45      irafx = agrif_irhox() 
     46      irafy = agrif_irhoy() 
     47 
     48      nx = nlci ; ny = nlcj 
    3949 
    4050   !       IF(jperio /=1 .AND. jperio/=4 .AND. jperio/=6 ) THEN 
     
    4454   !          nbghostcellsfine_tot_y= nbghostcellsfine_y +1 
    4555   !       ELSE 
    46    !         nx = (nbcellsx)+2*nbghostcellsfine_x  
     56   !         nx = (nbcellsx)+2*nbghostcellsfine_x 
    4757   !         ny = (nbcellsy)+2*nbghostcellsfine+2 
    4858   !         nbghostcellsfine_tot_x= 1 
     
    5363   !       nx = nbcellsx+irafx 
    5464   !       ny = nbcellsy+irafy 
    55         
    56   WRITE(*,*) ' ' 
    57   WRITE(*,*)'Size of the High resolution grid: ',nx,' x ',ny 
    58   WRITE(*,*) ' ' 
    59         
    60        call agrif_init_lonlat() 
    61        ln_perio=.FALSE.  
    62        if( jperio ==1 .OR. jperio==2 .OR. jperio==4) ln_perio=.TRUE. 
    63  
    64        where (glamt < -180) glamt = glamt +360. 
    65        if (ln_perio) then 
    66          glamt(1,:)=glamt(nx-1,:) 
    67          glamt(nx,:)=glamt(2,:) 
    68        endif 
    69   
    70        where (glamu < -180) glamu = glamu +360. 
    71        if (ln_perio) then 
    72          glamu(1,:)=glamu(nx-1,:) 
    73          glamu(nx,:)=glamu(2,:) 
    74        endif 
    75  
    76       where (glamv < -180) glamv = glamv +360. 
    77        if (ln_perio) then 
    78          glamv(1,:)=glamv(nx-1,:) 
    79          glamv(nx,:)=glamv(2,:) 
    80        endif 
    81  
    82       where (glamf < -180) glamf = glamf +360. 
    83        if (ln_perio) then 
    84          glamf(1,:)=glamf(nx-1,:) 
    85          glamf(nx,:)=glamf(2,:) 
    86        endif 
    87  
    88        call agrif_init_scales() 
    89  
    90          
    91 END SUBROUTINE Agrif_InitValues_cont   
    92  
    93  
    94 subroutine agrif_declare_var() 
    95 use par_oce 
    96 use dom_oce 
    97 use agrif_profiles 
    98 use agrif_parameters 
    99  
    100    IMPLICIT NONE 
    101     
    102 integer :: ind1, ind2, ind3 
    103 integer nx,ny 
    104 integer nbghostcellsfine_tot_x, nbghostcellsfine_tot_y 
    105 INTEGER :: irafx 
    106 !!---------------------------------------------------------------------- 
    107  
    108    ! 1. Declaration of the type of variable which have to be interpolated 
    109    !--------------------------------------------------------------------- 
    110  nx=nlci ; ny=nlcj 
    111  
    112 !ind2 = nbghostcellsfine_tot_x + 1 
    113 !ind3 = nbghostcellsfine_tot_y + 1 
    114 ind2 = 2 + nbghostcells 
    115 ind3 = ind2 
    116 nbghostcellsfine_tot_x=nbghostcells+1 
    117 nbghostcellsfine_tot_y=nbghostcells+1 
    118  
    119 irafx = Agrif_irhox() 
    120  
    121 CALL agrif_nemo_init  ! specific namelist part if needed 
    122  
    123 CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),glamt_id) 
    124 CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),glamu_id) 
    125 CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),glamv_id) 
    126 CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),glamf_id) 
    127  
    128 CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),gphit_id) 
    129 CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),gphiu_id) 
    130 CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),gphiv_id) 
    131 CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),gphif_id) 
    132  
    133 ! Horizontal scale factors 
    134  
    135 CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),e1t_id) 
    136 CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),e1u_id) 
    137 CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),e1v_id) 
    138 CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),e1f_id) 
    139  
    140 CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),e2t_id) 
    141 CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),e2u_id) 
    142 CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),e2v_id) 
    143 CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),e2f_id) 
    144  
    145 ! Bathymetry 
    146  
    147 CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),bathy_id) 
    148  
    149 ! Vertical scale factors 
    150 CALL agrif_declare_variable((/2,2,0/),(/ind2,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/nx,ny,jpk/),e3t_id) 
    151 CALL agrif_declare_variable((/2,2,0/),(/ind2,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/nx,ny,jpk/),e3t_copy_id) 
    152 CALL agrif_declare_variable((/2,2,0/),(/ind2,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/nx,ny,jpk+1/),e3t_connect_id) 
    153  
    154 CALL agrif_declare_variable((/1,2,0/),(/ind2-1,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/nx,ny,jpk/),e3u_id) 
    155 CALL agrif_declare_variable((/2,1,0/),(/ind2,ind3-1,0/),(/'x','y','N'/),(/1,1,1/),(/nx,ny,jpk/),e3v_id) 
    156  
    157 ! Bottom level 
    158  
    159 CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),bottom_level_id) 
    160  
    161 CALL Agrif_Set_bcinterp(glamt_id,interp=AGRIF_linear) 
    162 CALL Agrif_Set_interp(glamt_id,interp=AGRIF_linear) 
    163 CALL Agrif_Set_bc( glamt_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    164  
    165 CALL Agrif_Set_bcinterp(glamu_id,interp=AGRIF_linear) 
    166 CALL Agrif_Set_interp(glamu_id,interp=AGRIF_linear) 
    167 CALL Agrif_Set_bc( glamu_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    168  
    169 CALL Agrif_Set_bcinterp(glamv_id,interp=AGRIF_linear) 
    170 CALL Agrif_Set_interp(glamv_id,interp=AGRIF_linear) 
    171 CALL Agrif_Set_bc( glamv_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    172  
    173 CALL Agrif_Set_bcinterp(glamf_id,interp=AGRIF_linear) 
    174 CALL Agrif_Set_interp(glamf_id,interp=AGRIF_linear) 
    175 CALL Agrif_Set_bc( glamf_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    176  
    177 CALL Agrif_Set_bcinterp(gphit_id,interp=AGRIF_linear) 
    178 CALL Agrif_Set_interp(gphit_id,interp=AGRIF_linear) 
    179 CALL Agrif_Set_bc( gphit_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    180  
    181 CALL Agrif_Set_bcinterp(gphiu_id,interp=AGRIF_linear) 
    182 CALL Agrif_Set_interp(gphiu_id,interp=AGRIF_linear) 
    183 CALL Agrif_Set_bc( gphiu_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    184  
    185 CALL Agrif_Set_bcinterp(gphiv_id,interp=AGRIF_linear) 
    186 CALL Agrif_Set_interp(gphiv_id,interp=AGRIF_linear) 
    187 CALL Agrif_Set_bc( gphiv_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    188  
    189 CALL Agrif_Set_bcinterp(gphif_id,interp=AGRIF_linear) 
    190 CALL Agrif_Set_interp(gphif_id,interp=AGRIF_linear) 
    191 CALL Agrif_Set_bc( gphif_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    192  
    193 ! 
    194  
    195 CALL Agrif_Set_bcinterp(e1t_id,interp=AGRIF_ppm) 
    196 CALL Agrif_Set_interp(e1t_id,interp=AGRIF_ppm) 
    197 CALL Agrif_Set_bc( e1t_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    198  
    199 CALL Agrif_Set_bcinterp(e1u_id, interp1=Agrif_linear, interp2=AGRIF_ppm) 
    200 CALL Agrif_Set_interp(e1u_id, interp1=Agrif_linear, interp2=AGRIF_ppm) 
    201 CALL Agrif_Set_bc( e1u_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    202  
    203 CALL Agrif_Set_bcinterp(e1v_id,interp1=AGRIF_ppm, interp2=Agrif_linear) 
    204 CALL Agrif_Set_interp(e1v_id, interp1=AGRIF_ppm, interp2=Agrif_linear) 
    205 CALL Agrif_Set_bc( e1v_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    206  
    207 CALL Agrif_Set_bcinterp(e1f_id,interp=AGRIF_linear) 
    208 CALL Agrif_Set_interp(e1f_id,interp=AGRIF_linear) 
    209 CALL Agrif_Set_bc( e1f_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    210  
    211 CALL Agrif_Set_bcinterp(e2t_id,interp=AGRIF_ppm) 
    212 CALL Agrif_Set_interp(e2t_id,interp=AGRIF_ppm) 
    213 CALL Agrif_Set_bc( e2t_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    214  
    215 CALL Agrif_Set_bcinterp(e2u_id,interp1=Agrif_linear, interp2=AGRIF_ppm) 
    216 CALL Agrif_Set_interp(e2u_id,interp1=Agrif_linear, interp2=AGRIF_ppm) 
    217 CALL Agrif_Set_bc( e2u_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    218  
    219 CALL Agrif_Set_bcinterp(e2v_id,interp1=AGRIF_ppm, interp2=Agrif_linear) 
    220 CALL Agrif_Set_interp(e2v_id,interp1=AGRIF_ppm, interp2=Agrif_linear) 
    221 CALL Agrif_Set_bc( e2v_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    222  
    223 CALL Agrif_Set_bcinterp(e2f_id,interp=AGRIF_linear) 
    224 CALL Agrif_Set_interp(e2f_id,interp=AGRIF_linear) 
    225 CALL Agrif_Set_bc( e2f_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    226  
    227 CALL Agrif_Set_bcinterp(bathy_id,interp=AGRIF_linear) 
    228 CALL Agrif_Set_interp(bathy_id,interp=AGRIF_linear) 
    229 CALL Agrif_Set_bc( bathy_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    230  
    231 ! Vertical scale factors 
    232 CALL Agrif_Set_bcinterp(e3t_id,interp=AGRIF_ppm) 
    233 CALL Agrif_Set_interp(e3t_id,interp=AGRIF_ppm) 
    234 CALL Agrif_Set_bc( e3t_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    235 CALL Agrif_Set_Updatetype( e3t_id, update = AGRIF_Update_Average) 
    236  
    237 CALL Agrif_Set_bcinterp(e3t_copy_id,interp=AGRIF_constant) 
    238 CALL Agrif_Set_interp(e3t_copy_id,interp=AGRIF_constant) 
    239 CALL Agrif_Set_bc( e3t_copy_id, (/-npt_copy*irafx-1,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/)) 
    240  
    241 CALL Agrif_Set_bcinterp(e3t_connect_id,interp=AGRIF_ppm) 
    242 CALL Agrif_Set_interp(e3t_connect_id,interp=AGRIF_ppm) 
    243 CALL Agrif_Set_bc( e3t_connect_id, (/-(npt_copy+npt_connect)*irafx-1,-npt_copy*irafx-2/)) 
    244  
    245 CALL Agrif_Set_bcinterp(e3u_id, interp1=Agrif_linear, interp2=AGRIF_ppm) 
    246 CALL Agrif_Set_interp(e3u_id, interp1=Agrif_linear, interp2=AGRIF_ppm) 
    247 CALL Agrif_Set_bc( e3u_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    248 CALL Agrif_Set_Updatetype(e3u_id,update1 = Agrif_Update_Copy, update2 = Agrif_Update_Average) 
    249  
    250 CALL Agrif_Set_bcinterp(e3v_id,interp1=AGRIF_ppm, interp2=Agrif_linear) 
    251 CALL Agrif_Set_interp(e3v_id, interp1=AGRIF_ppm, interp2=Agrif_linear) 
    252 CALL Agrif_Set_bc( e3v_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
    253 CALL Agrif_Set_Updatetype(e3v_id,update1 = Agrif_Update_Average, update2 = Agrif_Update_Copy) 
    254     
    255 ! Bottom level 
    256 CALL Agrif_Set_bcinterp(bottom_level_id,interp=AGRIF_constant) 
    257 CALL Agrif_Set_interp(bottom_level_id,interp=AGRIF_constant) 
    258 CALL Agrif_Set_bc( bottom_level_id, (/-npt_copy*irafx-1,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/)) 
    259 CALL Agrif_Set_Updatetype( bottom_level_id, update = AGRIF_Update_Max) 
    260  
    261 end subroutine agrif_declare_var 
    262  
    263  
    264 subroutine agrif_init_lonlat() 
    265 use agrif_profiles 
    266 use agrif_util 
    267 external :: init_glamt, init_glamu, init_glamv, init_glamf 
    268 external :: init_gphit, init_gphiu, init_gphiv, init_gphif 
    269  
    270 call Agrif_Init_variable(glamt_id, procname = init_glamt) 
    271 call Agrif_Init_variable(glamu_id, procname = init_glamu) 
    272 call Agrif_Init_variable(glamv_id, procname = init_glamv) 
    273 call Agrif_Init_variable(glamf_id, procname = init_glamf) 
    274  
    275 call Agrif_Init_variable(gphit_id, procname = init_gphit) 
    276 call Agrif_Init_variable(gphiu_id, procname = init_gphiu) 
    277 call Agrif_Init_variable(gphiv_id, procname = init_gphiv) 
    278 call Agrif_Init_variable(gphif_id, procname = init_gphif) 
    279  
    280 end subroutine agrif_init_lonlat 
    281  
    282 subroutine agrif_init_scales() 
    283 use agrif_profiles 
    284 use agrif_util 
    285 external :: init_e1t, init_e1u, init_e1v, init_e1f 
    286 external :: init_e2t, init_e2u, init_e2v, init_e2f 
    287  
    288 call Agrif_Init_variable(e1t_id, procname = init_e1t) 
    289 call Agrif_Init_variable(e1u_id, procname = init_e1u) 
    290 call Agrif_Init_variable(e1v_id, procname = init_e1v) 
    291 call Agrif_Init_variable(e1f_id, procname = init_e1f) 
    292  
    293 call Agrif_Init_variable(e2t_id, procname = init_e2t) 
    294 call Agrif_Init_variable(e2u_id, procname = init_e2u) 
    295 call Agrif_Init_variable(e2v_id, procname = init_e2v) 
    296 call Agrif_Init_variable(e2f_id, procname = init_e2f) 
    297  
    298 end subroutine agrif_init_scales 
    299  
    300  
    301  
    302    SUBROUTINE init_glamt( ptab, i1, i2, j1, j2, before, nb,ndir) 
    303    use dom_oce 
     65 
     66      WRITE(*,*) ' ' 
     67      WRITE(*,*)'Size of the High resolution grid: ',nx,' x ',ny 
     68      WRITE(*,*) ' ' 
     69 
     70      CALL agrif_init_lonlat() 
     71      ln_perio = .FALSE. 
     72      IF( jperio == 1 .OR. jperio == 2 .OR. jperio == 4 ) ln_perio=.TRUE. 
     73 
     74      WHERE (glamt < -180) glamt = glamt +360. 
     75      WHERE (glamt > +180) glamt = glamt -360. 
     76 
     77      CALL lbc_lnk( 'glamt', glamt, 'T', 1._wp) 
     78      CALL lbc_lnk( 'gphit', gphit, 'T', 1._wp) 
     79 
     80      WHERE (glamu < -180) glamu = glamu +360. 
     81      WHERE (glamu > +180) glamu = glamu -360. 
     82      CALL lbc_lnk( 'glamu', glamu, 'U', 1._wp) 
     83      CALL lbc_lnk( 'gphiu', gphiu, 'U', 1._wp) 
     84 
     85      WHERE (glamv < -180) glamv = glamv +360. 
     86      WHERE (glamv > +180) glamv = glamv -360. 
     87      CALL lbc_lnk( 'glamv', glamv, 'V', 1._wp) 
     88      CALL lbc_lnk( 'gphiv', gphiv, 'V', 1._wp) 
     89 
     90      WHERE (glamf < -180) glamf = glamf +360. 
     91      WHERE (glamf > +180) glamf = glamf -360. 
     92      CALL lbc_lnk( 'glamf', glamf, 'F', 1._wp) 
     93      CALL lbc_lnk( 'gphif', gphif, 'F', 1._wp) 
     94 
     95      ! Correct South and North 
     96      IF ((nbondj == -1).OR.(nbondj == 2)) THEN 
     97         glamt(:,1) = glamt(:,2) 
     98         gphit(:,1) = gphit(:,2) 
     99         glamu(:,1) = glamu(:,2) 
     100         gphiu(:,1) = gphiu(:,2) 
     101         glamv(:,1) = glamv(:,2) 
     102         gphiv(:,1) = gphiv(:,2) 
     103      ENDIF 
     104      IF ((nbondj == 1).OR.(nbondj == 2)) THEN 
     105         glamt(:,jpj) = glamt(:,jpj-1) 
     106         gphit(:,jpj) = gphit(:,jpj-1) 
     107         glamu(:,jpj) = glamu(:,jpj-1) 
     108         gphiu(:,jpj) = gphiu(:,jpj-1) 
     109         glamv(:,jpj) = glamv(:,jpj-1) 
     110         gphiv(:,jpj) = gphiv(:,jpj-1) 
     111         glamf(:,jpj) = glamf(:,jpj-1) 
     112         gphif(:,jpj) = gphif(:,jpj-1) 
     113      ENDIF 
     114 
     115      ! Correct West and East 
     116      IF( jperio /= 1 ) THEN 
     117         IF((nbondi == -1) .OR. (nbondi == 2) ) THEN 
     118            glamt(1,:) = glamt(2,:) 
     119            gphit(1,:) = gphit(2,:) 
     120            glamu(1,:) = glamu(2,:) 
     121            gphiu(1,:) = gphiu(2,:) 
     122            glamv(1,:) = glamv(2,:) 
     123            gphiv(1,:) = gphiv(2,:) 
     124         ENDIF 
     125         IF( (nbondi == 1) .OR. (nbondi == 2) ) THEN 
     126            glamt(jpi,:) = glamt(jpi-1,:) 
     127            gphit(jpi,:) = gphit(jpi-1,:) 
     128            glamu(jpi,:) = glamu(jpi-1,:) 
     129            gphiu(jpi,:) = gphiu(jpi-1,:) 
     130            glamv(jpi,:) = glamv(jpi-1,:) 
     131            gphiv(jpi,:) = gphiv(jpi-1,:) 
     132            glamf(jpi,:) = glamf(jpi-1,:) 
     133            gphif(jpi,:) = gphif(jpi-1,:) 
     134         ENDIF 
     135      ENDIF 
     136 
     137      CALL agrif_init_scales() 
     138 
     139      ! Correct South and North 
     140      IF( (nbondj == -1) .OR. (nbondj == 2) ) THEN 
     141         e1t(:,1) = e1t(:,2) 
     142         e2t(:,1) = e2t(:,2) 
     143         e1u(:,1) = e1u(:,2) 
     144         e2u(:,1) = e2u(:,2) 
     145         e1v(:,1) = e1v(:,2) 
     146         e2v(:,1) = e2v(:,2) 
     147      ENDIF 
     148      IF( (nbondj == 1) .OR. (nbondj == 2) ) THEN 
     149         e1t(:,jpj) = e1t(:,jpj-1) 
     150         e2t(:,jpj) = e2t(:,jpj-1) 
     151         e1u(:,jpj) = e1u(:,jpj-1) 
     152         e2u(:,jpj) = e2u(:,jpj-1) 
     153         e1v(:,jpj) = e1v(:,jpj-1) 
     154         e2v(:,jpj) = e2v(:,jpj-1) 
     155         e1f(:,jpj) = e1f(:,jpj-1) 
     156         e2f(:,jpj) = e2f(:,jpj-1) 
     157      ENDIF 
     158 
     159      ! Correct West and East 
     160      IF( jperio /= 1 ) THEN 
     161         IF( (nbondj == -1) .OR. (nbondj == 2) ) THEN 
     162            e1t(1,:) = e1t(2,:) 
     163            e2t(1,:) = e2t(2,:) 
     164            e1u(1,:) = e1u(2,:) 
     165            e2u(1,:) = e2u(2,:) 
     166            e1v(1,:) = e1v(2,:) 
     167            e2v(1,:) = e2v(2,:) 
     168         ENDIF 
     169         IF( (nbondj == 1) .OR. (nbondj == 2) ) THEN 
     170            e1t(jpi,:) = e1t(jpi-1,:) 
     171            e2t(jpi,:) = e2t(jpi-1,:) 
     172            e1u(jpi,:) = e1u(jpi-1,:) 
     173            e2u(jpi,:) = e2u(jpi-1,:) 
     174            e1v(jpi,:) = e1v(jpi-1,:) 
     175            e2v(jpi,:) = e2v(jpi-1,:) 
     176            e1f(jpi,:) = e1f(jpi-1,:) 
     177            e2f(jpi,:) = e2f(jpi-1,:) 
     178         ENDIF 
     179      ENDIF 
     180 
     181   END SUBROUTINE Agrif_InitValues_cont 
     182 
     183 
     184   SUBROUTINE agrif_declare_var() 
     185      !!---------------------------------------------------------------------- 
     186      !!                 *** ROUTINE Agrif_InitValues_cont *** 
     187      !! 
     188      !! ** Purpose :: Declaration of variables to be interpolated 
     189      !!---------------------------------------------------------------------- 
     190      USE par_oce 
     191      USE dom_oce 
     192      USE agrif_profiles 
     193      USE agrif_parameters 
     194 
     195      IMPLICIT NONE 
     196 
     197      INTEGER :: ind1, ind2, ind3 
     198      INTEGER :: nx, ny 
     199      INTEGER ::nbghostcellsfine_tot_x, nbghostcellsfine_tot_y 
     200      INTEGER :: irafx 
     201 
     202      EXTERNAL :: nemo_mapping 
     203 
     204      ! 1. Declaration of the type of variable which have to be interpolated 
     205      !--------------------------------------------------------------------- 
     206 
     207      nx=nlci ; ny=nlcj 
     208 
     209      ind2 = 2 + nbghostcells_x 
     210      ind3 = 2 + nbghostcells_y_s 
     211      nbghostcellsfine_tot_x=nbghostcells_x+1 
     212      nbghostcellsfine_tot_y=max(nbghostcells_y_s,nbghostcells_y_n)+1 
     213 
     214      irafx = Agrif_irhox() 
     215 
     216      ! In case of East-West periodicity, prevent AGRIF interpolation at east and west boundaries 
     217      ! The procnames will not be CALLed at these boundaries 
     218      if (jperio == 1) THEN 
     219        CALL Agrif_Set_NearCommonBorderX(.TRUE.) 
     220        CALL Agrif_Set_DistantCommonBorderX(.TRUE.) 
     221      endif 
     222      if (.not.ln_bry_south) THEN 
     223        CALL Agrif_Set_NearCommonBorderY(.TRUE.) 
     224      endif 
     225 
     226      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),glamt_id) 
     227      CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),glamu_id) 
     228      CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),glamv_id) 
     229      CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),glamf_id) 
     230 
     231      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),gphit_id) 
     232      CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),gphiu_id) 
     233      CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),gphiv_id) 
     234      CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),gphif_id) 
     235 
     236      ! Horizontal scale factors 
     237 
     238      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),e1t_id) 
     239      CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),e1u_id) 
     240      CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),e1v_id) 
     241      CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),e1f_id) 
     242 
     243      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),e2t_id) 
     244      CALL agrif_declare_variable((/1,2/),(/ind2-1,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),e2u_id) 
     245      CALL agrif_declare_variable((/2,1/),(/ind2,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),e2v_id) 
     246      CALL agrif_declare_variable((/1,1/),(/ind2-1,ind3-1/),(/'x','y'/),(/1,1/),(/nx,ny/),e2f_id) 
     247 
     248      ! Bathymetry 
     249 
     250      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),bathy_id) 
     251 
     252      ! Vertical scale factors 
     253      CALL agrif_declare_variable((/2,2,0/),(/ind2,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/nx,ny,jpk/),e3t_id) 
     254      CALL agrif_declare_variable((/2,2,0/),(/ind2,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/nx,ny,jpk/),e3t_copy_id) 
     255      CALL agrif_declare_variable((/2,2,0/),(/ind2,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/nx,ny,jpk+1/),e3t_connect_id) 
     256 
     257      CALL agrif_declare_variable((/1,2,0/),(/ind2-1,ind3,0/),(/'x','y','N'/),(/1,1,1/),(/nx,ny,jpk/),e3u_id) 
     258      CALL agrif_declare_variable((/2,1,0/),(/ind2,ind3-1,0/),(/'x','y','N'/),(/1,1,1/),(/nx,ny,jpk/),e3v_id) 
     259 
     260      ! Bottom level 
     261 
     262      CALL agrif_declare_variable((/2,2/),(/ind2,ind3/),(/'x','y'/),(/1,1/),(/nx,ny/),bottom_level_id) 
     263 
     264      CALL Agrif_Set_bcinterp(glamt_id,interp=AGRIF_linear) 
     265      CALL Agrif_Set_interp(glamt_id,interp=AGRIF_linear) 
     266      CALL Agrif_Set_bc( glamt_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
     267 
     268      CALL Agrif_Set_bcinterp(glamu_id,interp=AGRIF_linear) 
     269      CALL Agrif_Set_interp(glamu_id,interp=AGRIF_linear) 
     270      CALL Agrif_Set_bc( glamu_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     271 
     272      CALL Agrif_Set_bcinterp(glamv_id,interp=AGRIF_linear) 
     273      CALL Agrif_Set_interp(glamv_id,interp=AGRIF_linear) 
     274      CALL Agrif_Set_bc( glamv_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     275 
     276      CALL Agrif_Set_bcinterp(glamf_id,interp=AGRIF_linear) 
     277      CALL Agrif_Set_interp(glamf_id,interp=AGRIF_linear) 
     278      CALL Agrif_Set_bc( glamf_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     279 
     280      CALL Agrif_Set_bcinterp(gphit_id,interp=AGRIF_linear) 
     281      CALL Agrif_Set_interp(gphit_id,interp=AGRIF_linear) 
     282      CALL Agrif_Set_bc( gphit_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
     283 
     284      CALL Agrif_Set_bcinterp(gphiu_id,interp=AGRIF_linear) 
     285      CALL Agrif_Set_interp(gphiu_id,interp=AGRIF_linear) 
     286      CALL Agrif_Set_bc( gphiu_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     287 
     288      CALL Agrif_Set_bcinterp(gphiv_id,interp=AGRIF_linear) 
     289      CALL Agrif_Set_interp(gphiv_id,interp=AGRIF_linear) 
     290      CALL Agrif_Set_bc( gphiv_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     291 
     292      CALL Agrif_Set_bcinterp(gphif_id,interp=AGRIF_linear) 
     293      CALL Agrif_Set_interp(gphif_id,interp=AGRIF_linear) 
     294      CALL Agrif_Set_bc( gphif_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     295 
     296      ! 
     297 
     298      CALL Agrif_Set_bcinterp(e1t_id,interp=AGRIF_ppm) 
     299      CALL Agrif_Set_interp(e1t_id,interp=AGRIF_ppm) 
     300      CALL Agrif_Set_bc( e1t_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
     301 
     302      CALL Agrif_Set_bcinterp(e1u_id, interp1=Agrif_linear, interp2=AGRIF_ppm) 
     303      CALL Agrif_Set_interp(e1u_id, interp1=Agrif_linear, interp2=AGRIF_ppm) 
     304      CALL Agrif_Set_bc( e1u_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     305 
     306      CALL Agrif_Set_bcinterp(e1v_id,interp1=AGRIF_ppm, interp2=Agrif_linear) 
     307      CALL Agrif_Set_interp(e1v_id, interp1=AGRIF_ppm, interp2=Agrif_linear) 
     308      CALL Agrif_Set_bc( e1v_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     309 
     310      CALL Agrif_Set_bcinterp(e1f_id,interp=AGRIF_linear) 
     311      CALL Agrif_Set_interp(e1f_id,interp=AGRIF_linear) 
     312      CALL Agrif_Set_bc( e1f_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     313 
     314      CALL Agrif_Set_bcinterp(e2t_id,interp=AGRIF_ppm) 
     315      CALL Agrif_Set_interp(e2t_id,interp=AGRIF_ppm) 
     316      CALL Agrif_Set_bc( e2t_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
     317 
     318      CALL Agrif_Set_bcinterp(e2u_id,interp1=Agrif_linear, interp2=AGRIF_ppm) 
     319      CALL Agrif_Set_interp(e2u_id,interp1=Agrif_linear, interp2=AGRIF_ppm) 
     320      CALL Agrif_Set_bc( e2u_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     321 
     322      CALL Agrif_Set_bcinterp(e2v_id,interp1=AGRIF_ppm, interp2=Agrif_linear) 
     323      CALL Agrif_Set_interp(e2v_id,interp1=AGRIF_ppm, interp2=Agrif_linear) 
     324      CALL Agrif_Set_bc( e2v_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     325 
     326      CALL Agrif_Set_bcinterp(e2f_id,interp=AGRIF_linear) 
     327      CALL Agrif_Set_interp(e2f_id,interp=AGRIF_linear) 
     328      CALL Agrif_Set_bc( e2f_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)/) ) 
     329 
     330      CALL Agrif_Set_bcinterp(bathy_id,interp=AGRIF_linear) 
     331      CALL Agrif_Set_interp(bathy_id,interp=AGRIF_linear) 
     332      CALL Agrif_Set_bc( bathy_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
     333 
     334      ! Vertical scale factors 
     335      CALL Agrif_Set_bcinterp(e3t_id,interp=AGRIF_ppm) 
     336      CALL Agrif_Set_interp(e3t_id,interp=AGRIF_ppm) 
     337      CALL Agrif_Set_bc( e3t_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
     338      CALL Agrif_Set_Updatetype( e3t_id, update = AGRIF_Update_Average) 
     339 
     340      CALL Agrif_Set_bcinterp(e3t_copy_id,interp=AGRIF_constant) 
     341      CALL Agrif_Set_interp(e3t_copy_id,interp=AGRIF_constant) 
     342      CALL Agrif_Set_bc( e3t_copy_id, (/-npt_copy*irafx,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/)) 
     343 
     344      CALL Agrif_Set_bcinterp(e3t_connect_id,interp=AGRIF_ppm) 
     345      CALL Agrif_Set_interp(e3t_connect_id,interp=AGRIF_ppm) 
     346      CALL Agrif_Set_bc( e3t_connect_id, (/-(npt_copy+npt_connect)*irafx-1,-npt_copy*irafx/)) 
     347 
     348      CALL Agrif_Set_bcinterp(e3u_id, interp1=Agrif_linear, interp2=AGRIF_ppm) 
     349      CALL Agrif_Set_interp(e3u_id, interp1=Agrif_linear, interp2=AGRIF_ppm) 
     350      CALL Agrif_Set_bc( e3u_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
     351      CALL Agrif_Set_Updatetype(e3u_id,update1 = Agrif_Update_Copy, update2 = Agrif_Update_Average) 
     352 
     353      CALL Agrif_Set_bcinterp(e3v_id,interp1=AGRIF_ppm, interp2=Agrif_linear) 
     354      CALL Agrif_Set_interp(e3v_id, interp1=AGRIF_ppm, interp2=Agrif_linear) 
     355      CALL Agrif_Set_bc( e3v_id, (/0,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/) ) 
     356      CALL Agrif_Set_Updatetype(e3v_id,update1 = Agrif_Update_Average, update2 = Agrif_Update_Copy) 
     357 
     358      ! Bottom level 
     359      CALL Agrif_Set_bcinterp(bottom_level_id,interp=AGRIF_constant) 
     360      CALL Agrif_Set_interp(bottom_level_id,interp=AGRIF_constant) 
     361      CALL Agrif_Set_bc( bottom_level_id, (/-npt_copy*irafx,max(nbghostcellsfine_tot_x,nbghostcellsfine_tot_y)-1/)) 
     362      CALL Agrif_Set_Updatetype( bottom_level_id, update = AGRIF_Update_Max) 
     363 
     364      CALL Agrif_Set_ExternalMapping(nemo_mapping) 
     365 
     366   END SUBROUTINE agrif_declare_var 
     367 
     368   SUBROUTINE nemo_mapping(ndim,ptx,pty,bounds,bounds_chunks,correction_required,nb_chunks) 
     369      USE dom_oce 
     370      INTEGER :: ndim 
     371      INTEGER :: ptx, pty 
     372      INTEGER,DIMENSION(ndim,2,2) :: bounds 
     373      INTEGER,DIMENSION(:,:,:,:),allocatable :: bounds_chunks 
     374      LOGICAL,DIMENSION(:),allocatable :: correction_required 
     375      INTEGER :: nb_chunks 
     376      INTEGER :: i 
     377 
     378      IF (agrif_debug_interp) THEN 
     379         DO i = 1, ndim 
     380             print *,'direction = ',i,bounds(i,1,2),bounds(i,2,2) 
     381         END DO 
     382      ENDIF 
     383 
     384      IF( bounds(2,2,2) > jpjglo ) THEN 
     385         IF( bounds(2,1,2) <= jpjglo ) THEN 
     386            nb_chunks = 2 
     387            ALLOCATE(bounds_chunks(nb_chunks,ndim,2,2)) 
     388            ALLOCATE(correction_required(nb_chunks)) 
     389            DO i = 1, nb_chunks 
     390               bounds_chunks(i,:,:,:) = bounds 
     391            END DO 
     392 
     393         ! FIRST CHUNCK (for j<=jpjglo) 
     394 
     395            ! Original indices 
     396            bounds_chunks(1,1,1,1) = bounds(1,1,2) 
     397            bounds_chunks(1,1,2,1) = bounds(1,2,2) 
     398            bounds_chunks(1,2,1,1) = bounds(2,1,2) 
     399            bounds_chunks(1,2,2,1) = jpjglo 
     400 
     401            bounds_chunks(1,1,1,2) = bounds(1,1,2) 
     402            bounds_chunks(1,1,2,2) = bounds(1,2,2) 
     403            bounds_chunks(1,2,1,2) = bounds(2,1,2) 
     404            bounds_chunks(1,2,2,2) = jpjglo 
     405 
     406            ! Correction required or not 
     407            correction_required(1)=.FALSE. 
     408 
     409         ! SECOND CHUNCK (for j>jpjglo) 
     410 
     411            !Original indices 
     412            bounds_chunks(2,1,1,1) = bounds(1,1,2) 
     413            bounds_chunks(2,1,2,1) = bounds(1,2,2) 
     414            bounds_chunks(2,2,1,1) = jpjglo-2 
     415            bounds_chunks(2,2,2,1) = bounds(2,2,2) 
     416 
     417           ! Where to find them 
     418           ! We use the relation TAB(ji,jj)=TAB(jpiglo-ji+2,jpjglo-2-(jj-jpjglo)) 
     419 
     420            IF (ptx == 2) THEN ! T, V points 
     421               bounds_chunks(2,1,1,2) = jpiglo-bounds(1,2,2)+2 
     422               bounds_chunks(2,1,2,2) = jpiglo-bounds(1,1,2)+2 
     423            ELSE ! U, F points 
     424               bounds_chunks(2,1,1,2) = jpiglo-bounds(1,2,2)+1 
     425               bounds_chunks(2,1,2,2) = jpiglo-bounds(1,1,2)+1 
     426            ENDIF 
     427 
     428            IF (pty == 2) THEN ! T, U points 
     429               bounds_chunks(2,2,1,2) = jpjglo-2-(bounds(2,2,2) -jpjglo) 
     430               bounds_chunks(2,2,2,2) = jpjglo-2-(jpjglo-2      -jpjglo) 
     431            ELSE ! V, F points 
     432               bounds_chunks(2,2,1,2) = jpjglo-3-(bounds(2,2,2) -jpjglo) 
     433               bounds_chunks(2,2,2,2) = jpjglo-3-(jpjglo-2      -jpjglo) 
     434            ENDIF 
     435       
     436            ! Correction required or not 
     437            correction_required(2)=.TRUE. 
     438 
     439         ELSE 
     440            
     441            nb_chunks = 1 
     442            ALLOCATE(bounds_chunks(nb_chunks,ndim,2,2)) 
     443            ALLOCATE(correction_required(nb_chunks)) 
     444            DO i=1,nb_chunks 
     445                bounds_chunks(i,:,:,:) = bounds 
     446            END DO 
     447 
     448            bounds_chunks(1,1,1,1) = bounds(1,1,2) 
     449            bounds_chunks(1,1,2,1) = bounds(1,2,2) 
     450            bounds_chunks(1,2,1,1) = bounds(2,1,2) 
     451            bounds_chunks(1,2,2,1) = bounds(2,2,2) 
     452 
     453            bounds_chunks(1,1,1,2) = jpiglo-bounds(1,2,2)+2 
     454            bounds_chunks(1,1,2,2) = jpiglo-bounds(1,1,2)+2 
     455 
     456            bounds_chunks(1,2,1,2) = jpjglo-2-(bounds(2,2,2)-jpjglo) 
     457            bounds_chunks(1,2,2,2) = jpjglo-2-(bounds(2,1,2)-jpjglo) 
     458 
     459            IF (ptx == 2) THEN ! T, V points 
     460               bounds_chunks(1,1,1,2) = jpiglo-bounds(1,2,2)+2 
     461               bounds_chunks(1,1,2,2) = jpiglo-bounds(1,1,2)+2 
     462            ELSE ! U, F points 
     463               bounds_chunks(1,1,1,2) = jpiglo-bounds(1,2,2)+1 
     464               bounds_chunks(1,1,2,2) = jpiglo-bounds(1,1,2)+1 
     465            ENDIF 
     466 
     467            IF (pty == 2) THEN ! T, U points 
     468               bounds_chunks(1,2,1,2) = jpjglo-2-(bounds(2,2,2) -jpjglo) 
     469               bounds_chunks(1,2,2,2) = jpjglo-2-(bounds(2,1,2) -jpjglo) 
     470            ELSE ! V, F points 
     471               bounds_chunks(1,2,1,2) = jpjglo-3-(bounds(2,2,2) -jpjglo) 
     472               bounds_chunks(1,2,2,2) = jpjglo-3-(bounds(2,1,2) -jpjglo) 
     473            ENDIF 
     474 
     475            correction_required(1)=.TRUE. 
     476 
     477         ENDIF  ! bounds(2,1,2) <= jpjglo 
     478 
     479      ELSE IF (bounds(1,1,2) < 1) THEN 
     480          
     481         IF (bounds(1,2,2) > 0) THEN 
     482            nb_chunks = 2 
     483            ALLOCATE(correction_required(nb_chunks)) 
     484            correction_required=.FALSE. 
     485            ALLOCATE(bounds_chunks(nb_chunks,ndim,2,2)) 
     486            DO i=1,nb_chunks 
     487               bounds_chunks(i,:,:,:) = bounds 
     488            END DO 
     489 
     490            bounds_chunks(1,1,1,2) = bounds(1,1,2)+jpiglo-2 
     491            bounds_chunks(1,1,2,2) = 1+jpiglo-2 
     492 
     493            bounds_chunks(1,1,1,1) = bounds(1,1,2) 
     494            bounds_chunks(1,1,2,1) = 1 
     495 
     496            bounds_chunks(2,1,1,2) = 2 
     497            bounds_chunks(2,1,2,2) = bounds(1,2,2) 
     498 
     499            bounds_chunks(2,1,1,1) = 2 
     500            bounds_chunks(2,1,2,1) = bounds(1,2,2) 
     501         ELSE 
     502            nb_chunks = 1 
     503            ALLOCATE(correction_required(nb_chunks)) 
     504            correction_required=.FALSE. 
     505            ALLOCATE(bounds_chunks(nb_chunks,ndim,2,2)) 
     506            DO i=1,nb_chunks 
     507               bounds_chunks(i,:,:,:) = bounds 
     508            END DO 
     509            bounds_chunks(1,1,1,2) = bounds(1,1,2)+jpiglo-2 
     510            bounds_chunks(1,1,2,2) = bounds(1,2,2)+jpiglo-2 
     511 
     512            bounds_chunks(1,1,1,1) = bounds(1,1,2) 
     513            bounds_chunks(1,1,2,1) = bounds(1,2,2) 
     514         ENDIF 
     515       
     516      ELSE 
     517       
     518         nb_chunks=1 
     519         ALLOCATE(correction_required(nb_chunks)) 
     520         correction_required=.FALSE. 
     521         ALLOCATE(bounds_chunks(nb_chunks,ndim,2,2)) 
     522         DO i=1,nb_chunks 
     523            bounds_chunks(i,:,:,:) = bounds 
     524         END DO 
     525         bounds_chunks(1,1,1,2) = bounds(1,1,2) 
     526         bounds_chunks(1,1,2,2) = bounds(1,2,2) 
     527         bounds_chunks(1,2,1,2) = bounds(2,1,2) 
     528         bounds_chunks(1,2,2,2) = bounds(2,2,2) 
     529 
     530         bounds_chunks(1,1,1,1) = bounds(1,1,2) 
     531         bounds_chunks(1,1,2,1) = bounds(1,2,2) 
     532         bounds_chunks(1,2,1,1) = bounds(2,1,2) 
     533         bounds_chunks(1,2,2,1) = bounds(2,2,2) 
     534 
     535      ENDIF 
     536 
     537   END SUBROUTINE nemo_mapping 
     538 
     539   FUNCTION agrif_external_switch_index(ptx,pty,i1,isens) 
     540      USE dom_oce 
     541      INTEGER :: ptx, pty, i1, isens 
     542      INTEGER :: agrif_external_switch_index 
     543 
     544      IF( isens == 1 )  THEN 
     545         IF( ptx == 2 ) THEN ! T, V points 
     546            agrif_external_switch_index = jpiglo-i1+2 
     547         ELSE ! U, F points 
     548            agrif_external_switch_index = jpiglo-i1+1 
     549         ENDIF 
     550      ELSE IF (isens ==2) THEN 
     551         IF (pty == 2) THEN ! T, U points 
     552            agrif_external_switch_index = jpjglo-2-(i1 -jpjglo) 
     553         ELSE ! V, F points 
     554            agrif_external_switch_index = jpjglo-3-(i1 -jpjglo) 
     555         ENDIF 
     556      ENDIF 
     557 
     558   END FUNCTION agrif_external_switch_index 
     559 
     560   SUBROUTINE correct_field(tab2d,i1,i2,j1,j2) 
     561      USE dom_oce 
     562      INTEGER :: i1,i2,j1,j2 
     563      REAL,DIMENSION(i1:i2,j1:j2) :: tab2d 
     564 
     565      INTEGER :: i,j 
     566      REAL,DIMENSION(i1:i2,j1:j2) :: tab2dtemp 
     567 
     568      tab2dtemp = tab2d 
     569 
     570      DO j=j1,j2 
     571         DO i=i1,i2 
     572        tab2d(i,j)=tab2dtemp(i2-(i-i1),j2-(j-j1)) 
     573         END DO 
     574      ENDDO 
     575 
     576   END SUBROUTINE correct_field 
     577 
     578   SUBROUTINE agrif_init_lonlat() 
     579      USE agrif_profiles 
     580      USE agrif_util 
     581      USE dom_oce 
     582      
     583      EXTERNAL :: init_glamt, init_glamu, init_glamv, init_glamf 
     584      EXTERNAL :: init_gphit, init_gphiu, init_gphiv, init_gphif 
     585      EXTERNAL :: longitude_linear_interp 
     586 
     587      INTEGER :: ji,jj,i1,i2,j1,j2 
     588      REAL, DIMENSION(jpi,jpj) :: tab2dtemp 
     589      INTEGER :: ind2,ind3 
     590      INTEGER :: irhox, irhoy 
     591 
     592      irhox = agrif_irhox() 
     593      irhoy = agrif_irhoy() 
     594      CALL Agrif_Set_external_linear_interp(longitude_linear_interp) 
     595 
     596      CALL Agrif_Init_variable(glamt_id, procname = init_glamt) 
     597      CALL Agrif_Init_variable(glamu_id, procname = init_glamu) 
     598      CALL Agrif_Init_variable(glamv_id, procname = init_glamv) 
     599      CALL Agrif_Init_variable(glamf_id, procname = init_glamf) 
     600      CALL Agrif_UnSet_external_linear_interp() 
     601 
     602      CALL Agrif_Init_variable(gphit_id, procname = init_gphit) 
     603      CALL Agrif_Init_variable(gphiu_id, procname = init_gphiu) 
     604      CALL Agrif_Init_variable(gphiv_id, procname = init_gphiv) 
     605      CALL Agrif_Init_variable(gphif_id, procname = init_gphif) 
     606 
     607   END SUBROUTINE agrif_init_lonlat 
     608 
     609   REAL FUNCTION longitude_linear_interp(x1,x2,coeff) 
     610      REAL :: x1, x2, coeff 
     611      REAL :: val_interp 
     612 
     613      IF( (x1*x2 <= -50*50) ) THEN 
     614         IF( x1 < 0 ) THEN 
     615            val_interp = coeff *(x1+360.) + (1.-coeff) *x2 
     616         ELSE 
     617            val_interp = coeff *x1 + (1.-coeff) *(x2+360.) 
     618         ENDIF 
     619         IF ((val_interp) >=180.) val_interp = val_interp - 360. 
     620      ELSE 
     621         val_interp = coeff * x1 + (1.-coeff) * x2 
     622      ENDIF 
     623      longitude_linear_interp = val_interp 
     624 
     625   END FUNCTION longitude_linear_interp 
     626 
     627   SUBROUTINE agrif_init_scales() 
     628      USE agrif_profiles 
     629      USE agrif_util 
     630      USE dom_oce 
     631      USE lbclnk 
     632      LOGICAL :: ln_perio 
     633      INTEGER nx,ny 
     634 
     635      EXTERNAL :: init_e1t, init_e1u, init_e1v, init_e1f 
     636      EXTERNAL :: init_e2t, init_e2u, init_e2v, init_e2f 
     637 
     638      nx = nlci ; ny = nlcj 
     639 
     640      ln_perio=.FALSE. 
     641      if( jperio ==1 .OR. jperio==2 .OR. jperio==4) ln_perio=.TRUE. 
     642 
     643      CALL Agrif_Init_variable(e1t_id, procname = init_e1t) 
     644      CALL Agrif_Init_variable(e1u_id, procname = init_e1u) 
     645      CALL Agrif_Init_variable(e1v_id, procname = init_e1v) 
     646      CALL Agrif_Init_variable(e1f_id, procname = init_e1f) 
     647 
     648      CALL Agrif_Init_variable(e2t_id, procname = init_e2t) 
     649      CALL Agrif_Init_variable(e2u_id, procname = init_e2u) 
     650      CALL Agrif_Init_variable(e2v_id, procname = init_e2v) 
     651      CALL Agrif_Init_variable(e2f_id, procname = init_e2f) 
     652 
     653      CALL lbc_lnk( 'e1t', e1t, 'T', 1._wp) 
     654      CALL lbc_lnk( 'e2t', e2t, 'T', 1._wp) 
     655      CALL lbc_lnk( 'e1u', e1u, 'U', 1._wp) 
     656      CALL lbc_lnk( 'e2u', e2u, 'U', 1._wp) 
     657      CALL lbc_lnk( 'e1v', e1v, 'V', 1._wp) 
     658      CALL lbc_lnk( 'e2v', e2v, 'V', 1._wp) 
     659      CALL lbc_lnk( 'e1f', e1f, 'F', 1._wp) 
     660      CALL lbc_lnk( 'e2f', e2f, 'F', 1._wp) 
     661 
     662   END SUBROUTINE agrif_init_scales 
     663 
     664   SUBROUTINE init_glamt( ptab, i1, i2, j1, j2,  before, nb,ndir) 
     665      USE dom_oce 
    304666      !!---------------------------------------------------------------------- 
    305667      !!                  ***  ROUTINE interpsshn  *** 
    306       !!----------------------------------------------------------------------   
     668      !!---------------------------------------------------------------------- 
     669      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     670      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     671      LOGICAL                         , INTENT(in   ) ::   before 
     672      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     673 
     674      ! 
     675      !!---------------------------------------------------------------------- 
     676      ! 
     677      IF( before) THEN 
     678         ptab(i1:i2,j1:j2) = glamt(i1:i2,j1:j2) 
     679      ELSE 
     680          glamt(i1:i2,j1:j2) = ptab(i1:i2,j1:j2) 
     681      ENDIF 
     682      ! 
     683   END SUBROUTINE init_glamt 
     684 
     685   SUBROUTINE init_glamu( ptab, i1, i2, j1, j2, before, nb,ndir) 
     686      USE dom_oce 
     687      !!---------------------------------------------------------------------- 
     688      !!                  ***  ROUTINE interpsshn  *** 
     689      !!---------------------------------------------------------------------- 
    307690      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    308691      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     
    311694      LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    312695      ! 
    313       !!----------------------------------------------------------------------   
    314       ! 
    315          western_side  = (nb == 1).AND.(ndir == 1) 
    316          eastern_side  = (nb == 1).AND.(ndir == 2) 
    317          southern_side = (nb == 2).AND.(ndir == 1) 
    318          northern_side = (nb == 2).AND.(ndir == 2) 
    319       IF( before) THEN 
    320          ptab(i1:i2,j1:j2) = glamt(i1:i2,j1:j2) 
    321       ELSE 
    322          glamt(i1:i2,j1:j2)=ptab 
    323       ENDIF 
    324       ! 
    325    END SUBROUTINE init_glamt 
    326  
    327     SUBROUTINE init_glamu( ptab, i1, i2, j1, j2, before, nb,ndir) 
    328     use dom_oce 
     696      !!---------------------------------------------------------------------- 
     697      ! 
     698      IF( before) THEN 
     699         ptab(i1:i2,j1:j2) = glamu(i1:i2,j1:j2) 
     700      ELSE 
     701          glamu(i1:i2,j1:j2) = ptab(i1:i2,j1:j2) 
     702      ENDIF 
     703      ! 
     704   END SUBROUTINE init_glamu 
     705 
     706   SUBROUTINE init_glamv( ptab, i1, i2, j1, j2, before, nb,ndir) 
     707      USE dom_oce 
    329708      !!---------------------------------------------------------------------- 
    330709      !!                  ***  ROUTINE interpsshn  *** 
    331       !!----------------------------------------------------------------------   
    332       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    333       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    334       LOGICAL                         , INTENT(in   ) ::   before 
    335       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    336       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    337       ! 
    338       !!----------------------------------------------------------------------   
    339       ! 
    340          western_side  = (nb == 1).AND.(ndir == 1) 
    341          eastern_side  = (nb == 1).AND.(ndir == 2) 
    342          southern_side = (nb == 2).AND.(ndir == 1) 
    343          northern_side = (nb == 2).AND.(ndir == 2) 
    344       IF( before) THEN 
    345          ptab(i1:i2,j1:j2) = glamu(i1:i2,j1:j2) 
    346       ELSE 
    347          glamu(i1:i2,j1:j2)=ptab 
    348       ENDIF 
    349       ! 
    350    END SUBROUTINE init_glamu 
    351  
    352    SUBROUTINE init_glamv( ptab, i1, i2, j1, j2, before, nb,ndir) 
    353    use dom_oce 
     710      !!---------------------------------------------------------------------- 
     711      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     712      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     713      LOGICAL                         , INTENT(in   ) ::   before 
     714      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     715      ! 
     716      !!---------------------------------------------------------------------- 
     717      ! 
     718      IF( before) THEN 
     719         ptab(i1:i2,j1:j2) = glamv(i1:i2,j1:j2) 
     720      ELSE 
     721          glamv(i1:i2,j1:j2) = ptab(i1:i2,j1:j2) 
     722      ENDIF 
     723      ! 
     724   END SUBROUTINE init_glamv 
     725 
     726   SUBROUTINE init_glamf( ptab, i1, i2, j1, j2,  before, nb,ndir) 
     727      USE dom_oce 
     728      !!---------------------------------------------------------------------- 
     729      !!                  ***  ROUTINE init_glamf  *** 
     730      !!---------------------------------------------------------------------- 
     731      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     732      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     733      LOGICAL                         , INTENT(in   ) ::   before 
     734      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     735      ! 
     736      !!---------------------------------------------------------------------- 
     737      ! 
     738      IF( before) THEN 
     739         ptab(i1:i2,j1:j2) = glamf(i1:i2,j1:j2) 
     740      ELSE 
     741          glamf(i1:i2,j1:j2) = ptab(i1:i2,j1:j2) 
     742      ENDIF 
     743      ! 
     744   END SUBROUTINE init_glamf 
     745 
     746   SUBROUTINE init_gphit( ptab, i1, i2, j1, j2, before, nb,ndir) 
     747      USE dom_oce 
     748      !!---------------------------------------------------------------------- 
     749      !!                  ***  ROUTINE init_gphit  *** 
     750      !!---------------------------------------------------------------------- 
     751      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     752      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     753      LOGICAL                         , INTENT(in   ) ::   before 
     754      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     755      ! 
     756      !!---------------------------------------------------------------------- 
     757      ! 
     758      IF( before) THEN 
     759         ptab(i1:i2,j1:j2) = gphit(i1:i2,j1:j2) 
     760      ELSE 
     761         gphit(i1:i2,j1:j2)=ptab 
     762      ENDIF 
     763      ! 
     764   END SUBROUTINE init_gphit 
     765 
     766   SUBROUTINE init_gphiu( ptab, i1, i2, j1, j2, before, nb,ndir) 
     767      USE dom_oce 
     768      !!---------------------------------------------------------------------- 
     769      !!                  ***  ROUTINE init_gphiu  *** 
     770      !!---------------------------------------------------------------------- 
     771      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     772      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     773      LOGICAL                         , INTENT(in   ) ::   before 
     774      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     775      ! 
     776      !!---------------------------------------------------------------------- 
     777      ! 
     778      IF( before) THEN 
     779         ptab(i1:i2,j1:j2) = gphiu(i1:i2,j1:j2) 
     780      ELSE 
     781         gphiu(i1:i2,j1:j2)=ptab 
     782      ENDIF 
     783      ! 
     784   END SUBROUTINE init_gphiu 
     785 
     786   SUBROUTINE init_gphiv( ptab, i1, i2, j1, j2, before, nb,ndir) 
     787      USE dom_oce 
     788      !!---------------------------------------------------------------------- 
     789      !!                  ***  ROUTINE init_gphiv  *** 
     790      !!---------------------------------------------------------------------- 
     791      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     792      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     793      LOGICAL                         , INTENT(in   ) ::   before 
     794      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     795      ! 
     796      !!---------------------------------------------------------------------- 
     797 
     798      IF( before) THEN 
     799         ptab(i1:i2,j1:j2) = gphiv(i1:i2,j1:j2) 
     800      ELSE 
     801         gphiv(i1:i2,j1:j2)=ptab 
     802      ENDIF 
     803      ! 
     804   END SUBROUTINE init_gphiv 
     805 
     806 
     807   SUBROUTINE init_gphif( ptab, i1, i2, j1, j2, before, nb,ndir) 
     808      USE dom_oce 
     809      !!---------------------------------------------------------------------- 
     810      !!                  ***  ROUTINE init_gphif  *** 
     811      !!---------------------------------------------------------------------- 
     812      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     813      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     814      LOGICAL                         , INTENT(in   ) ::   before 
     815      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     816      ! 
     817      !!---------------------------------------------------------------------- 
     818      ! 
     819      IF( before) THEN 
     820         ptab(i1:i2,j1:j2) = gphif(i1:i2,j1:j2) 
     821      ELSE 
     822         gphif(i1:i2,j1:j2)=ptab 
     823      ENDIF 
     824      ! 
     825   END SUBROUTINE init_gphif 
     826 
     827 
     828   SUBROUTINE init_e1t( ptab, i1, i2, j1, j2, before, nb,ndir) 
     829      USE dom_oce 
     830      USE agrif_parameters 
     831      !!---------------------------------------------------------------------- 
     832      !!                  ***  ROUTINE init_e1t  *** 
     833      !!---------------------------------------------------------------------- 
     834      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     835      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     836      LOGICAL                         , INTENT(in   ) ::   before 
     837      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     838      ! 
     839      !!---------------------------------------------------------------------- 
     840      ! 
     841      INTEGER :: jj 
     842 
     843      IF( before) THEN 
     844        ! May need to extend at south boundary 
     845          IF (j1<1) THEN 
     846            IF (.NOT.agrif_child(ln_bry_south)) THEN 
     847              IF ((nbondj == -1).OR.(nbondj == 2)) THEN 
     848                DO jj=1,j2 
     849                  ptab(i1:i2,jj)=e1t(i1:i2,jj) 
     850                ENDDO 
     851                DO jj=j1,0 
     852                  ptab(i1:i2,jj)=e1t(i1:i2,1) 
     853                ENDDO 
     854              ENDIF 
     855            ELSE 
     856              stop "OUT OF BOUNDS" 
     857            ENDIF 
     858          ELSE 
     859             ptab(i1:i2,j1:j2) = e1t(i1:i2,j1:j2) 
     860          ENDIF 
     861      ELSE 
     862         e1t(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
     863      ENDIF 
     864      ! 
     865   END SUBROUTINE init_e1t 
     866 
     867   SUBROUTINE init_e1u( ptab, i1, i2, j1, j2, before, nb,ndir) 
     868      USE dom_oce 
     869      USE agrif_parameters 
     870      !!---------------------------------------------------------------------- 
     871      !!                  ***  ROUTINE init_e1u  *** 
     872      !!---------------------------------------------------------------------- 
     873      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     874      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     875      LOGICAL                         , INTENT(in   ) ::   before 
     876      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     877      ! 
     878      !!---------------------------------------------------------------------- 
     879      ! 
     880      INTEGER :: jj 
     881 
     882      IF( before) THEN 
     883          IF (j1<1) THEN 
     884            IF (.NOT.agrif_child(ln_bry_south)) THEN 
     885              IF ((nbondj == -1).OR.(nbondj == 2)) THEN 
     886                DO jj=1,j2 
     887                  ptab(i1:i2,jj)=e1u(i1:i2,jj) 
     888                ENDDO 
     889                DO jj=j1,0 
     890                  ptab(i1:i2,jj)=e1u(i1:i2,1) 
     891                ENDDO 
     892              ENDIF 
     893            ELSE 
     894              stop "OUT OF BOUNDS" 
     895            ENDIF 
     896          ELSE 
     897             ptab(i1:i2,j1:j2) = e1u(i1:i2,j1:j2) 
     898          ENDIF 
     899      ELSE 
     900         e1u(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
     901      ENDIF 
     902      ! 
     903   END SUBROUTINE init_e1u 
     904 
     905   SUBROUTINE init_e1v( ptab, i1, i2, j1, j2, before, nb,ndir) 
     906      USE dom_oce 
     907      !!---------------------------------------------------------------------- 
     908      !!                  ***  ROUTINE init_e1v  *** 
     909      !!---------------------------------------------------------------------- 
     910      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     911      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     912      LOGICAL                         , INTENT(in   ) ::   before 
     913      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     914      ! 
     915      !!---------------------------------------------------------------------- 
     916      ! 
     917      IF( before) THEN 
     918         ptab(i1:i2,j1:j2) = e1v(i1:i2,j1:j2) 
     919      ELSE 
     920         e1v(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
     921      ENDIF 
     922      ! 
     923   END SUBROUTINE init_e1v 
     924 
     925   SUBROUTINE init_e1f( ptab, i1, i2, j1, j2, before, nb,ndir) 
     926      USE dom_oce 
     927      !!---------------------------------------------------------------------- 
     928      !!                  ***  ROUTINE init_e1f  *** 
     929      !!---------------------------------------------------------------------- 
     930      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     931      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     932      LOGICAL                         , INTENT(in   ) ::   before 
     933      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     934      ! 
     935      !!---------------------------------------------------------------------- 
     936      ! 
     937      IF( before) THEN 
     938         ptab(i1:i2,j1:j2) = e1f(i1:i2,j1:j2) 
     939      ELSE 
     940         e1f(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
     941      ENDIF 
     942      ! 
     943   END SUBROUTINE init_e1f 
     944 
     945   SUBROUTINE init_e2t( ptab, i1, i2, j1, j2, before, nb,ndir) 
     946      USE dom_oce 
     947      USE agrif_parameters 
     948      !!---------------------------------------------------------------------- 
     949      !!                  ***  ROUTINE init_e2t  *** 
     950      !!---------------------------------------------------------------------- 
     951      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     952      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     953      LOGICAL                         , INTENT(in   ) ::   before 
     954      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     955      ! 
     956      !!---------------------------------------------------------------------- 
     957      ! 
     958      INTEGER :: jj 
     959 
     960      IF( before) THEN 
     961          IF (j1<1) THEN 
     962            IF (.NOT.agrif_child(ln_bry_south)) THEN 
     963              IF ((nbondj == -1).OR.(nbondj == 2)) THEN 
     964                DO jj=1,j2 
     965                  ptab(i1:i2,jj)=e2t(i1:i2,jj) 
     966                ENDDO 
     967                DO jj=j1,0 
     968                  ptab(i1:i2,jj)=e2t(i1:i2,1) 
     969                ENDDO 
     970              ENDIF 
     971            ELSE 
     972              stop "OUT OF BOUNDS" 
     973            ENDIF 
     974          ELSE 
     975             ptab(i1:i2,j1:j2) = e2t(i1:i2,j1:j2) 
     976          ENDIF 
     977      ELSE 
     978         e2t(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
     979      ENDIF 
     980      ! 
     981   END SUBROUTINE init_e2t 
     982 
     983   SUBROUTINE init_e2u( ptab, i1, i2, j1, j2, before, nb,ndir) 
     984   USE dom_oce 
     985   USE agrif_parameters 
    354986      !!---------------------------------------------------------------------- 
    355987      !!                  ***  ROUTINE interpsshn  *** 
    356       !!----------------------------------------------------------------------   
    357       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    358       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    359       LOGICAL                         , INTENT(in   ) ::   before 
    360       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    361       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    362       ! 
    363       !!----------------------------------------------------------------------   
    364       ! 
    365          western_side  = (nb == 1).AND.(ndir == 1) 
    366          eastern_side  = (nb == 1).AND.(ndir == 2) 
    367          southern_side = (nb == 2).AND.(ndir == 1) 
    368          northern_side = (nb == 2).AND.(ndir == 2) 
    369       IF( before) THEN 
    370          ptab(i1:i2,j1:j2) = glamv(i1:i2,j1:j2) 
    371       ELSE 
    372          glamv(i1:i2,j1:j2)=ptab 
    373       ENDIF 
    374       ! 
    375    END SUBROUTINE init_glamv 
    376  
    377    SUBROUTINE init_glamf( ptab, i1, i2, j1, j2, before, nb,ndir) 
    378    use dom_oce 
     988      !!---------------------------------------------------------------------- 
     989      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     990      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     991      LOGICAL                         , INTENT(in   ) ::   before 
     992      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     993      ! 
     994      !!---------------------------------------------------------------------- 
     995      ! 
     996      INTEGER :: jj 
     997 
     998      IF( before) THEN 
     999          IF (j1<1) THEN 
     1000            IF (.NOT.agrif_child(ln_bry_south)) THEN 
     1001              IF ((nbondj == -1).OR.(nbondj == 2)) THEN 
     1002                DO jj=1,j2 
     1003                  ptab(i1:i2,jj)=e2u(i1:i2,jj) 
     1004                ENDDO 
     1005                DO jj=j1,0 
     1006                  ptab(i1:i2,jj)=e2u(i1:i2,1) 
     1007                ENDDO 
     1008              ENDIF 
     1009            ELSE 
     1010              stop "OUT OF BOUNDS" 
     1011            ENDIF 
     1012          ELSE 
     1013             ptab(i1:i2,j1:j2) = e2u(i1:i2,j1:j2) 
     1014          ENDIF 
     1015      ELSE 
     1016         e2u(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
     1017      ENDIF 
     1018      ! 
     1019   END SUBROUTINE init_e2u 
     1020 
     1021   SUBROUTINE init_e2v( ptab, i1, i2, j1, j2, before, nb,ndir) 
     1022      USE dom_oce 
    3791023      !!---------------------------------------------------------------------- 
    3801024      !!                  ***  ROUTINE interpsshn  *** 
    381       !!----------------------------------------------------------------------   
    382       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    383       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    384       LOGICAL                         , INTENT(in   ) ::   before 
    385       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    386       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    387       ! 
    388       !!----------------------------------------------------------------------   
    389       ! 
    390          western_side  = (nb == 1).AND.(ndir == 1) 
    391          eastern_side  = (nb == 1).AND.(ndir == 2) 
    392          southern_side = (nb == 2).AND.(ndir == 1) 
    393          northern_side = (nb == 2).AND.(ndir == 2) 
    394       IF( before) THEN 
    395          ptab(i1:i2,j1:j2) = glamf(i1:i2,j1:j2) 
    396       ELSE 
    397          glamf(i1:i2,j1:j2)=ptab 
    398       ENDIF 
    399       ! 
    400    END SUBROUTINE init_glamf 
    401  
    402    SUBROUTINE init_gphit( ptab, i1, i2, j1, j2, before, nb,ndir) 
    403    use dom_oce 
     1025      !!---------------------------------------------------------------------- 
     1026      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     1027      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     1028      LOGICAL                         , INTENT(in   ) ::   before 
     1029      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     1030      ! 
     1031      !!---------------------------------------------------------------------- 
     1032      IF( before) THEN 
     1033         ptab(i1:i2,j1:j2) = e2v(i1:i2,j1:j2) 
     1034      ELSE 
     1035         e2v(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
     1036      ENDIF 
     1037      ! 
     1038   END SUBROUTINE init_e2v 
     1039 
     1040   SUBROUTINE init_e2f( ptab, i1, i2, j1, j2, before, nb,ndir) 
     1041   USE dom_oce 
    4041042      !!---------------------------------------------------------------------- 
    4051043      !!                  ***  ROUTINE interpsshn  *** 
    406       !!----------------------------------------------------------------------   
    407       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    408       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    409       LOGICAL                         , INTENT(in   ) ::   before 
    410       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    411       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    412       ! 
    413       !!----------------------------------------------------------------------   
    414       ! 
    415          western_side  = (nb == 1).AND.(ndir == 1) 
    416          eastern_side  = (nb == 1).AND.(ndir == 2) 
    417          southern_side = (nb == 2).AND.(ndir == 1) 
    418          northern_side = (nb == 2).AND.(ndir == 2) 
    419       IF( before) THEN 
    420          ptab(i1:i2,j1:j2) = gphit(i1:i2,j1:j2) 
    421       ELSE 
    422          gphit(i1:i2,j1:j2)=ptab 
    423       ENDIF 
    424       ! 
    425    END SUBROUTINE init_gphit 
    426  
    427     SUBROUTINE init_gphiu( ptab, i1, i2, j1, j2, before, nb,ndir) 
    428     use dom_oce 
    429       !!---------------------------------------------------------------------- 
    430       !!                  ***  ROUTINE interpsshn  *** 
    431       !!----------------------------------------------------------------------   
    432       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    433       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    434       LOGICAL                         , INTENT(in   ) ::   before 
    435       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    436       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    437       ! 
    438       !!----------------------------------------------------------------------   
    439       ! 
    440          western_side  = (nb == 1).AND.(ndir == 1) 
    441          eastern_side  = (nb == 1).AND.(ndir == 2) 
    442          southern_side = (nb == 2).AND.(ndir == 1) 
    443          northern_side = (nb == 2).AND.(ndir == 2) 
    444       IF( before) THEN 
    445          ptab(i1:i2,j1:j2) = gphiu(i1:i2,j1:j2) 
    446       ELSE 
    447          gphiu(i1:i2,j1:j2)=ptab 
    448       ENDIF 
    449       ! 
    450    END SUBROUTINE init_gphiu 
    451  
    452     SUBROUTINE init_gphiv( ptab, i1, i2, j1, j2, before, nb,ndir) 
    453     use dom_oce 
    454       !!---------------------------------------------------------------------- 
    455       !!                  ***  ROUTINE interpsshn  *** 
    456       !!----------------------------------------------------------------------   
    457       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    458       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    459       LOGICAL                         , INTENT(in   ) ::   before 
    460       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    461       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    462       ! 
    463       !!----------------------------------------------------------------------   
    464       ! 
    465          western_side  = (nb == 1).AND.(ndir == 1) 
    466          eastern_side  = (nb == 1).AND.(ndir == 2) 
    467          southern_side = (nb == 2).AND.(ndir == 1) 
    468          northern_side = (nb == 2).AND.(ndir == 2) 
    469       IF( before) THEN 
    470          ptab(i1:i2,j1:j2) = gphiv(i1:i2,j1:j2) 
    471       ELSE 
    472          gphiv(i1:i2,j1:j2)=ptab 
    473       ENDIF 
    474       ! 
    475    END SUBROUTINE init_gphiv 
    476  
    477  
    478    SUBROUTINE init_gphif( ptab, i1, i2, j1, j2, before, nb,ndir) 
    479    use dom_oce 
    480       !!---------------------------------------------------------------------- 
    481       !!                  ***  ROUTINE interpsshn  *** 
    482       !!----------------------------------------------------------------------   
    483       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    484       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    485       LOGICAL                         , INTENT(in   ) ::   before 
    486       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    487       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    488       ! 
    489       !!----------------------------------------------------------------------   
    490       ! 
    491          western_side  = (nb == 1).AND.(ndir == 1) 
    492          eastern_side  = (nb == 1).AND.(ndir == 2) 
    493          southern_side = (nb == 2).AND.(ndir == 1) 
    494          northern_side = (nb == 2).AND.(ndir == 2) 
    495       IF( before) THEN 
    496          ptab(i1:i2,j1:j2) = gphif(i1:i2,j1:j2) 
    497       ELSE 
    498          gphif(i1:i2,j1:j2)=ptab 
    499       ENDIF 
    500       ! 
    501    END SUBROUTINE init_gphif 
    502  
    503  
    504    SUBROUTINE init_e1t( ptab, i1, i2, j1, j2, before, nb,ndir) 
    505    use dom_oce 
    506       !!---------------------------------------------------------------------- 
    507       !!                  ***  ROUTINE interpsshn  *** 
    508       !!----------------------------------------------------------------------   
    509       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    510       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    511       LOGICAL                         , INTENT(in   ) ::   before 
    512       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    513       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    514       ! 
    515       !!----------------------------------------------------------------------   
    516       ! 
    517          western_side  = (nb == 1).AND.(ndir == 1) 
    518          eastern_side  = (nb == 1).AND.(ndir == 2) 
    519          southern_side = (nb == 2).AND.(ndir == 1) 
    520          northern_side = (nb == 2).AND.(ndir == 2) 
    521       IF( before) THEN 
    522          ptab(i1:i2,j1:j2) = e1t(i1:i2,j1:j2) 
    523       ELSE 
    524          e1t(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
    525       ENDIF 
    526       ! 
    527    END SUBROUTINE init_e1t 
    528  
    529    SUBROUTINE init_e1u( ptab, i1, i2, j1, j2, before, nb,ndir) 
    530    use dom_oce 
    531       !!---------------------------------------------------------------------- 
    532       !!                  ***  ROUTINE interpsshn  *** 
    533       !!----------------------------------------------------------------------   
    534       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    535       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    536       LOGICAL                         , INTENT(in   ) ::   before 
    537       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    538       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    539       ! 
    540       !!----------------------------------------------------------------------   
    541       ! 
    542          western_side  = (nb == 1).AND.(ndir == 1) 
    543          eastern_side  = (nb == 1).AND.(ndir == 2) 
    544          southern_side = (nb == 2).AND.(ndir == 1) 
    545          northern_side = (nb == 2).AND.(ndir == 2) 
    546       IF( before) THEN 
    547          ptab(i1:i2,j1:j2) = e1u(i1:i2,j1:j2) 
    548       ELSE 
    549          e1u(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
    550       ENDIF 
    551       ! 
    552    END SUBROUTINE init_e1u 
    553  
    554    SUBROUTINE init_e1v( ptab, i1, i2, j1, j2, before, nb,ndir) 
    555    use dom_oce 
    556       !!---------------------------------------------------------------------- 
    557       !!                  ***  ROUTINE interpsshn  *** 
    558       !!----------------------------------------------------------------------   
    559       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    560       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    561       LOGICAL                         , INTENT(in   ) ::   before 
    562       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    563       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    564       ! 
    565       !!----------------------------------------------------------------------   
    566       ! 
    567          western_side  = (nb == 1).AND.(ndir == 1) 
    568          eastern_side  = (nb == 1).AND.(ndir == 2) 
    569          southern_side = (nb == 2).AND.(ndir == 1) 
    570          northern_side = (nb == 2).AND.(ndir == 2) 
    571       IF( before) THEN 
    572          ptab(i1:i2,j1:j2) = e1v(i1:i2,j1:j2) 
    573       ELSE 
    574          e1v(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
    575       ENDIF 
    576       ! 
    577    END SUBROUTINE init_e1v 
    578  
    579    SUBROUTINE init_e1f( ptab, i1, i2, j1, j2, before, nb,ndir) 
    580    use dom_oce 
    581       !!---------------------------------------------------------------------- 
    582       !!                  ***  ROUTINE interpsshn  *** 
    583       !!----------------------------------------------------------------------   
    584       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    585       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    586       LOGICAL                         , INTENT(in   ) ::   before 
    587       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    588       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    589       ! 
    590       !!----------------------------------------------------------------------   
    591       ! 
    592          western_side  = (nb == 1).AND.(ndir == 1) 
    593          eastern_side  = (nb == 1).AND.(ndir == 2) 
    594          southern_side = (nb == 2).AND.(ndir == 1) 
    595          northern_side = (nb == 2).AND.(ndir == 2) 
    596       IF( before) THEN 
    597          ptab(i1:i2,j1:j2) = e1f(i1:i2,j1:j2) 
    598       ELSE 
    599          e1f(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
    600       ENDIF 
    601       ! 
    602    END SUBROUTINE init_e1f 
    603  
    604   SUBROUTINE init_e2t( ptab, i1, i2, j1, j2, before, nb,ndir) 
    605    use dom_oce 
    606       !!---------------------------------------------------------------------- 
    607       !!                  ***  ROUTINE interpsshn  *** 
    608       !!----------------------------------------------------------------------   
    609       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    610       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    611       LOGICAL                         , INTENT(in   ) ::   before 
    612       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    613       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    614       ! 
    615       !!----------------------------------------------------------------------   
    616       ! 
    617          western_side  = (nb == 1).AND.(ndir == 1) 
    618          eastern_side  = (nb == 1).AND.(ndir == 2) 
    619          southern_side = (nb == 2).AND.(ndir == 1) 
    620          northern_side = (nb == 2).AND.(ndir == 2) 
    621       IF( before) THEN 
    622          ptab(i1:i2,j1:j2) = e2t(i1:i2,j1:j2) 
    623       ELSE 
    624          e2t(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
    625       ENDIF 
    626       ! 
    627    END SUBROUTINE init_e2t 
    628  
    629    SUBROUTINE init_e2u( ptab, i1, i2, j1, j2, before, nb,ndir) 
    630    use dom_oce 
    631       !!---------------------------------------------------------------------- 
    632       !!                  ***  ROUTINE interpsshn  *** 
    633       !!----------------------------------------------------------------------   
    634       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    635       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    636       LOGICAL                         , INTENT(in   ) ::   before 
    637       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    638       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    639       ! 
    640       !!----------------------------------------------------------------------   
    641       ! 
    642          western_side  = (nb == 1).AND.(ndir == 1) 
    643          eastern_side  = (nb == 1).AND.(ndir == 2) 
    644          southern_side = (nb == 2).AND.(ndir == 1) 
    645          northern_side = (nb == 2).AND.(ndir == 2) 
    646       IF( before) THEN 
    647          ptab(i1:i2,j1:j2) = e2u(i1:i2,j1:j2) 
    648       ELSE 
    649          e2u(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
    650       ENDIF 
    651       ! 
    652    END SUBROUTINE init_e2u 
    653  
    654    SUBROUTINE init_e2v( ptab, i1, i2, j1, j2, before, nb,ndir) 
    655    use dom_oce 
    656       !!---------------------------------------------------------------------- 
    657       !!                  ***  ROUTINE interpsshn  *** 
    658       !!----------------------------------------------------------------------   
    659       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    660       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    661       LOGICAL                         , INTENT(in   ) ::   before 
    662       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    663       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    664       ! 
    665       !!----------------------------------------------------------------------   
    666       ! 
    667          western_side  = (nb == 1).AND.(ndir == 1) 
    668          eastern_side  = (nb == 1).AND.(ndir == 2) 
    669          southern_side = (nb == 2).AND.(ndir == 1) 
    670          northern_side = (nb == 2).AND.(ndir == 2) 
    671       IF( before) THEN 
    672          ptab(i1:i2,j1:j2) = e2v(i1:i2,j1:j2) 
    673       ELSE 
    674          e2v(i1:i2,j1:j2)=ptab/Agrif_rhoy() 
    675       ENDIF 
    676       ! 
    677    END SUBROUTINE init_e2v 
    678  
    679    SUBROUTINE init_e2f( ptab, i1, i2, j1, j2, before, nb,ndir) 
    680    use dom_oce 
    681       !!---------------------------------------------------------------------- 
    682       !!                  ***  ROUTINE interpsshn  *** 
    683       !!----------------------------------------------------------------------   
    684       INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
    685       REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
    686       LOGICAL                         , INTENT(in   ) ::   before 
    687       INTEGER                         , INTENT(in   ) ::   nb , ndir 
    688       LOGICAL  ::   western_side, eastern_side,northern_side,southern_side 
    689       ! 
    690       !!----------------------------------------------------------------------   
    691       ! 
    692          western_side  = (nb == 1).AND.(ndir == 1) 
    693          eastern_side  = (nb == 1).AND.(ndir == 2) 
    694          southern_side = (nb == 2).AND.(ndir == 1) 
    695          northern_side = (nb == 2).AND.(ndir == 2) 
     1044      !!---------------------------------------------------------------------- 
     1045      INTEGER                         , INTENT(in   ) ::   i1, i2, j1, j2 
     1046      REAL, DIMENSION(i1:i2,j1:j2), INTENT(inout) ::   ptab 
     1047      LOGICAL                         , INTENT(in   ) ::   before 
     1048      INTEGER                         , INTENT(in   ) ::   nb , ndir 
     1049      ! 
     1050      !!---------------------------------------------------------------------- 
     1051      ! 
    6961052      IF( before) THEN 
    6971053         ptab(i1:i2,j1:j2) = e2f(i1:i2,j1:j2) 
     
    7031059 
    7041060 
    705 SUBROUTINE agrif_nemo_init 
    706 USE agrif_parameters 
    707 USE in_out_manager 
    708 USE lib_mpp 
    709  
    710     
    711    !! 
    712    IMPLICIT NONE 
    713     
    714    INTEGER ::   ios 
    715     
    716    NAMELIST/namagrif/ nn_cln_update,ln_spc_dyn,rn_sponge_tra,rn_sponge_dyn,ln_chk_bathy,npt_connect, npt_copy 
     1061   SUBROUTINE agrif_nemo_init 
     1062      USE agrif_parameters 
     1063      USE dom_oce 
     1064      USE in_out_manager 
     1065      USE lib_mpp 
     1066      !! 
     1067      IMPLICIT NONE 
     1068 
     1069      INTEGER ::   ios 
     1070 
     1071      NAMELIST/namagrif/ nn_cln_update,ln_spc_dyn,rn_sponge_tra,rn_sponge_dyn,ln_chk_bathy,npt_connect,   & 
     1072      &  npt_copy, ln_bry_south 
    7171073 
    7181074      REWIND( numnam_ref )              ! Namelist namagrif in reference namelist : nesting parameters 
     
    7301086         WRITE(numout,*) '~~~~~~~' 
    7311087         WRITE(numout,*) '   Namelist namagrif : set nesting parameters' 
    732          WRITE(numout,*) '      npt_copy     = ', npt_copy 
    733          WRITE(numout,*) '      npt_connect  = ', npt_connect 
    734       ENDIF 
    735        
    736 END SUBROUTINE agrif_nemo_init 
    737     
    738  
    739 SUBROUTINE Agrif_detect( kg, ksizex ) 
     1088         WRITE(numout,*) '      npt_copy             = ', npt_copy 
     1089         WRITE(numout,*) '      npt_connect          = ', npt_connect 
     1090         WRITE(numout,*) '      ln_bry_south  = ', ln_bry_south 
     1091      ENDIF 
     1092 
     1093   ! Set the number of ghost cells according to periodicity 
     1094 
     1095      nbghostcells_x = nbghostcells 
     1096      nbghostcells_y_s = nbghostcells 
     1097      nbghostcells_y_n = nbghostcells 
     1098 
     1099      IF (.not.agrif_root()) THEN 
     1100        IF (jperio == 1) THEN 
     1101          nbghostcells_x = 0 
     1102        ENDIF 
     1103        IF (.NOT.ln_bry_south) THEN 
     1104          nbghostcells_y_s = 0 
     1105        ENDIF 
     1106      ENDIF 
     1107 
     1108   END SUBROUTINE agrif_nemo_init 
     1109 
     1110   SUBROUTINE Agrif_detect( kg, ksizex ) 
    7401111      !!---------------------------------------------------------------------- 
    7411112      !!                      *** ROUTINE Agrif_detect *** 
    7421113      !!---------------------------------------------------------------------- 
    743    INTEGER, DIMENSION(2) :: ksizex 
    744    INTEGER, DIMENSION(ksizex(1),ksizex(2)) :: kg  
    745       !!---------------------------------------------------------------------- 
    746    ! 
    747    RETURN 
    748    ! 
    749 END SUBROUTINE Agrif_detect 
    750 SUBROUTINE agrif_before_regridding 
    751 END SUBROUTINE agrif_before_regridding 
    752  
    753 SUBROUTINE Agrif_InvLoc( indloc, nprocloc, i, indglob ) 
    754       !!---------------------------------------------------------------------- 
    755       !!                     *** ROUTINE Agrif_InvLoc *** 
    756       !!---------------------------------------------------------------------- 
    757    USE dom_oce 
    758    !! 
    759    IMPLICIT NONE 
    760    ! 
    761    INTEGER :: indglob, indloc, nprocloc, i 
    762       !!---------------------------------------------------------------------- 
    763    ! 
    764    SELECT CASE( i ) 
    765    CASE(1)   ;   indglob = indloc + nimppt(nprocloc+1) - 1 
    766    CASE(2)   ;   indglob = indloc + njmppt(nprocloc+1) - 1 
    767    CASE DEFAULT 
    768       indglob = indloc 
    769    END SELECT 
    770    ! 
    771 END SUBROUTINE Agrif_InvLoc 
    772  
    773 SUBROUTINE Agrif_get_proc_info( imin, imax, jmin, jmax ) 
     1114      INTEGER, DIMENSION(2) :: ksizex 
     1115      INTEGER, DIMENSION(ksizex(1),ksizex(2)) :: kg 
     1116      !!---------------------------------------------------------------------- 
     1117      ! 
     1118      RETURN 
     1119      ! 
     1120   END SUBROUTINE Agrif_detect 
     1121 
     1122   SUBROUTINE agrif_before_regridding 
     1123   END SUBROUTINE agrif_before_regridding 
     1124 
     1125# if defined  key_mpp_mpi 
     1126   SUBROUTINE Agrif_InvLoc( indloc, nprocloc, i, indglob ) 
     1127         !!---------------------------------------------------------------------- 
     1128         !!                     *** ROUTINE Agrif_InvLoc *** 
     1129         !!---------------------------------------------------------------------- 
     1130      USE dom_oce 
     1131      !! 
     1132      IMPLICIT NONE 
     1133      ! 
     1134      INTEGER :: indglob, indloc, nprocloc, i 
     1135         !!---------------------------------------------------------------------- 
     1136      ! 
     1137      SELECT CASE( i ) 
     1138      CASE(1)   ;   indglob = indloc + nimppt(nprocloc+1) - 1 
     1139      CASE(2)   ;   indglob = indloc + njmppt(nprocloc+1) - 1 
     1140      CASE DEFAULT 
     1141         indglob = indloc 
     1142      END SELECT 
     1143      ! 
     1144   END SUBROUTINE Agrif_InvLoc 
     1145 
     1146   SUBROUTINE Agrif_get_proc_info( imin, imax, jmin, jmax ) 
    7741147      !!---------------------------------------------------------------------- 
    7751148      !!                 *** ROUTINE Agrif_get_proc_info *** 
    7761149      !!---------------------------------------------------------------------- 
    777    USE par_oce 
    778    USE dom_oce  
    779    !! 
    780    IMPLICIT NONE 
    781    ! 
    782    INTEGER, INTENT(out) :: imin, imax 
    783    INTEGER, INTENT(out) :: jmin, jmax 
    784       !!---------------------------------------------------------------------- 
    785    ! 
    786    imin = nimppt(Agrif_Procrank+1)  ! ????? 
    787    jmin = njmppt(Agrif_Procrank+1)  ! ????? 
    788    imax = imin + jpi - 1 
    789    jmax = jmin + jpj - 1 
    790    !  
    791 END SUBROUTINE Agrif_get_proc_info 
    792  
    793 SUBROUTINE Agrif_estimate_parallel_cost(imin, imax,jmin, jmax, nbprocs, grid_cost) 
     1150      USE par_oce 
     1151      USE dom_oce 
     1152      !! 
     1153      IMPLICIT NONE 
     1154      ! 
     1155      INTEGER, INTENT(out) :: imin, imax 
     1156      INTEGER, INTENT(out) :: jmin, jmax 
     1157         !!---------------------------------------------------------------------- 
     1158      ! 
     1159      imin = nimppt(Agrif_Procrank+1)  ! ????? 
     1160      jmin = njmppt(Agrif_Procrank+1)  ! ????? 
     1161      imax = imin + jpi - 1 
     1162      jmax = jmin + jpj - 1 
     1163      ! 
     1164   END SUBROUTINE Agrif_get_proc_info 
     1165 
     1166   SUBROUTINE Agrif_estimate_parallel_cost(imin, imax,jmin, jmax, nbprocs, grid_cost) 
    7941167      !!---------------------------------------------------------------------- 
    7951168      !!                 *** ROUTINE Agrif_estimate_parallel_cost *** 
    7961169      !!---------------------------------------------------------------------- 
    797    USE par_oce 
    798    !! 
    799    IMPLICIT NONE 
    800    ! 
    801    INTEGER,  INTENT(in)  :: imin, imax 
    802    INTEGER,  INTENT(in)  :: jmin, jmax 
    803    INTEGER,  INTENT(in)  :: nbprocs 
    804    REAL(wp), INTENT(out) :: grid_cost 
    805       !!---------------------------------------------------------------------- 
    806    ! 
    807    grid_cost = REAL(imax-imin+1,wp)*REAL(jmax-jmin+1,wp) / REAL(nbprocs,wp) 
    808    ! 
    809 END SUBROUTINE Agrif_estimate_parallel_cost 
    810  
     1170      USE par_oce 
     1171      !! 
     1172      IMPLICIT NONE 
     1173      ! 
     1174      INTEGER,  INTENT(in)  :: imin, imax 
     1175      INTEGER,  INTENT(in)  :: jmin, jmax 
     1176      INTEGER,  INTENT(in)  :: nbprocs 
     1177      REAL(wp), INTENT(out) :: grid_cost 
     1178         !!---------------------------------------------------------------------- 
     1179      ! 
     1180      grid_cost = REAL(imax-imin+1,wp)*REAL(jmax-jmin+1,wp) / REAL(nbprocs,wp) 
     1181      ! 
     1182   END SUBROUTINE Agrif_estimate_parallel_cost 
     1183# endif 
    8111184#endif 
Note: See TracChangeset for help on using the changeset viewer.