From fd2a721246060e2731af86ac70c741cbbc1ed41a Mon Sep 17 00:00:00 2001
From: Philippe WAUTELET <philippe.wautelet@aero.obs-mip.fr>
Date: Mon, 24 Apr 2023 13:43:39 +0200
Subject: [PATCH] Philippe 24/04/2023: OpenACC: continuation of: merge ZSOLVER/
 version into MNH/

---
 src/LIB/SURCOUCHE/src/mode_exchange2_ll.f90 |  22 ++-
 src/MNH/advec_4th_order_aux.f90             |  65 ++------
 src/MNH/advection_uvw_cen.f90               |  23 ++-
 src/MNH/advecuvw_4th.f90                    |  34 +++-
 src/MNH/advecuvw_rk.f90                     | 147 +++++++++++-------
 src/MNH/advecuvw_weno_k.f90                 |  42 +++--
 src/MNH/contrav.f90                         | 163 +++++++++++---------
 src/MNH/get_halo.f90                        |  10 +-
 src/MNH/modeln.f90                          |  36 ++++-
 9 files changed, 340 insertions(+), 202 deletions(-)

diff --git a/src/LIB/SURCOUCHE/src/mode_exchange2_ll.f90 b/src/LIB/SURCOUCHE/src/mode_exchange2_ll.f90
index 93d537416..e124cb81e 100644
--- a/src/LIB/SURCOUCHE/src/mode_exchange2_ll.f90
+++ b/src/LIB/SURCOUCHE/src/mode_exchange2_ll.f90
@@ -120,6 +120,8 @@ implicit none
 !
   TYPE(HALO2LIST_ll), POINTER :: TZHALO2LIST
   INTEGER :: JJ ! loop counter
+
+  REAL , POINTER , CONTIGUOUS , DIMENSION(:,:) :: ZWEST,ZEAST,ZSOUTH,ZNORTH
 !
 !-------------------------------------------------------------------------------
 !
@@ -132,16 +134,26 @@ implicit none
 !
 !*       1.1   Allocate the current HALO2_ll
 !
-    ALLOCATE(TZHALO2LIST%HALO2) 
+    ALLOCATE(TZHALO2LIST%HALO2)
     ALLOCATE(TZHALO2LIST%HALO2%WEST(KDIMY, KDIMZ))
     ALLOCATE(TZHALO2LIST%HALO2%EAST(KDIMY, KDIMZ))
     ALLOCATE(TZHALO2LIST%HALO2%SOUTH(KDIMX, KDIMZ))
     ALLOCATE(TZHALO2LIST%HALO2%NORTH(KDIMX, KDIMZ))
+    ZWEST => TZHALO2LIST%HALO2%WEST
+    ZEAST => TZHALO2LIST%HALO2%EAST
+    ZSOUTH => TZHALO2LIST%HALO2%SOUTH
+    ZNORTH => TZHALO2LIST%HALO2%NORTH
+    !Remark: OpenACC data delete done in DEL_HALO2_ll
+    !$acc enter data create(ZWEST,ZEAST,ZSOUTH,ZNORTH)
+
+    !$acc kernels
+    ZWEST=0.
+    ZEAST=0.
+    ZSOUTH=0.
+    ZNORTH=0.
+    !$acc end kernels
+
     ALLOCATE(TZHALO2LIST%NEXT)
-    TZHALO2LIST%HALO2%WEST=0.
-    TZHALO2LIST%HALO2%EAST=0.
-    TZHALO2LIST%HALO2%SOUTH=0.
-    TZHALO2LIST%HALO2%NORTH=0.
 !
 !*       1.2   Go to the next HALO2_ll, or terminate the list
 !
diff --git a/src/MNH/advec_4th_order_aux.f90 b/src/MNH/advec_4th_order_aux.f90
index 238a9f843..9838c4762 100644
--- a/src/MNH/advec_4th_order_aux.f90
+++ b/src/MNH/advec_4th_order_aux.f90
@@ -159,13 +159,8 @@ INTEGER:: ILUOUT,IRESP   ! for prints
 ! JUAN ACC
 LOGICAL                                               :: GWEST , GEAST
 LOGICAL                                               :: GSOUTH , GNORTH
-#ifndef MNH_OPENACC
-REAL, DIMENSION(:,:), ALLOCATABLE :: ZHALO2_WEST,  ZHALO2_EAST
-REAL, DIMENSION(:,:), ALLOCATABLE :: ZHALO2_SOUTH, ZHALO2_NORTH
-#else
 REAL, DIMENSION(:,:), pointer, contiguous :: ZHALO2_WEST,  ZHALO2_EAST
 REAL, DIMENSION(:,:), pointer, contiguous :: ZHALO2_SOUTH, ZHALO2_NORTH
-#endif
 !
 
 !$acc data present( PMEANX, PMEANY, PFIELDT )
@@ -175,34 +170,11 @@ IF (MPPDB_INITIALIZED) THEN
   CALL MPPDB_CHECK(PFIELDT,"ADVEC_4TH_ORDER_ALGO beg:PFIELDT")
 END IF
 
-#ifndef MNH_OPENACC
-allocate( zhalo2_west ( size( pfieldt, 2 ), size( pfieldt, 3 ) ) )
-allocate( zhalo2_east ( size( pfieldt, 2 ), size( pfieldt, 3 ) ) )
-allocate( zhalo2_south( size( pfieldt, 2 ), size( pfieldt, 3 ) ) )
-allocate( zhalo2_north( size( pfieldt, 2 ), size( pfieldt, 3 ) ) )
-#else
-!Pin positions in the pools of MNH memory
-CALL MNH_MEM_POSITION_PIN( 'ADVEC_4TH_ORDER_ALGO' )
-
-CALL MNH_MEM_GET( zhalo2_west,  size( pfieldt, 2 ), size( pfieldt, 3 ) )
-CALL MNH_MEM_GET( zhalo2_east,  size( pfieldt, 2 ), size( pfieldt, 3 ) )
-CALL MNH_MEM_GET( zhalo2_south, size( pfieldt, 2 ), size( pfieldt, 3 ) )
-CALL MNH_MEM_GET( zhalo2_north, size( pfieldt, 2 ), size( pfieldt, 3 ) )
-
-!$acc data present ( zhalo2_west, zhalo2_east, zhalo2_south, zhalo2_north )
-#endif
-
 !-------------------------------------------------------------------------------
 !
 !*       0.3.     COMPUTES THE DOMAIN DIMENSIONS
 !                 ------------------------------
 !
-#ifdef MNH_OPENACC
-CALL INIT_ON_HOST_AND_DEVICE(ZHALO2_WEST,-1e99,'ADVEC_4TH_ORDER_ALGO::ZHALO2_WEST')
-CALL INIT_ON_HOST_AND_DEVICE(ZHALO2_EAST,-1e99,'ADVEC_4TH_ORDER_ALGO::ZHALO2_EAST')
-CALL INIT_ON_HOST_AND_DEVICE(ZHALO2_SOUTH,-1e99,'ADVEC_4TH_ORDER_ALGO::ZHALO2_SOUTH')
-CALL INIT_ON_HOST_AND_DEVICE(ZHALO2_NORTH,-1e99,'ADVEC_4TH_ORDER_ALGO::ZHALO2_NORTH')
-#endif
 !
 CALL GET_INDICE_ll(IIB,IJB,IIE,IJE)
 !
@@ -236,11 +208,10 @@ CASE ('CYCL')          ! In that case one must have HLBCX(1) == HLBCX(2)
 #ifdef MNH_OPENACC
 call Print_msg( NVERB_WARNING, 'GEN', 'ADVEC_4TH_ORDER_ALGO', 'OpenACC: HLBCX(1) AND CYCL not yet tested' )
 #endif
-ZHALO2_WEST(:,:) = TPHALO2%WEST(:,:)
-ZHALO2_EAST(:,:) = TPHALO2%EAST(:,:)
-!$acc update device (ZHALO2_WEST,ZHALO2_EAST)
+ZHALO2_WEST => TPHALO2%WEST
+ZHALO2_EAST => TPHALO2%EAST
 !
-!$acc kernels present(PMEANX)
+!$acc kernels present(PMEANX,ZHALO2_WEST,ZHALO2_EAST)
     IW=IIB+1
     IE=IIE
 !
@@ -281,11 +252,10 @@ ZHALO2_EAST(:,:) = TPHALO2%EAST(:,:)
 !
 CASE ('OPEN','WALL','NEST')
 !
-ZHALO2_WEST(:,:) = TPHALO2%WEST(:,:)
-ZHALO2_EAST(:,:) = TPHALO2%EAST(:,:)
-!$acc update device (ZHALO2_WEST,ZHALO2_EAST)
+ZHALO2_WEST => TPHALO2%WEST
+ZHALO2_EAST => TPHALO2%EAST
 !
-!$acc kernels present(PMEANX)
+!$acc kernels present(PMEANX,ZHALO2_WEST,ZHALO2_EAST)
   IF (GWEST) THEN
     IF(KGRID == 2) THEN
       IW=IIB+2          ! special case of C grid
@@ -364,11 +334,10 @@ IF ( .NOT. L2D ) THEN
 #ifdef MNH_OPENACC
 call Print_msg( NVERB_WARNING, 'GEN', 'ADVEC_4TH_ORDER_ALGO', 'OpenACC: HLBCX(2) AND CYCL not yet tested' )
 #endif
-ZHALO2_SOUTH(:,:) = TPHALO2%SOUTH(:,:) 
-ZHALO2_NORTH(:,:) = TPHALO2%NORTH(:,:)
-!$acc update device (ZHALO2_SOUTH,ZHALO2_NORTH)
+ZHALO2_SOUTH => TPHALO2%SOUTH 
+ZHALO2_NORTH => TPHALO2%NORTH
 !
-!$acc kernels present(PMEANY)
+!$acc kernels present(PMEANY,ZHALO2_SOUTH,ZHALO2_NORTH)
 !
 !
       IS=IJB+1
@@ -408,12 +377,11 @@ ZHALO2_NORTH(:,:) = TPHALO2%NORTH(:,:)
 !*       2.2    NON CYCLIC CASE IN THE Y DIRECTION
 !
   CASE ('OPEN','WALL','NEST')
+!     
+ZHALO2_SOUTH => TPHALO2%SOUTH
+ZHALO2_NORTH => TPHALO2%NORTH
 !
-ZHALO2_SOUTH(:,:) = TPHALO2%SOUTH(:,:)
-ZHALO2_NORTH(:,:) = TPHALO2%NORTH(:,:)
-!$acc update device (ZHALO2_SOUTH,ZHALO2_NORTH)
-!
-!$acc kernels present(PMEANY)
+!$acc kernels present(PMEANY,ZHALO2_SOUTH,ZHALO2_NORTH)
     IF (GSOUTH) THEN
       IF(KGRID == 3) THEN
         IS=IJB+2          ! special case of C grid
@@ -491,13 +459,6 @@ END IF
 
 !$acc end data
 
-#ifdef MNH_OPENACC
-!Release all memory allocated with MNH_MEM_GET calls since last call to MNH_MEM_POSITION_PIN
-CALL MNH_MEM_RELEASE( 'ADVEC_4TH_ORDER_ALGO' )
-#endif
-
-!$acc end data
-
 !-------------------------------------------------------------------------------
 !
 END SUBROUTINE ADVEC_4TH_ORDER_ALGO
diff --git a/src/MNH/advection_uvw_cen.f90 b/src/MNH/advection_uvw_cen.f90
index ddd55cc8b..340947ef1 100644
--- a/src/MNH/advection_uvw_cen.f90
+++ b/src/MNH/advection_uvw_cen.f90
@@ -16,7 +16,11 @@ INTERFACE
                            PUT, PVT, PWT,                          &
                            PRHODJ, PDXX, PDYY, PDZZ, PDZX, PDZY,   &
                            PRUS,PRVS, PRWS,                        &
+#ifndef MNH_OPENACC
                            TPHALO2MLIST                            )
+#else
+                           TPHALO2_UT,TPHALO2_VT,TPHALO2_WT )
+#endif
 !
 USE MODD_ARGSLIST_ll, ONLY : HALO2LIST_ll
 !
@@ -37,7 +41,11 @@ REAL, DIMENSION(:,:,:),   INTENT(INOUT) :: PRUS , PRVS  , PRWS
                                                   ! Sources terms 
 !
 ! halo lists for 4th order advection
+#ifndef MNH_OPENACC
 TYPE(HALO2LIST_ll), POINTER :: TPHALO2MLIST ! momentum variables
+#else
+TYPE(HALO2LIST_ll), POINTER :: TPHALO2_UT,TPHALO2_VT,TPHALO2_WT
+#endif
 !
 END SUBROUTINE ADVECTION_UVW_CEN
 !
@@ -53,7 +61,11 @@ END MODULE MODI_ADVECTION_UVW_CEN
                            PUT, PVT, PWT,                          &
                            PRHODJ, PDXX, PDYY, PDZZ, PDZX, PDZY,   &
                            PRUS,PRVS, PRWS,                        &
+#ifndef MNH_OPENACC
                            TPHALO2MLIST                            )
+#else
+                           TPHALO2_UT,TPHALO2_VT,TPHALO2_WT        )
+#endif
 !     ##########################################################################
 !
 !!****  *ADVECTION * - routine to call the specialized advection routines
@@ -141,7 +153,11 @@ REAL, DIMENSION(:,:,:),   INTENT(INOUT) :: PRUS , PRVS  , PRWS
                                                   ! Sources terms 
 !
 ! halo lists for 4th order advection
+#ifndef MNH_OPENACC
 TYPE(HALO2LIST_ll), POINTER :: TPHALO2MLIST ! momentum variables
+#else
+TYPE(HALO2LIST_ll), POINTER :: TPHALO2_UT,TPHALO2_VT,TPHALO2_WT
+#endif
 !
 !
 !*       0.2   declarations of local variables
@@ -378,7 +394,12 @@ IF (HUVW_ADV_SCHEME=='CEN2ND' ) THEN
 ELSEIF (HUVW_ADV_SCHEME=='CEN4TH') THEN
 ! 
    CALL ADVECUVW_4TH ( HLBCX, HLBCY, ZRUCT, ZRVCT, ZRWCT,            &
-                       PUT, PVT, PWT, ZRUS, ZRVS, ZRWS, TPHALO2MLIST )
+                       PUT, PVT, PWT, ZRUS, ZRVS, ZRWS,              &
+#ifndef MNH_OPENACC
+                       TPHALO2MLIST                                  )
+#else
+                       TPHALO2_UT,TPHALO2_VT,TPHALO2_WT              )
+#endif
 !
 END IF
 !
diff --git a/src/MNH/advecuvw_4th.f90 b/src/MNH/advecuvw_4th.f90
index f68752969..2cca13537 100644
--- a/src/MNH/advecuvw_4th.f90
+++ b/src/MNH/advecuvw_4th.f90
@@ -10,7 +10,12 @@
 INTERFACE
 !
       SUBROUTINE ADVECUVW_4TH ( HLBCX, HLBCY, PRUCT, PRVCT, PRWCT,           &
-                                PUT, PVT, PWT, PRUS, PRVS, PRWS, TPHALO2LIST )              
+                                PUT, PVT, PWT, PRUS, PRVS, PRWS,             &
+#ifndef MNH_OPENACC
+                                TPHALO2LIST )
+#else
+                                TPHALO2_UT,TPHALO2_VT,TPHALO2_WT )
+#endif
 !
 USE MODD_ARGSLIST_ll, ONLY : HALO2LIST_ll
 !
@@ -25,7 +30,11 @@ REAL, DIMENSION(:,:,:),   INTENT(IN) :: PUT, PVT, PWT        ! U,V,W at t
 !
 REAL, DIMENSION(:,:,:),   INTENT(INOUT) :: PRUS, PRVS, PRWS     ! Source terms
 !
+#ifndef MNH_OPENACC
 TYPE(HALO2LIST_ll), POINTER :: TPHALO2LIST ! list for diffusion
+#else
+TYPE(HALO2LIST_ll), POINTER :: TPHALO2_UT,TPHALO2_VT,TPHALO2_WT
+#endif
 !
 END SUBROUTINE ADVECUVW_4TH
 !
@@ -36,7 +45,12 @@ END MODULE MODI_ADVECUVW_4TH
 !
 !     ######################################################################
       SUBROUTINE ADVECUVW_4TH ( HLBCX, HLBCY, PRUCT, PRVCT, PRWCT,           &
-                                PUT, PVT, PWT, PRUS, PRVS, PRWS, TPHALO2LIST )              
+                                PUT, PVT, PWT, PRUS, PRVS, PRWS,             &
+#ifndef MNH_OPENACC
+                                TPHALO2LIST )
+#else
+                                TPHALO2_UT,TPHALO2_VT,TPHALO2_WT )
+#endif
 !     ######################################################################
 !
 !!****  *ADVECUVW_4TH * - routine to compute the 4th order centered
@@ -141,7 +155,11 @@ REAL, DIMENSION(:,:,:),   INTENT(IN) :: PUT, PVT, PWT        ! Variables at t
 !
 REAL, DIMENSION(:,:,:),   INTENT(INOUT) :: PRUS, PRVS, PRWS     ! Source terms
 !
+#ifndef MNH_OPENACC
 TYPE(HALO2LIST_ll), POINTER :: TPHALO2LIST ! list for diffusion
+#else
+TYPE(HALO2LIST_ll), POINTER :: TPHALO2_UT,TPHALO2_VT,TPHALO2_WT
+#endif
 !
 !*       0.2   Declarations of local variables :
 !
@@ -216,7 +234,11 @@ CALL MNH_MEM_GET( ZTEMP4, IIU, IJU, IKU )
 !
 IGRID = 2
 !!$IF(NHALO == 1) THEN
+#ifndef MNH_OPENACC
   TZHALO2LIST => TPHALO2LIST
+#else
+  TZHALO2LIST => TPHALO2_UT
+#endif
   CALL ADVEC_4TH_ORDER_ALGO(HLBCX, HLBCY, PUT, IGRID, ZMEANX, ZMEANY, &
                             TZHALO2LIST%HALO2 )
 !!$ELSE
@@ -269,7 +291,11 @@ PRUS(:,:,:) = PRUS(:,:,:) - ZTEMP4
 !
 IGRID = 3
 !!$IF(NHALO == 1) THEN
+#ifndef MNH_OPENACC
   TZHALO2LIST => TZHALO2LIST%NEXT
+#else
+  TZHALO2LIST => TPHALO2_VT
+#endif
   CALL ADVEC_4TH_ORDER_ALGO(HLBCX, HLBCY, PVT, IGRID, ZMEANX, ZMEANY, &
                             TZHALO2LIST%HALO2 )
 !!$ELSE
@@ -324,7 +350,11 @@ CALL MPPDB_CHECK(PRUCT,"ADVECUVW_4TH 02: PRUCT")
 IGRID = 4
 !
 !!$IF(NHALO == 1) THEN
+#ifndef MNH_OPENACC
   TZHALO2LIST => TZHALO2LIST%NEXT
+#else
+  TZHALO2LIST => TPHALO2_WT
+#endif
   CALL ADVEC_4TH_ORDER_ALGO(HLBCX, HLBCY, PWT, IGRID, ZMEANX, ZMEANY, &
                             TZHALO2LIST%HALO2 )
 !!$ELSE
diff --git a/src/MNH/advecuvw_rk.f90 b/src/MNH/advecuvw_rk.f90
index c3de8116b..c46759984 100644
--- a/src/MNH/advecuvw_rk.f90
+++ b/src/MNH/advecuvw_rk.f90
@@ -113,7 +113,8 @@ END MODULE MODI_ADVECUVW_RK
 !
 USE MODD_ARGSLIST_ll, ONLY: LIST_ll, HALO2LIST_ll
 USE MODD_CONF,        ONLY: NHALO
-USE MODD_IBM_PARAM_n, ONLY: LIBM, CIBM_ADV, XIBM_LS, XIBM_EPSI
+USE MODD_IBM_PARAM_n, ONLY: LIBM, XIBM_LS, XIBM_EPSI
+USE MODD_IBM_PARAM_n, ONLY: MODD_CIBM_ADV => CIBM_ADV
 USE MODD_PARAMETERS,  ONLY: JPVEXT
 USE MODD_SUB_MODEL_n, ONLY: XT_IBM_FORC
 !
@@ -196,6 +197,10 @@ REAL, DIMENSION(:),   POINTER, CONTIGUOUS :: ZBUTS! Butcher array coefficients
 !JUAN
 TYPE(LIST_ll), POINTER      :: TZFIELDMT_ll ! list of fields to exchange
 TYPE(HALO2LIST_ll), POINTER :: TZHALO2MT_ll ! momentum variables
+#ifdef MNH_OPENACC
+TYPE(HALO2LIST_ll), SAVE , POINTER :: TZHALO2_UT,TZHALO2_VT,TZHALO2_WT
+LOGICAL , SAVE :: GFIRST_CALL_ADVECUVW_RK = .TRUE.
+#endif
 INTEGER                     :: INBVAR
 INTEGER :: IIU, IJU, IKU ! array sizes
 !JUAN
@@ -213,7 +218,9 @@ TYPE(LIST_ll), POINTER      :: TZFIELDS4_ll ! list of fields to exchange
 LOGICAL :: GIBM !Intermediate variable used to work around a Cray compiler bug (CCE 13.0.0)
 REAL    :: ZIBM_EPSI !Intermediate variable used to work around a Cray compiler bug (CCE 13.0.0)
 REAL          :: ZTIME1,ZTIME2
-INTEGER :: JJI,JJJ,JJK
+CHARACTER(LEN=6) :: CIBM_ADV
+LOGICAL :: GIBM_FREEZE,GIBM_LOWORD,GIBM_FORCIN
+INTEGER :: JII,JJI,JKI
 !-------------------------------------------------------------------------------
 
 IF (MPPDB_INITIALIZED) THEN
@@ -237,9 +244,14 @@ END IF
 
 GIBM = LIBM
 ZIBM_EPSI = XIBM_EPSI
+CIBM_ADV = MODD_CIBM_ADV
+
+GIBM_FREEZE = ( GIBM .AND. CIBM_ADV=='FREEZE' )
+GIBM_LOWORD = ( GIBM .AND. CIBM_ADV=='LOWORD' )
+GIBM_FORCIN = ( GIBM .AND. CIBM_ADV=='FORCIN' )
 
 #ifdef MNH_OPENACC
-if ( LIBM ) call Print_msg( NVERB_FATAL, 'GEN', 'ADVECUVW_RK', 'OpenACC: LIBM=T not yet implemented' )
+if ( GIBM ) call Print_msg( NVERB_FATAL, 'GEN', 'ADVECUVW_RK', 'OpenACC: LIBM=T not yet implemented' )
 #endif
 !
 !*       0.     INITIALIZATION
@@ -408,7 +420,7 @@ IF ( GIBM ) THEN
 !$acc end kernels
 END IF
 !
-IF (GIBM .AND. CIBM_ADV=='FREEZE') THEN
+IF (GIBM_FREEZE) THEN
 !$acc kernels present(ZIBM)
   WHERE (XIBM_LS(:,:,:,2).GT.-ZIBM_EPSI) ZIBM(:,:,:,1) = 0.
   WHERE (XIBM_LS(:,:,:,3).GT.-ZIBM_EPSI) ZIBM(:,:,:,2) = 0.
@@ -416,19 +428,19 @@ IF (GIBM .AND. CIBM_ADV=='FREEZE') THEN
 !$acc end kernels
 ENDIF
 !
-!$acc kernels present_cr(PU,PV,PW,PRUS_ADV,PRVS_ADV,PRWS_ADV,ZUT,ZVT,ZWT)
-PRUS_ADV = 0.
-PRVS_ADV = 0.
-PRWS_ADV = 0.
+!$acc kernels present(PU,PV,PW,PRUS_ADV,PRVS_ADV,PRWS_ADV,ZUT,ZVT,ZWT)
+PRUS_ADV(:,:,:) = 0.
+PRVS_ADV(:,:,:) = 0.
+PRWS_ADV(:,:,:) = 0.
 !
 !-------------------------------------------------------------------------------
 !
 !*       2.     Wind guess before RK loop
 !        --------------------------------
 !
-ZUT = PU
-ZVT = PV
-ZWT = PW
+ZUT(:,:,:) = PU(:,:,:)
+ZVT(:,:,:) = PV(:,:,:)
+ZWT(:,:,:) = PW(:,:,:)
 !$acc end kernels
 !
 #ifndef MNH_OPENACC
@@ -446,7 +458,9 @@ CALL ADD3DFIELD_ll( TZFIELDMT_ll, ZUT, 'ADVECUVW_RK::ZUT' )
 CALL ADD3DFIELD_ll( TZFIELDMT_ll, ZVT, 'ADVECUVW_RK::ZVT' )
 CALL ADD3DFIELD_ll( TZFIELDMT_ll, ZWT, 'ADVECUVW_RK::ZWT' )
 INBVAR = 3
+#ifndef MNH_OPENACC
 CALL INIT_HALO2_ll(TZHALO2MT_ll,INBVAR,SIZE(PUT,1),SIZE(PUT,2),SIZE(PWT,3))
+#endif
 !
 !$acc kernels present_cr(ZRUS,ZRVS,ZRWS)
 ZRUS(:, :, :, : ) = 0.
@@ -475,19 +489,24 @@ RKLOOP: DO JS = 1, ISPL
   CALL UPDATE_HALO_ll(TZFIELDMT_ll,IINFO_ll)
   CALL UPDATE_HALO2_ll(TZFIELDMT_ll, TZHALO2MT_ll, IINFO_ll)
 #else
-! acc update self(ZUT,ZVT,ZWT)
-  CALL GET_HALO_D(ZUT,HNAME='ZUT')
-  CALL GET_HALO_D(ZVT,HNAME='ZVT')
-  CALL GET_HALO_D(ZWT,HNAME='ZWT')
-  CALL UPDATE_HALO2_ll(TZFIELDMT_ll, TZHALO2MT_ll, IINFO_ll)
-! acc update device(ZUT,ZVT,ZWT)       
+  IF (GFIRST_CALL_ADVECUVW_RK) THEN
+    GFIRST_CALL_ADVECUVW_RK = .FALSE.
+    NULLIFY(TZHALO2_UT,TZHALO2_VT,TZHALO2_WT)
+    CALL INIT_HALO2_ll(TZHALO2_UT,1,IIU,IJU,IKU)
+    CALL INIT_HALO2_ll(TZHALO2_VT,1,IIU,IJU,IKU)
+    CALL INIT_HALO2_ll(TZHALO2_WT,1,IIU,IJU,IKU)
+  END IF
+
+  CALL GET_HALO2_DF(ZUT,TZHALO2_UT,HNAME='ADVECUVW_RK::ZUT')
+  CALL GET_HALO2_DF(ZVT,TZHALO2_VT,HNAME='ADVECUVW_RK::ZVT')
+  CALL GET_HALO2_DF(ZWT,TZHALO2_WT,HNAME='ADVECUVW_RK::ZWT')
 #endif
 !
 !*       4.     Advection with WENO
 !        --------------------------
 !
-  IF (GIBM .AND. CIBM_ADV=='LOWORD') THEN
-!$acc kernels present_cr(ZIBM,ZRUS,ZRVS,ZRWS)
+  IF (GIBM_LOWORD) THEN
+!$acc kernels present(ZIBM)
     ZIBM(:,:,:,1)=ZRUS(:,:,:,JS)
     ZIBM(:,:,:,2)=ZRVS(:,:,:,JS)
     ZIBM(:,:,:,3)=ZRWS(:,:,:,JS)
@@ -498,42 +517,51 @@ RKLOOP: DO JS = 1, ISPL
     CALL ADVECUVW_WENO_K (HLBCX, HLBCY, KWENO_ORDER, ZUT, ZVT, ZWT,     &
                         PRUCT, PRVCT, PRWCT,                            &
                         ZRUS(:,:,:,JS), ZRVS(:,:,:,JS), ZRWS(:,:,:,JS), &
+#ifndef MNH_OPENACC
                         TZHALO2MT_ll                                    )
+#else
+                        TZHALO2_UT,TZHALO2_VT,TZHALO2_WT                )
+#endif
   ELSE IF ((HUVW_ADV_SCHEME=='CEN4TH') .AND. (HTEMP_SCHEME=='RKC4')) THEN
     CALL ADVECUVW_4TH (HLBCX, HLBCY, PRUCT, PRVCT, PRWCT,               &
                        ZUT, ZVT, ZWT,                                   &
                        ZRUS(:,:,:,JS), ZRVS(:,:,:,JS), ZRWS(:,:,:,JS),  &
-                       TZHALO2MT_ll )
+#ifndef MNH_OPENACC
+                       TZHALO2MT_ll                                     )
+#else
+                       TZHALO2_UT,TZHALO2_VT,TZHALO2_WT                 )
+#endif
   ENDIF
 !
-  IF (GIBM .AND. CIBM_ADV=='LOWORD') THEN
+  IF (GIBM_LOWORD) THEN
     IF (HUVW_ADV_SCHEME=='WENO_K') THEN
       CALL ADVECUVW_WENO_K (HLBCX, HLBCY,           3, ZUT, ZVT, ZWT,    &
                           PRUCT, PRVCT, PRWCT,                           &
                           ZIBM(:,:,:,1),  ZIBM(:,:,:,2),  ZIBM(:,:,:,3) ,&
-                        TZHALO2MT_ll                                    )
+#ifndef MNH_OPENACC
+                          TZHALO2MT_ll                                   )
+#else
+                          TZHALO2_UT,TZHALO2_VT,TZHALO2_WT               )
+#endif
     ENDIF
     IF (HUVW_ADV_SCHEME=='CEN4TH') THEN
-       !$acc update self(ZIBM)
        CALL ADVECUVW_2ND (ZUT, ZVT, ZWT, PRUCT, PRVCT, PRWCT,            &
                           ZIBM(:,:,:,1),  ZIBM(:,:,:,2),  ZIBM(:,:,:,3))
-       !$acc update device(ZIBM)
     ENDIF
-    !$acc kernels present_cr(ZIBM,ZRUS,ZRVS,ZRWS)
-    !$mnh_expand_where(JJI=1:IIU,JJJ=1:IJU,JJK=1:IKU)
-    WHERE(XIBM_LS(:,:,:,2).GT.-ZIBM_EPSI)
-       ZRUS(:,:,:,JS)=ZIBM(:,:,:,1)
-    END WHERE
-    WHERE(XIBM_LS(:,:,:,3).GT.-ZIBM_EPSI)
-       ZRVS(:,:,:,JS)=ZIBM(:,:,:,2)
-    END WHERE
-    WHERE(XIBM_LS(:,:,:,4).GT.-ZIBM_EPSI)
-       ZRWS(:,:,:,JS)=ZIBM(:,:,:,3)
-    END WHERE
-    !$mnh_end_expand_where()
-    ZIBM(:,:,:,:)=1.
-    !$acc end kernels
-  ENDIF
+! JE: Memory Leak with nvhpc22.2 with WHERE , even when GIBM FALSE /!\ -> replace DO CONCURRENT
+    DO CONCURRENT ( JII=1:IIU , JJI=1:IJU , JKI=1:IKU )
+      IF (XIBM_LS(JII,JJI,JKI,2).GT.-ZIBM_EPSI) THEN
+        ZRUS(JII,JJI,JKI,JS)=ZIBM(JII,JJI,JKI,1)
+      END IF
+      IF (XIBM_LS(JII,JJI,JKI,3).GT.-ZIBM_EPSI) THEN
+        ZRVS(JII,JJI,JKI,JS)=ZIBM(JII,JJI,JKI,2)
+      END IF
+      IF (XIBM_LS(JII,JJI,JKI,4).GT.-ZIBM_EPSI) THEN
+        ZRWS(JII,JJI,JKI,JS)=ZIBM(JII,JJI,JKI,3)
+      END IF
+    END DO
+    ZIBM(:,:,:,:) = 1.
+  END IF
 !
   write ( ynum, '( I3 )' ) js
 #ifndef MNH_OPENACC
@@ -550,13 +578,22 @@ RKLOOP: DO JS = 1, ISPL
 ! acc update device(ZRUS(:,:,:,JS),ZRVS(:,:,:,JS),ZRWS(:,:,:,JS))
 #endif
 !
-  IF (GIBM .AND. CIBM_ADV=='FREEZE') THEN
-    WHERE(XIBM_LS(:,:,:,2).GT.-ZIBM_EPSI) ZRUS(:,:,:,JS)=ZUT(:,:,:)*PMXM_RHODJ(:,:,:)/PTSTEP
-    WHERE(XIBM_LS(:,:,:,3).GT.-ZIBM_EPSI) ZRVS(:,:,:,JS)=ZVT(:,:,:)*PMYM_RHODJ(:,:,:)/PTSTEP
-    WHERE(XIBM_LS(:,:,:,4).GT.-ZIBM_EPSI) ZRWS(:,:,:,JS)=ZWT(:,:,:)*PMZM_RHODJ(:,:,:)/PTSTEP
-  ENDIF
+  IF (GIBM_FREEZE) THEN
+  ! JE: Memory Leak with nvhpc22.2 with WHERE , even when GIBM FALSE /!\ -> replace DO CONCURRENT
+    DO CONCURRENT ( JII=1:IIU , JJI=1:IJU , JKI=1:IKU )
+      IF (XIBM_LS(JII,JJI,JKI,2).GT.-ZIBM_EPSI) THEN
+        ZRUS(JII,JJI,JKI,JS)=ZUT(JII,JJI,JKI)*PMXM_RHODJ(JII,JJI,JKI)/PTSTEP
+      END IF
+      IF (XIBM_LS(JII,JJI,JKI,3).GT.-ZIBM_EPSI) THEN
+        ZRVS(JII,JJI,JKI,JS)=ZVT(JII,JJI,JKI)*PMYM_RHODJ(JII,JJI,JKI)/PTSTEP
+      END IF
+      IF (XIBM_LS(JII,JJI,JKI,4).GT.-ZIBM_EPSI) THEN
+        ZRWS(JII,JJI,JKI,JS)=ZWT(JII,JJI,JKI)*PMZM_RHODJ(JII,JJI,JKI)/PTSTEP
+      END IF
+    END DO
+  END IF
 
-  IF (GIBM .AND. CIBM_ADV=='FORCIN') THEN
+  IF (GIBM_FORCIN) THEN
     CALL SECOND_MNH(ZTIME1)
     CALL IBM_FORCING_ADV(ZRUS(:,:,:,JS),ZRVS(:,:,:,JS),ZRWS(:,:,:,JS))
     CALL SECOND_MNH(ZTIME2)
@@ -565,17 +602,21 @@ RKLOOP: DO JS = 1, ISPL
 !
 ! Guesses at the end of the RK loop
 !
-!$acc kernels present_cr(PRUS_ADV,PRVS_ADV,PRWS_ADV,ZBUTS) present_cr(ZRUS,ZRVS,ZRWS,ZIBM)
+
   IF ( .NOT. GIBM ) THEN
+!$acc kernels present(PRUS_ADV,PRVS_ADV,PRWS_ADV,ZBUTS) present(ZRUS,ZRVS,ZRWS)
     PRUS_ADV(:,:,:) = PRUS_ADV(:,:,:) + ZBUTS(JS) * ZRUS(:,:,:,JS)
     PRVS_ADV(:,:,:) = PRVS_ADV(:,:,:) + ZBUTS(JS) * ZRVS(:,:,:,JS)
     PRWS_ADV(:,:,:) = PRWS_ADV(:,:,:) + ZBUTS(JS) * ZRWS(:,:,:,JS)
-  ELSE
+!$acc end kernels
+ ELSE
+!$acc kernels present(PRUS_ADV,PRVS_ADV,PRWS_ADV,ZBUTS) present(ZRUS,ZRVS,ZRWS,ZIBM)
     PRUS_ADV(:,:,:) = PRUS_ADV(:,:,:) + ZBUTS(JS) * ZRUS(:,:,:,JS) * ZIBM(:,:,:,1)
     PRVS_ADV(:,:,:) = PRVS_ADV(:,:,:) + ZBUTS(JS) * ZRVS(:,:,:,JS) * ZIBM(:,:,:,2)
     PRWS_ADV(:,:,:) = PRWS_ADV(:,:,:) + ZBUTS(JS) * ZRWS(:,:,:,JS) * ZIBM(:,:,:,3)
-  END IF
 !$acc end kernels
+  END IF
+
 !
   IF ( JS < ISPL ) THEN
 !PW: note: 20211025: kernels split because performance problems if in 1 block with NVHPC 21.9
@@ -583,7 +624,7 @@ RKLOOP: DO JS = 1, ISPL
 ! !$acc & present(ZRUS,ZRVS,ZRWS,ZIBM) present(PRUS_OTHER,PRVS_OTHER,PRWS_OTHER) &
 ! !$acc & present(PMXM_RHODJ,PMYM_RHODJ,PMZM_RHODJ)
 !
-!$acc kernels present_cr( ZUT, ZVT, ZWT )
+!$acc kernels present( ZUT, ZVT, ZWT )
     ZUT(:,:,:) = PU(:,:,:)
     ZVT(:,:,:) = PV(:,:,:)
     ZWT(:,:,:) = PW(:,:,:)
@@ -594,7 +635,7 @@ RKLOOP: DO JS = 1, ISPL
 ! Intermediate guesses inside the RK loop
 !
       IF ( .NOT. GIBM ) THEN
-!$acc kernels present_cr(ZUT,ZVT,ZWT,ZRUS,ZRVS,ZRWS )
+!$acc kernels present( ZUT, ZVT, ZWT, ZRUS, ZRVS, ZRWS )
         ZUT(:,:,:) = ZUT(:,:,:) + ZBUT(JS,JI) *  PTSTEP *  &
           ( ZRUS(:,:,:,JI) + PRUS_OTHER(:,:,:) ) / PMXM_RHODJ(:,:,:)
         ZVT(:,:,:) = ZVT(:,:,:) + ZBUT(JS,JI) *  PTSTEP *  &
@@ -603,7 +644,7 @@ RKLOOP: DO JS = 1, ISPL
           ( ZRWS(:,:,:,JI) + PRWS_OTHER(:,:,:) ) / PMZM_RHODJ(:,:,:)
 !$acc end kernels
       ELSE
-!$acc kernels present_cr(ZUT,ZVT,ZWT,ZRUS,ZRVS,ZRWS,ZIBM )
+!$acc kernels present( ZUT, ZVT, ZWT, ZRUS, ZRVS, ZRWS, ZIBM )
         ZUT(:,:,:) = ZUT(:,:,:) + ZBUT(JS,JI) *  PTSTEP *  &
           ( ZRUS(:,:,:,JI) + PRUS_OTHER(:,:,:) ) / PMXM_RHODJ(:,:,:) * ZIBM(:,:,:,1)
         ZVT(:,:,:) = ZVT(:,:,:) + ZBUT(JS,JI) *  PTSTEP *  &
@@ -615,7 +656,7 @@ RKLOOP: DO JS = 1, ISPL
 !
     END DO
 ! !$acc end kernels
-!$acc update self(ZUT,ZVT,ZWT)
+! acc update self(ZUT,ZVT,ZWT)
   END IF
 !
 ! End of the RK loop
@@ -623,7 +664,9 @@ RKLOOP: DO JS = 1, ISPL
 !
 !
 CALL CLEANLIST_ll(TZFIELDMT_ll)
+#ifndef MNH_OPENACC
 CALL DEL_HALO2_ll(TZHALO2MT_ll)
+#endif
 !-------------------------------------------------------------------------------
 
 !$acc end data
diff --git a/src/MNH/advecuvw_weno_k.f90 b/src/MNH/advecuvw_weno_k.f90
index edd110e87..e438471e2 100644
--- a/src/MNH/advecuvw_weno_k.f90
+++ b/src/MNH/advecuvw_weno_k.f90
@@ -9,9 +9,13 @@
 !
 INTERFACE
 !
-      SUBROUTINE ADVECUVW_WENO_K(HLBCX, HLBCY, KWENO_ORDER, PUT, PVT, PWT,       &
-                             PRUCT, PRVCT, PRWCT, PRUS, PRVS, PRWS, TPHALO2LIST  &
-                             )
+      SUBROUTINE ADVECUVW_WENO_K( HLBCX, HLBCY, KWENO_ORDER, PUT, PVT, PWT, &
+                                  PRUCT, PRVCT, PRWCT, PRUS, PRVS, PRWS,    &
+#ifndef MNH_OPENACC
+                                  TPHALO2LIST )
+#else
+                                  TPHALO2_UT,TPHALO2_VT,TPHALO2_WT )
+#endif
 !
 USE MODD_ARGSLIST_ll, ONLY : HALO2LIST_ll
 !
@@ -28,7 +32,11 @@ REAL, DIMENSION(:,:,:),   INTENT(IN) :: PUT, PVT, PWT        ! U,V,W at t
 !
 REAL, DIMENSION(:,:,:),   INTENT(INOUT) :: PRUS, PRVS, PRWS     ! Source terms
 !
+#ifndef MNH_OPENACC
 TYPE(HALO2LIST_ll), POINTER :: TPHALO2LIST ! list for diffusion
+#else
+TYPE(HALO2LIST_ll), POINTER :: TPHALO2_UT,TPHALO2_VT,TPHALO2_WT
+#endif
 !
 !
 END SUBROUTINE ADVECUVW_WENO_K
@@ -37,11 +45,15 @@ END INTERFACE
 !
 END MODULE MODI_ADVECUVW_WENO_K
 !
-!     ##########################################################################
-      SUBROUTINE ADVECUVW_WENO_K(HLBCX, HLBCY, KWENO_ORDER, PUT, PVT, PWT,      &
-                             PRUCT, PRVCT, PRWCT, PRUS, PRVS, PRWS, TPHALO2LIST &
-                             )
-!     ##########################################################################
+!     #######################################################################
+      SUBROUTINE ADVECUVW_WENO_K( HLBCX, HLBCY, KWENO_ORDER, PUT, PVT, PWT, &
+                                  PRUCT, PRVCT, PRWCT, PRUS, PRVS, PRWS,    &
+#ifndef MNH_OPENACC
+                                  TPHALO2LIST )
+#else
+                                  TPHALO2_UT,TPHALO2_VT,TPHALO2_WT )
+#endif
+!     #######################################################################
 !
 !!    AUTHOR
 !!    ------
@@ -97,7 +109,11 @@ REAL, DIMENSION(:,:,:),   INTENT(IN) :: PUT, PVT, PWT     ! Variables at t
 !
 REAL, DIMENSION(:,:,:),   INTENT(INOUT) :: PRUS, PRVS, PRWS     ! Source terms
 !
+#ifndef MNH_OPENACC
 TYPE(HALO2LIST_ll), POINTER :: TPHALO2LIST ! list for diffusion
+#else
+TYPE(HALO2LIST_ll), POINTER :: TPHALO2_UT,TPHALO2_VT,TPHALO2_WT
+#endif
 !
 !*       0.2   Declarations of local variables :
 !
@@ -156,9 +172,15 @@ CALL INIT_ON_HOST_AND_DEVICE(ZWORK,2e90,'ADVECUVW_WENO_K::ZWORK')
 !------------------------- ADVECTION OF MOMENTUM ------------------------------
 !
 !
+#ifndef MNH_OPENACC
 TZHALO2_UT => TPHALO2LIST                   ! 1rst add3dfield in model_n
 TZHALO2_VT => TPHALO2LIST%NEXT              ! 2nd  add3dfield in model_n
 TZHALO2_WT => TPHALO2LIST%NEXT%NEXT         ! 3rst add3dfield in model_n
+#else
+TZHALO2_UT => TPHALO2_UT
+TZHALO2_VT => TPHALO2_VT
+TZHALO2_WT => TPHALO2_WT
+#endif
 !
 !      -------------------------------------------------------
 
@@ -238,7 +260,7 @@ CASE(1) ! WENO 1
   PRVS(:,:,:) = PRVS(:,:,:) - ZWORK(:,:,:)
 !$acc end kernels
 !
-  !PRVS = PRVS - DZF(UP_MZ(PVT,MYM(PRWCT)))
+  !PRVS = PRVS - DZF( UP_MZ(PVT,MYM(PRWCT)))
   CALL MYM_DEVICE(PRWCT,ZWORK)
   CALL UP_MZ_DEVICE(PVT,ZWORK,ZMEAN)
   CALL DZF_DEVICE( ZMEAN, ZWORK )
@@ -674,7 +696,7 @@ CASE(5) ! WENO 5
 !
 END SELECT
 !             ---------------------------------
-!$acc update self(PRUS,PRVS,PRWS)
+! acc update self(PRUS,PRVS,PRWS)
 
 !$acc end data
 
diff --git a/src/MNH/contrav.f90 b/src/MNH/contrav.f90
index 5918c3455..7629cbc99 100644
--- a/src/MNH/contrav.f90
+++ b/src/MNH/contrav.f90
@@ -234,20 +234,8 @@ END IF
 !              ------------------------------------
 !
 IF (LFLAT) THEN
-   !
-   PRWCT(:,:,:) = PRWT(:,:,:) / PDZZ(:,:,:)
-   !
-   IF (KADV_ORDER == 4 ) THEN
-      NULLIFY(TZFIELD_U)
-      NULLIFY(TZFIELD_V)
-      CALL ADD3DFIELD_ll( TZFIELD_U, PRUCT, 'CONTRAV::PRUCT' )
-      CALL ADD3DFIELD_ll( TZFIELD_V, PRVCT, 'CONTRAV::PRVCT' )
-      CALL UPDATE_HALO_ll(TZFIELD_U,IINFO_ll)
-      CALL UPDATE_HALO_ll(TZFIELD_V,IINFO_ll)
-   ENDIF
-   !
-   RETURN
-   !
+  PRWCT(:,:,:) = PRWT(:,:,:) / PDZZ(:,:,:)
+  RETURN
 END IF
 !
 !*       3.    Compute the vertical contravariant components (general case)
@@ -577,7 +565,8 @@ real                                :: ZTMP1, ZTMP2 ! Intermediate work variable
 REAL, DIMENSION(:,:),   POINTER, CONTIGUOUS :: ZU_EAST, ZV_NORTH, ZDZX_EAST, ZDZY_NORTH
 REAL, DIMENSION(:,:,:), POINTER, CONTIGUOUS :: Z1,Z2       ! Work arrays
 TYPE(LIST_ll),          POINTER     :: TZFIELD_U, TZFIELD_V, TZFIELD_DZX, TZFIELD_DZY
-TYPE(HALO2LIST_ll),     POINTER     :: TZHALO2_U, TZHALO2_V, TZHALO2_DZX, TZHALO2_DZY
+TYPE(HALO2LIST_ll), SAVE, POINTER  :: TZHALO2_U, TZHALO2_V, TZHALO2_DZX, TZHALO2_DZY
+LOGICAL , SAVE :: GFIRST_CALL_CONTRAV_DEVICE = .TRUE.
 !
 LOGICAL :: GWEST,GEAST,GSOUTH,GNORTH
 !
@@ -624,11 +613,31 @@ CALL GET_INDICE_ll( IIB,IJB,IIE,IJE)
 IKB=1+JPVEXT
 IKE=IKU - JPVEXT
 
+IF (KADV_ORDER == 4 ) THEN
+   IF( .NOT. LFLAT) THEN
+      IF ( GFIRST_CALL_CONTRAV_DEVICE ) THEN
+         GFIRST_CALL_CONTRAV_DEVICE = .FALSE.
+         NULLIFY(TZHALO2_U)
+         NULLIFY(TZHALO2_V)
+         NULLIFY(TZHALO2_DZX)
+         NULLIFY(TZHALO2_DZY)
+         CALL INIT_HALO2_ll(TZHALO2_U,1,IIU,IJU,IKU)
+         CALL INIT_HALO2_ll(TZHALO2_V,1,IIU,IJU,IKU)
+         CALL INIT_HALO2_ll(TZHALO2_DZX,1,IIU,IJU,IKU)
+         CALL INIT_HALO2_ll(TZHALO2_DZY,1,IIU,IJU,IKU)
+    END IF
+    ZU_EAST => TZHALO2_U%HALO2%EAST
+    ZDZX_EAST => TZHALO2_DZX%HALO2%EAST
+    ZV_NORTH => TZHALO2_V%HALO2%NORTH
+    ZDZY_NORTH => TZHALO2_DZY%HALO2%NORTH
+   END IF
+END IF
+
 CALL CONTRAV_DEVICE_DIM(&
          PRUT,PRVT,PRWT,&
          PDXX,PDYY,PDZZ,PDZX,PDZY,&
-         PRUCT,PRVCT,PRWCT,Z1,Z2)
-!         ZU_EAST,ZV_NORTH,ZDZX_EAST,ZDZY_NORTH)
+         PRUCT,PRVCT,PRWCT,Z1,Z2,&
+         ZU_EAST,ZV_NORTH,ZDZX_EAST,ZDZY_NORTH)
 
 !Release all memory allocated with MNH_MEM_GET calls since last call to MNH_MEM_POSITION_PIN
 CALL MNH_MEM_RELEASE( 'CONTRAV_DEVICE' )
@@ -645,8 +654,8 @@ CONTAINS
   SUBROUTINE CONTRAV_DEVICE_DIM(&
          PRUT,PRVT,PRWT,&
          PDXX,PDYY,PDZZ,PDZX,PDZY,&
-         PRUCT,PRVCT,PRWCT,Z1,Z2)
-!         ZU_EAST,ZV_NORTH,ZDZX_EAST,ZDZY_NORTH)
+         PRUCT,PRVCT,PRWCT,Z1,Z2,&
+         ZU_EAST,ZV_NORTH,ZDZX_EAST,ZDZY_NORTH)
 
     IMPLICIT NONE
 
@@ -654,12 +663,14 @@ CONTAINS
          PRUT,PRVT,PRWT,&
          PDXX,PDYY,PDZZ,PDZX,PDZY,&
          PRUCT,PRVCT,PRWCT,Z1,Z2
-!    REAL :: ZU_EAST(IJU,IKU),ZV_NORTH(IIU,IKU),ZDZX_EAST(IJU,IKU),ZDZY_NORTH(IIU,IKU)
+    REAL :: ZU_EAST(IJU,IKU),ZV_NORTH(IIU,IKU),ZDZX_EAST(IJU,IKU),ZDZY_NORTH(IIU,IKU)
+
 !
 !  Begin Compute
 !
 
-!$acc data present( PRUT, PRVT, PRWT, PDXX, PDYY, PDZZ, PDZX, PDZY, PRUCT, PRVCT, PRWCT, Z1, Z2 )
+!$acc data present( PRUT, PRVT, PRWT, PDXX, PDYY, PDZZ, PDZX, PDZY, PRUCT, PRVCT, PRWCT, Z1, Z2 ) &
+!$acc&     present( ZU_EAST, ZV_NORTH, ZDZX_EAST, ZDZY_NORTH )
 
 IF (GDATA_ON_DEVICE) THEN
 !PW TODO:remplacer (ailleurs aussi...) 1/PDXX... par PINV_PDXX (fait pour la turbulence...) cfr MNH/turb_hor_splt.f90
@@ -675,47 +686,52 @@ END IF
 !
 IF (KADV_ORDER == 4 ) THEN
  IF( .NOT. LFLAT) THEN
-  NULLIFY(TZFIELD_U)
-  NULLIFY(TZFIELD_V)
-  CALL ADD3DFIELD_ll( TZFIELD_U, PRUCT, 'CONTRAV::PRUCT' )
-  CALL ADD3DFIELD_ll( TZFIELD_V, PRVCT, 'CONTRAV::PRVCT' )
-  CALL UPDATE_HALO_ll(TZFIELD_U,IINFO_ll)
-  CALL UPDATE_HALO_ll(TZFIELD_V,IINFO_ll)
-!!$ IF( NHALO==1 ) THEN 
-  NULLIFY(TZFIELD_DZX)
-  NULLIFY(TZFIELD_DZY)
-  CALL ADD3DFIELD_ll( TZFIELD_DZX, PDZX, 'CONTRAV::PDZX' )
-  CALL ADD3DFIELD_ll( TZFIELD_DZY, PDZY, 'CONTRAV::PDZY' )
-  NULLIFY(TZHALO2_U)
-  NULLIFY(TZHALO2_V)
-  NULLIFY(TZHALO2_DZX)
-  NULLIFY(TZHALO2_DZY)
-  CALL INIT_HALO2_ll(TZHALO2_U,1,IIU,IJU,IKU)
-  CALL INIT_HALO2_ll(TZHALO2_V,1,IIU,IJU,IKU)
-  CALL INIT_HALO2_ll(TZHALO2_DZX,1,IIU,IJU,IKU)
-  CALL INIT_HALO2_ll(TZHALO2_DZY,1,IIU,IJU,IKU)
-  CALL UPDATE_HALO2_ll(TZFIELD_U, TZHALO2_U, IINFO_ll)
-  CALL UPDATE_HALO2_ll(TZFIELD_V, TZHALO2_V, IINFO_ll)
-  CALL UPDATE_HALO2_ll(TZFIELD_DZX, TZHALO2_DZX, IINFO_ll)
-  CALL UPDATE_HALO2_ll(TZFIELD_DZY, TZHALO2_DZY, IINFO_ll)
-!$acc update device(PRUCT,PRVCT)
-!!$ END IF
+!!$  NULLIFY(TZFIELD_U)
+!!$  NULLIFY(TZFIELD_V)
+!!$  CALL ADD3DFIELD_ll( TZFIELD_U, PRUCT, 'CONTRAV::PRUCT' )
+!!$  CALL ADD3DFIELD_ll( TZFIELD_V, PRVCT, 'CONTRAV::PRVCT' )
+!!$  CALL UPDATE_HALO_ll(TZFIELD_U,IINFO_ll)
+!!$  CALL UPDATE_HALO_ll(TZFIELD_V,IINFO_ll)
+!!$ !!$ IF( NHALO==1 ) THEN
+!!$  NULLIFY(TZFIELD_DZX)
+!!$  NULLIFY(TZFIELD_DZY)
+!!$  CALL ADD3DFIELD_ll( TZFIELD_DZX, PDZX, 'CONTRAV::PDZX' )
+!!$  CALL ADD3DFIELD_ll( TZFIELD_DZY, PDZY, 'CONTRAV::PDZY' )
+!!$    IF ( GFIRST_CALL_CONTRAV_DEVICE ) THEN
+!!$       GFIRST_CALL_CONTRAV_DEVICE = .FALSE.
+!!$       NULLIFY(TZHALO2_U)
+!!$       NULLIFY(TZHALO2_V)
+!!$       NULLIFY(TZHALO2_DZX)
+!!$       NULLIFY(TZHALO2_DZY)
+!!$       CALL INIT_HALO2_ll(TZHALO2_U,1,IIU,IJU,IKU)
+!!$       CALL INIT_HALO2_ll(TZHALO2_V,1,IIU,IJU,IKU)
+!!$       CALL INIT_HALO2_ll(TZHALO2_DZX,1,IIU,IJU,IKU)
+!!$       CALL INIT_HALO2_ll(TZHALO2_DZY,1,IIU,IJU,IKU)
+!!$    END IF
+!!$  CALL UPDATE_HALO2_ll(TZFIELD_U, TZHALO2_U, IINFO_ll)
+!!$  CALL UPDATE_HALO2_ll(TZFIELD_V, TZHALO2_V, IINFO_ll)
+!!$  CALL UPDATE_HALO2_ll(TZFIELD_DZX, TZHALO2_DZX, IINFO_ll)
+!!$  CALL UPDATE_HALO2_ll(TZFIELD_DZY, TZHALO2_DZY, IINFO_ll)
+  !
+  CALL GET_HALO2_DF(PRUCT,TZHALO2_U,'CONTRAV::PRUCT')
+  CALL GET_HALO2_DF(PRVCT,TZHALO2_V,'CONTRAV::PRVCT')
+  CALL GET_HALO2_DF(PDZX,TZHALO2_DZX,'CONTRAV::PDZX')
+  CALL GET_HALO2_DF(PDZY,TZHALO2_DZY,'CONTRAV::PDZY')
+!!$!$acc update device(PRUCT,PRVCT)
+!!$ !!$ END IF
 !
   !PW: necessary because pointers does not work with OpenACC (PGI 16.1)
-  !Pin positions in the pools of MNH memory
-  CALL MNH_MEM_POSITION_PIN( 'CONTRAV_DEVICE 2' )
-
-  CALL MNH_MEM_GET( zu_east,    IJU, IKU )
-  CALL MNH_MEM_GET( zv_north,   IIU, IKU )
-  CALL MNH_MEM_GET( zdzx_east,  IJU, IKU )
-  CALL MNH_MEM_GET( zdzy_north, IIU, IKU )
-  ZU_EAST(:,:)    = TZHALO2_U%HALO2%EAST
-  ZDZX_EAST(:,:)  = TZHALO2_DZX%HALO2%EAST
-  ZV_NORTH(:,:)   = TZHALO2_V%HALO2%NORTH
-  ZDZY_NORTH(:,:) = TZHALO2_DZY%HALO2%NORTH
-!$acc update device(ZU_EAST,ZV_NORTH,ZDZX_EAST,ZDZY_NORTH)
- END IF
-END IF
+!!$  ALLOCATE(ZU_EAST(IJU,IKU),ZV_NORTH(IIU,IKU),ZDZX_EAST(IJU,IKU),ZDZY_NORTH(IIU,IKU))
+!!$  !$acc enter data create( zu_east, zv_north, zdzx_east, zdzy_north )
+!!$  !$acc kernels
+!!$  ZU_EAST => TZHALO2_U%HALO2%EAST
+!!$  ZDZX_EAST => TZHALO2_DZX%HALO2%EAST
+!!$  ZV_NORTH => TZHALO2_V%HALO2%NORTH
+!!$  ZDZY_NORTH => TZHALO2_DZY%HALO2%NORTH
+!!$  !$acc end kernels
+!!$!$acc update device(ZU_EAST,ZV_NORTH,ZDZX_EAST,ZDZY_NORTH)
+ END IF ! NOT FLAT
+END IF ! KADV_ORDER == 4
 !
 !
 !*       2.    Compute the vertical contravariant components (flat case)
@@ -726,7 +742,7 @@ FLAT: IF (LFLAT) THEN
 !$acc kernels
     PRWCT(:,:,:) = PRWT(:,:,:) / PDZZ(:,:,:)
 !$acc end kernels
-!$acc update self(PRWCT)
+! acc update self(PRWCT)
   ELSE
     PRWCT(:,:,:) = PRWT(:,:,:) / PDZZ(:,:,:)
   END IF
@@ -914,20 +930,21 @@ END IF
 !$acc kernels
 PRWCT(:,:,1) = - PRWCT(:,:,3)     ! Mirror hypothesis
 !$acc end kernels
-!$acc update self(PRWCT)
+! acc update self(PRWCT)
 !
 IF (KADV_ORDER == 4 ) THEN
-  CALL MNH_MEM_RELEASE( 'CONTRAV_DEVICE 2' )
-  CALL CLEANLIST_ll(TZFIELD_U)
-  CALL CLEANLIST_ll(TZFIELD_V)
-!!$  IF (NHALO==1) THEN
-    CALL CLEANLIST_ll(TZFIELD_DZX)
-    CALL CLEANLIST_ll(TZFIELD_DZY)
-    CALL DEL_HALO2_ll(TZHALO2_U)
-    CALL DEL_HALO2_ll(TZHALO2_V)
-    CALL DEL_HALO2_ll(TZHALO2_DZX)
-    CALL DEL_HALO2_ll(TZHALO2_DZY)
-!!$  END IF
+!!$!$acc exit data delete( zu_east, zv_north, zdzx_east, zdzy_north )
+!!$  DEALLOCATE(ZU_EAST,ZV_NORTH,ZDZX_EAST,ZDZY_NORTH)
+!!$  CALL CLEANLIST_ll(TZFIELD_U)
+!!$  CALL CLEANLIST_ll(TZFIELD_V)
+!!$ !!$  IF (NHALO==1) THEN
+!!$    CALL CLEANLIST_ll(TZFIELD_DZX)
+!!$    CALL CLEANLIST_ll(TZFIELD_DZY)
+!!$    CALL DEL_HALO2_ll(TZHALO2_U)
+!!$    CALL DEL_HALO2_ll(TZHALO2_V)
+!!$    CALL DEL_HALO2_ll(TZHALO2_DZX)
+!!$    CALL DEL_HALO2_ll(TZHALO2_DZY)
+!!$ !!$  END IF
 END IF
 
 END IF FLAT
diff --git a/src/MNH/get_halo.f90 b/src/MNH/get_halo.f90
index 98277285d..bc3491ba4 100644
--- a/src/MNH/get_halo.f90
+++ b/src/MNH/get_halo.f90
@@ -2740,6 +2740,7 @@ END SUBROUTINE GET_HALO2_D
 !*       0.2   Declarations of local variables :
 !
   TYPE(HALO2LIST_ll), POINTER :: TZHALO2LIST
+  REAL , POINTER , CONTIGUOUS , DIMENSION(:,:) :: ZWEST,ZEAST,ZSOUTH,ZNORTH
 !
 !-------------------------------------------------------------------------------
 !
@@ -2749,7 +2750,14 @@ END SUBROUTINE GET_HALO2_D
 !
   DO WHILE(ASSOCIATED(TZHALO2LIST))
 !
-    TPHALO2LIST => TZHALO2LIST%NEXT
+   ZWEST => TZHALO2LIST%HALO2%WEST
+   ZEAST => TZHALO2LIST%HALO2%EAST
+   ZSOUTH => TZHALO2LIST%HALO2%SOUTH
+   ZNORTH => TZHALO2LIST%HALO2%NORTH
+    !Remark: OpenACC data create done in INIT_HALO2_ll
+   !$acc exit data delete(ZWEST,ZEAST,ZSOUTH,ZNORTH)
+
+   TPHALO2LIST => TZHALO2LIST%NEXT
     DEALLOCATE(TZHALO2LIST%HALO2%WEST)
     DEALLOCATE(TZHALO2LIST%HALO2%EAST)
     DEALLOCATE(TZHALO2LIST%HALO2%SOUTH)
diff --git a/src/MNH/modeln.f90 b/src/MNH/modeln.f90
index 37e3d4061..56b9978d1 100644
--- a/src/MNH/modeln.f90
+++ b/src/MNH/modeln.f90
@@ -563,6 +563,10 @@ REAL, DIMENSION(:,:,:), ALLOCATABLE :: ZJ
 !
 TYPE(LIST_ll), POINTER :: TZFIELDC_ll   ! list of fields to exchange
 TYPE(HALO2LIST_ll), POINTER :: TZHALO2C_ll   ! list of fields to exchange
+#ifdef MNH_OPENACC
+TYPE(HALO2LIST_ll), SAVE , POINTER :: TZHALO2_UT,TZHALO2_VT,TZHALO2_WT
+LOGICAL , SAVE :: GFIRST_CALL_MODELN = .TRUE.
+#endif
 LOGICAL :: GCLD                     ! conditionnal call for dust wet deposition
 LOGICAL :: GCLOUD_ONLY              ! conditionnal radiation computations for
                                 !      the only cloudy columns
@@ -1716,9 +1720,23 @@ IF ((CUVW_ADV_SCHEME(1:3)=='CEN') .AND. (CTEMP_SCHEME == 'LEFR')) THEN
       CALL ADD3DFIELD_ll( TZFIELDC_ll, XUT, 'MODEL_n::XUT' )
       CALL ADD3DFIELD_ll( TZFIELDC_ll, XVT, 'MODEL_n::XVT' )
       CALL ADD3DFIELD_ll( TZFIELDC_ll, XWT, 'MODEL_n::XWT' )
+#ifndef MNH_OPENACC
       CALL INIT_HALO2_ll(TZHALO2C_ll,3,IIU,IJU,IKU)
       CALL UPDATE_HALO_ll(TZFIELDC_ll,IINFO_ll)
       CALL UPDATE_HALO2_ll(TZFIELDC_ll, TZHALO2C_ll, IINFO_ll)
+#else
+  IF (GFIRST_CALL_MODELN) THEN
+    GFIRST_CALL_MODELN = .FALSE.
+    NULLIFY(TZHALO2_UT,TZHALO2_VT,TZHALO2_WT)
+    CALL INIT_HALO2_ll(TZHALO2_UT,1,IIU,IJU,IKU)
+    CALL INIT_HALO2_ll(TZHALO2_VT,1,IIU,IJU,IKU)
+    CALL INIT_HALO2_ll(TZHALO2_WT,1,IIU,IJU,IKU)
+  END IF
+
+  CALL GET_HALO2_DF(XUT,TZHALO2_UT,HNAME='XUT')
+  CALL GET_HALO2_DF(XVT,TZHALO2_VT,HNAME='XVT')
+  CALL GET_HALO2_DF(XWT,TZHALO2_WT,HNAME='XWT')
+#endif
 !$acc update device(XUT, XVT, XWT)
   END IF
 !$acc data copyin(XUM, XVM, XWM)    &
@@ -1730,13 +1748,19 @@ IF ((CUVW_ADV_SCHEME(1:3)=='CEN') .AND. (CTEMP_SCHEME == 'LEFR')) THEN
                            XUT, XVT, XWT,                          &
                            XRHODJ, XDXX, XDYY, XDZZ, XDZX, XDZY,   &
                            XRUS,XRVS, XRWS,                        &
+#ifndef MNH_OPENACC
                            TZHALO2C_ll                             )
+#else
+                           TZHALO2_UT,TZHALO2_VT,TZHALO2_WT        )
+#endif
 !$acc end data
   IF (CUVW_ADV_SCHEME=='CEN4TH') THEN
     CALL CLEANLIST_ll(TZFIELDC_ll)
     NULLIFY(TZFIELDC_ll)
+#ifndef MNH_OPENACC
     CALL  DEL_HALO2_ll(TZHALO2C_ll)
     NULLIFY(TZHALO2C_ll)
+#endif
   END IF
 ELSE
 
@@ -1825,9 +1849,9 @@ ZPABST(:,:,:) = XPABST(:,:,:)
 IF(.NOT. L1D) THEN
 !
 !$acc kernels
-  XRUS_PRES = XRUS
-  XRVS_PRES = XRVS
-  XRWS_PRES = XRWS
+  XRUS_PRES(:,:,:) = XRUS(:,:,:)
+  XRVS_PRES(:,:,:) = XRVS(:,:,:)
+  XRWS_PRES(:,:,:) = XRWS(:,:,:)
 !$acc end kernels
 !
 !$acc data present( XRHODJ, XDXX, XDYY, XDZZ, XDZX, XDZY )                    &
@@ -1853,9 +1877,9 @@ IF(.NOT. L1D) THEN
 !$acc update self( XRUS, XRVS, XRWS, XPABST )
 !
 !$acc kernels
-  XRUS_PRES = XRUS - XRUS_PRES + ZRUS
-  XRVS_PRES = XRVS - XRVS_PRES + ZRVS
-  XRWS_PRES = XRWS - XRWS_PRES + ZRWS
+  XRUS_PRES(:,:,:) = XRUS(:,:,:) - XRUS_PRES(:,:,:) + ZRUS(:,:,:)
+  XRVS_PRES(:,:,:) = XRVS(:,:,:) - XRVS_PRES(:,:,:) + ZRVS(:,:,:)
+  XRWS_PRES(:,:,:) = XRWS(:,:,:) - XRWS_PRES(:,:,:) + ZRWS(:,:,:)
 !$acc end kernels
 !$acc update self( XRUS_PRES, XRVS_PRES, XRWS_PRES )
   CALL MPPDB_CHECK3DM("after pressurez:XRU/V/WS",PRECISION,XRUS,XRVS,XRWS)
-- 
GitLab