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 14622 – NEMO

Changeset 14622


Ignore:
Timestamp:
2021-03-21T19:36:08+01:00 (3 years ago)
Author:
ldebreu
Message:

AGFdomcfg: 1) define some AGRIF variables in real8 2) remove the wrong indices translation in modarrays #2638

Location:
vendors/AGRIF/dev_r14608_AGRIF_domcfg/AGRIF_FILES
Files:
10 edited

Legend:

Unmodified
Added
Removed
  • vendors/AGRIF/dev_r14608_AGRIF_domcfg/AGRIF_FILES/modarrays.F90

    r14107 r14622  
    195195 
    196196    if (present(pvariable)) then 
     197      lower=-1 
     198      upper=-1 
    197199     if (variable%root_var%interptab(index) == 'N') then 
    198200        lower = pvariable%lb(index) 
    199201        upper = pvariable%ub(index) 
     202      else 
     203       lower = variable % lb(index) 
     204       upper = variable % ub(index) 
    200205      endif 
    201206    else 
     
    703708    integer, dimension(6), intent(out)          :: lb_child     !< Lower bound on the child grid 
    704709    integer, dimension(6), intent(out)          :: lb_parent    !< Lower bound on the parent grid 
    705     real, dimension(6),    intent(out)          :: s_child      !< Child  grid position (s_root = 0) 
    706     real, dimension(6),    intent(out)          :: s_parent     !< Parent grid position (s_root = 0) 
    707     real, dimension(6),    intent(out)          :: ds_child     !< Child  grid dx (ds_root = 1) 
    708     real, dimension(6),    intent(out)          :: ds_parent    !< Parent grid dx (ds_root = 1) 
     710    real(kind=8), dimension(6),    intent(out)  :: s_child      !< Child  grid position (s_root = 0) 
     711    real(kind=8), dimension(6),    intent(out)  :: s_parent     !< Parent grid position (s_root = 0) 
     712    real(kind=8), dimension(6),    intent(out)  :: ds_child     !< Child  grid dx (ds_root = 1) 
     713    real(kind=8), dimension(6),    intent(out)  :: ds_parent    !< Parent grid dx (ds_root = 1) 
    709714    integer,               intent(out)          :: nbdim        !< Number of dimensions 
    710715    logical,               intent(in)           :: interp       !< .true. if preprocess for interpolation, \n 
     
    739744            ds_parent(n) = Agrif_Parent_Gr % Agrif_dx(1) 
    740745            ! Take into account potential difference of first points 
    741             s_parent(n) = s_parent(n) + (lb_parent(n)-lb_child(n))*ds_parent(n) 
     746           ! s_parent(n) = s_parent(n) + (lb_parent(n)-lb_child(n))*ds_parent(n) 
    742747! 
    743748            if ( root_var % posvar(n) == 1 ) then 
     
    745750            else 
    746751                ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(1) - 1 
    747                 s_child(n)  = s_child(n)  + 0.5*ds_child(n) 
    748                 s_parent(n) = s_parent(n) + 0.5*ds_parent(n) 
     752                s_child(n)  = s_child(n)  + 0.5d0*ds_child(n) 
     753                s_parent(n) = s_parent(n) + 0.5d0*ds_parent(n) 
    749754            endif 
    750755! 
     
    759764            ds_parent(n) = Agrif_Parent_Gr % Agrif_dx(2) 
    760765            ! Take into account potential difference of first points 
    761             s_parent(n) = s_parent(n) + (lb_parent(n)-lb_child(n))*ds_parent(n) 
     766           ! s_parent(n) = s_parent(n) + (lb_parent(n)-lb_child(n))*ds_parent(n) 
    762767! 
    763768            if (root_var % posvar(n)==1) then 
     
    765770            else 
    766771                ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(2) - 1 
    767                 s_child(n)  = s_child(n)  + 0.5*ds_child(n) 
    768                 s_parent(n) = s_parent(n) + 0.5*ds_parent(n) 
     772                s_child(n)  = s_child(n)  + 0.5d0*ds_child(n) 
     773                s_parent(n) = s_parent(n) + 0.5d0*ds_parent(n) 
    769774            endif 
    770775! 
     
    779784            ds_parent(n) = Agrif_Parent_Gr % Agrif_dx(3) 
    780785            ! Take into account potential difference of first points 
    781             s_parent(n) = s_parent(n) + (lb_parent(n)-lb_child(n))*ds_parent(n) 
     786           ! s_parent(n) = s_parent(n) + (lb_parent(n)-lb_child(n))*ds_parent(n) 
    782787! 
    783788            if (root_var % posvar(n)==1) then 
     
    805810!           No interpolation but only a copy of the values of the grid variable 
    806811            lb_parent(n) = lb_child(n) 
    807             s_child(n)   = 0. 
    808             s_parent(n)  = 0. 
    809             ds_child(n)  = 1. 
    810             ds_parent(n) = 1. 
     812            s_child(n)   = 0.d0 
     813            s_parent(n)  = 0.d0 
     814            ds_child(n)  = 1.d0 
     815            ds_parent(n) = 1.d0 
    811816! 
    812817        end select 
  • vendors/AGRIF/dev_r14608_AGRIF_domcfg/AGRIF_FILES/modbc.F90

    r14107 r14622  
    6161    integer, dimension(6)  :: loctab_child      ! Indicates if the child grid has a common border 
    6262                                                !    with the root grid 
    63     real, dimension(6)     :: s_child, s_parent   ! Positions of the parent and child grids 
    64     real, dimension(6)     :: ds_child, ds_parent ! Space steps of the parent and child grids 
     63    real(kind=8), dimension(6)     :: s_child, s_parent   ! Positions of the parent and child grids 
     64    real(kind=8), dimension(6)     :: ds_child, ds_parent ! Space steps of the parent and child grids 
    6565! 
    6666    call PreProcessToInterpOrUpdate( parent,   child,       & 
     
    145145    INTEGER, DIMENSION(nbdim)   :: posvartab_Child      !< Position of the grid variable (1 or 2) 
    146146    INTEGER, DIMENSION(nbdim)   :: loctab_Child         !< Indicates if the child grid has a common border with the root grid 
    147     REAL   , DIMENSION(nbdim)   :: s_Child,  s_Parent   !< Positions of the parent and child grids 
    148     REAL   , DIMENSION(nbdim)   :: ds_Child, ds_Parent  !< Space steps of the parent and child grids 
     147    REAL(kind=8)   , DIMENSION(nbdim)   :: s_Child,  s_Parent   !< Positions of the parent and child grids 
     148    REAL(kind=8)   , DIMENSION(nbdim)   :: ds_Child, ds_Parent  !< Space steps of the parent and child grids 
    149149    INTEGER                             :: nbdim        !< Number of dimensions of the grid variable 
    150150    procedure()                         :: procname     !< Data recovery procedure 
  • vendors/AGRIF/dev_r14608_AGRIF_domcfg/AGRIF_FILES/modbcfunction.F90

    r14107 r14622  
    3030    use Agrif_Update 
    3131    use Agrif_Save 
     32    use Agrif_Arrays 
    3233! 
    3334    implicit none 
     
    431432end subroutine Agrif_Bc_variable 
    432433!=================================================================================================== 
     434!=================================================================================================== 
     435!  subroutine Agrif_Find_Nearest : find nearest point on the parent grid 
     436!--------------------------------------------------------------------------------------------------- 
     437subroutine Agrif_Find_Nearest ( tabvarsindic, fineloc, parentloc ) 
     438!--------------------------------------------------------------------------------------------------- 
     439    integer,        intent(in) :: tabvarsindic     !< indice of the variable in tabvars 
     440    integer,dimension(:), intent(in)    :: fineloc 
     441    integer,dimension(:,:), intent(out) :: parentloc 
     442 
     443! 
     444 
     445    integer :: indic 
     446 
     447    type(Agrif_Variable), pointer :: parent_var 
     448    type(Agrif_Variable), pointer :: child_var 
     449    integer :: i 
     450    integer, dimension(6)           :: nb_child     !< Number of cells on the child grid 
     451    integer, dimension(6)           :: ub_child     !< Upper bound on the child grid 
     452    integer, dimension(6)           :: lb_child     !< Lower bound on the child grid 
     453    integer, dimension(6)           :: lb_parent    !< Lower bound on the parent grid 
     454    real, dimension(6)              :: s_child      !< Child  grid position (s_root = 0) 
     455    real, dimension(6)              :: s_parent     !< Parent grid position (s_root = 0) 
     456    real, dimension(6)              :: ds_child     !< Child  grid dx (ds_root = 1) 
     457    real, dimension(6)              :: ds_parent    !< Parent grid dx (ds_root = 1) 
     458    integer                         :: nbdim        !< Number of dimensions 
     459    real, dimension(6) :: xfineloc 
     460! 
     461    if ( Agrif_Curgrid%level <= 0 ) return 
     462! 
     463! 
     464! 
     465        child_var  => Agrif_Search_Variable(Agrif_Curgrid,tabvarsindic) 
     466        parent_var => child_var % parent_var 
     467! 
     468 
     469    call PreProcessToInterpOrUpdate( parent_var,   child_var,       & 
     470                                     nb_child, ub_child,    & 
     471                                     lb_child, lb_parent,   & 
     472                                      s_child,  s_parent,   & 
     473                                     ds_child, ds_parent, nbdim, interp=.true.) 
     474 
     475    do i=1,nbdim 
     476      xfineloc(i) = s_Child(i) + (fineloc(i) - lb_Child(i)) * ds_Child(i) 
     477 
     478      parentloc(1,i) = lb_parent(i) + agrif_int((xfineloc(i)-s_Parent(i))/ds_Parent(i)) 
     479      parentloc(2,i) = lb_parent(i) + agrif_ceiling((xfineloc(i)-s_Parent(i))/ds_Parent(i)) 
     480    enddo  
     481! 
     482!--------------------------------------------------------------------------------------------------- 
     483end subroutine Agrif_Find_Nearest 
     484!=================================================================================================== 
    433485! 
    434486!=================================================================================================== 
  • vendors/AGRIF/dev_r14608_AGRIF_domcfg/AGRIF_FILES/modcluster.F90

    r5656 r14622  
    5454    TYPE(Agrif_LRectangle), pointer  :: parcours 
    5555    TYPE(Agrif_Grid)      , pointer  :: newgrid 
    56     REAL                             :: g_eps 
     56    REAL(kind=8)                     :: g_eps 
    5757    INTEGER                          :: i 
    5858! 
     
    131131    TYPE(Agrif_PGrid), pointer  :: parcours 
    132132! 
    133     REAL                  :: g_eps, newgrid_eps, eps 
    134     REAL   , DIMENSION(3) :: newmin, newmax 
    135     REAL   , DIMENSION(3) :: gmin, gmax 
    136     REAL   , DIMENSION(3) :: xmin 
     133    REAL(kind=8)                  :: g_eps, newgrid_eps, eps 
     134    REAL(kind=8)   , DIMENSION(3) :: newmin, newmax 
     135    REAL(kind=8)   , DIMENSION(3) :: gmin, gmax 
     136    REAL(kind=8)   , DIMENSION(3) :: xmin 
    137137    INTEGER, DIMENSION(3) :: igmin, inewmin 
    138138    INTEGER, DIMENSION(3) :: inewmax 
  • vendors/AGRIF/dev_r14608_AGRIF_domcfg/AGRIF_FILES/modgrids.F90

    r14107 r14622  
    4444    type(Agrif_Variable_i), dimension(:), allocatable :: tabvars_i  !< List of integer   grid variables 
    4545! 
    46     real   , dimension(3)              :: Agrif_x   !< global x, y and z position 
    47     real   , dimension(3)              :: Agrif_dx  !< global space step in the x, y and z direction 
     46    real(kind=8), dimension(3)         :: Agrif_x   !< global x, y and z position 
     47    real(kind=8)   , dimension(3)      :: Agrif_dx  !< global space step in the x, y and z direction 
    4848    real   , dimension(3)              :: Agrif_dt  !< global time  step in the x, y and z direction 
    4949    integer, dimension(3)              :: nb = 1    !< number of cells in the x, y and z direction 
  • vendors/AGRIF/dev_r14608_AGRIF_domcfg/AGRIF_FILES/modinterp.F90

    r14107 r14622  
    6666    integer, dimension(6) :: ub_child 
    6767    integer, dimension(6) :: lb_parent 
    68     real   , dimension(6) :: s_child,   s_parent 
    69     real   , dimension(6) :: ds_child, ds_parent 
     68    real(kind=8)   , dimension(6) :: s_child,   s_parent 
     69    real(kind=8)   , dimension(6) :: ds_child, ds_parent 
    7070    integer, dimension(child % root_var % nbdim,2,2) :: childarray 
    7171! 
     
    115115    INTEGER, DIMENSION(nbdim), INTENT(in)   :: pttab_Parent !< Index of the first point inside the domain 
    116116                                                            !<    for the parent grid variable 
    117     REAL,    DIMENSION(nbdim), INTENT(in)   :: s_Child,s_Parent   !< Positions of the parent and child grids 
    118     REAL,    DIMENSION(nbdim), INTENT(in)   :: ds_Child,ds_Parent !< Space steps of the parent and child grids 
     117    REAL(kind=8),    DIMENSION(nbdim), INTENT(in)   :: s_Child,s_Parent   !< Positions of the parent and child grids 
     118    REAL(kind=8),    DIMENSION(nbdim), INTENT(in)   :: ds_Child,ds_Parent !< Space steps of the parent and child grids 
    119119    TYPE(Agrif_Variable),      pointer      :: restore            !< Indicates points where interpolation 
    120120    LOGICAL,                   INTENT(in)   :: torestore          !< Indicates if the array restore is used 
     
    139139#endif 
    140140    LOGICAL, DIMENSION(nbdim)     :: noraftab 
    141     REAL   , DIMENSION(nbdim)     :: s_Child_temp,s_Parent_temp 
     141    REAL(kind=8)   , DIMENSION(nbdim)     :: s_Child_temp,s_Parent_temp 
    142142    INTEGER, DIMENSION(nbdim)     :: lowerbound, upperbound, coords 
    143143    INTEGER, DIMENSION(nbdim,2,2), INTENT(OUT) :: childarray 
     
    608608    Agrif_CurChildgrid=>Agrif_Curgrid 
    609609    call Agrif_ChildGrid_to_ParentGrid() 
     610 
     611 
    610612    do i=1,nb_chunks 
    611613    if (agrif_debug_interp) then 
    612614    print *,'PROCNAME POUR CHUCNK ',i 
    613615    endif 
    614  
     616     
    615617    if (member_chuncks(i)) then 
    616618        select case (nbdim) 
     
    13031305    INTEGER,                   intent(in)  :: nbdim 
    13041306    INTEGER, DIMENSION(nbdim), intent(out) :: indmin, indmax 
    1305     REAL,    DIMENSION(nbdim), intent(out) :: s_Parent_temp, s_child_temp 
    1306     REAL,    DIMENSION(nbdim), intent(in)  :: s_Child, ds_child 
    1307     REAL,    DIMENSION(nbdim), intent(in)  :: s_Parent,ds_Parent 
     1307    REAL(kind=8),    DIMENSION(nbdim), intent(out) :: s_Parent_temp, s_child_temp 
     1308    REAL(kind=8),    DIMENSION(nbdim), intent(in)  :: s_Child, ds_child 
     1309    REAL(kind=8),    DIMENSION(nbdim), intent(in)  :: s_Parent,ds_Parent 
    13081310    INTEGER, DIMENSION(nbdim), intent(in)  :: pttruetab, cetruetab 
    13091311    INTEGER, DIMENSION(nbdim), intent(in)  :: pttab_Child, pttab_Parent 
     
    13121314! 
    13131315    INTEGER :: i 
    1314     REAL,DIMENSION(nbdim) :: dim_newmin, dim_newmax 
     1316    REAL(kind=8),DIMENSION(nbdim) :: dim_newmin, dim_newmax 
    13151317! 
    13161318    dim_newmin = s_Child + (pttruetab - pttab_Child) * ds_Child 
     
    13801382    integer,            intent(in)  :: indmin, indmax 
    13811383    integer,            intent(in)  :: pttab_child, petab_child 
    1382     real,               intent(in)  :: s_child, s_parent 
    1383     real,               intent(in)  :: ds_child, ds_parent 
     1384    real(kind=8),               intent(in)  :: s_child, s_parent 
     1385    real(kind=8),               intent(in)  :: ds_child, ds_parent 
    13841386    real, dimension(            & 
    13851387        indmin:indmax           & 
     
    14151417    integer, dimension(2),              intent(in)  :: indmin, indmax 
    14161418    integer, dimension(2),              intent(in)  :: pttab_child, petab_child 
    1417     real,    dimension(2),              intent(in)  :: s_child, s_parent 
    1418     real,    dimension(2),              intent(in)  :: ds_child, ds_parent 
     1419    real(kind=8),    dimension(2),              intent(in)  :: s_child, s_parent 
     1420    real(kind=8),    dimension(2),              intent(in)  :: ds_child, ds_parent 
    14191421    real,    dimension(                 & 
    14201422        indmin(1):indmax(1),            & 
     
    15341536    integer, dimension(3),              intent(in)  :: indmin, indmax 
    15351537    integer, dimension(3),              intent(in)  :: pttab_child, petab_child 
    1536     real,    dimension(3),              intent(in)  :: s_child, s_parent 
    1537     real,    dimension(3),              intent(in)  :: ds_child, ds_parent 
     1538    real(kind=8),    dimension(3),              intent(in)  :: s_child, s_parent 
     1539    real(kind=8),    dimension(3),              intent(in)  :: ds_child, ds_parent 
    15381540    real,    dimension(                 & 
    15391541        indmin(1):indmax(1),            & 
     
    16411643    integer, dimension(4),              intent(in)  :: indmin, indmax 
    16421644    integer, dimension(4),              intent(in)  :: pttab_child, petab_child 
    1643     real,    dimension(4),              intent(in)  :: s_child, s_parent 
    1644     real,    dimension(4),              intent(in)  :: ds_child, ds_parent 
     1645    real(kind=8),    dimension(4),              intent(in)  :: s_child, s_parent 
     1646    real(kind=8),    dimension(4),              intent(in)  :: ds_child, ds_parent 
    16451647    real,    dimension(                 & 
    16461648        indmin(1):indmax(1),            & 
     
    17071709    integer, dimension(5),              intent(in)  :: indmin, indmax 
    17081710    integer, dimension(5),              intent(in)  :: pttab_child, petab_child 
    1709     real,    dimension(5),              intent(in)  :: s_child, s_parent 
    1710     real,    dimension(5),              intent(in)  :: ds_child, ds_parent 
     1711    real(kind=8),    dimension(5),              intent(in)  :: s_child, s_parent 
     1712    real(kind=8),    dimension(5),              intent(in)  :: ds_child, ds_parent 
    17111713    real,    dimension(                 & 
    17121714        indmin(1):indmax(1),            & 
     
    17801782    integer, dimension(6),                  intent(in)  :: indmin, indmax 
    17811783    integer, dimension(6),                  intent(in)  :: pttab_child, petab_child 
    1782     real,    dimension(6),                  intent(in)  :: s_child, s_parent 
    1783     real,    dimension(6),                  intent(in)  :: ds_child, ds_parent 
     1784    real(kind=8),    dimension(6),                  intent(in)  :: s_child, s_parent 
     1785    real(kind=8),    dimension(6),                  intent(in)  :: ds_child, ds_parent 
    17841786    real,    dimension(                 & 
    17851787        indmin(1):indmax(1),            & 
     
    18591861    REAL, DIMENSION(indmin:indmax),           INTENT(IN)    :: parenttab 
    18601862    REAL, DIMENSION(pttab_child:petab_child), INTENT(OUT)   :: childtab 
    1861     REAL                                                    :: s_parent, s_child 
    1862     REAL                                                    :: ds_parent,ds_child 
     1863    REAL(kind=8)                                            :: s_parent, s_child 
     1864    REAL(kind=8)                                            :: ds_parent,ds_child 
    18631865! 
    18641866    if ( (indmin == indmax) .and. (pttab_child == petab_child) ) then 
  • vendors/AGRIF/dev_r14608_AGRIF_domcfg/AGRIF_FILES/modsauv.F90

    r12420 r14622  
    251251! 
    252252    type(Agrif_PGrid), pointer  :: parcours ! Pointer for the recursive procedure 
    253     real    :: g_eps, eps, oldgrid_eps 
     253    real(kind=8)    :: g_eps, eps, oldgrid_eps 
    254254    integer :: out 
    255255    integer :: iii 
     
    332332! 
    333333    type(Agrif_PGrid), pointer  :: parcours ! Pointer for the recursive procedure 
    334     real    :: g_eps,eps,oldgrid_eps 
     334    real(kind=8)    :: g_eps,eps,oldgrid_eps 
    335335    integer :: out 
    336336    integer :: iii 
     
    416416    integer, dimension(6) :: nbtabold  ! Number of cells in each direction 
    417417    integer, dimension(6) :: nbtabnew  ! Number of cells in each direction 
    418     real,    dimension(6) :: snew,sold 
    419     real,    dimension(6) :: dsnew,dsold 
    420     real    :: eps 
     418    real(kind=8),    dimension(6) :: snew,sold 
     419    real(kind=8),    dimension(6) :: dsnew,dsold 
     420    real(kind=8)    :: eps 
    421421    integer :: n 
    422422! 
     
    532532    integer, dimension(nbdim),     intent(in)    :: pttabold 
    533533    integer, dimension(nbdim),     intent(in)    :: petabold 
    534     real,    dimension(nbdim),     intent(in)    :: snew, sold 
    535     real,    dimension(nbdim),     intent(in)    :: dsnew,dsold 
     534    real(kind=8),    dimension(nbdim),     intent(in)    :: snew, sold 
     535    real(kind=8),    dimension(nbdim),     intent(in)    :: dsnew,dsold 
    536536    integer,                       intent(in)    :: nbdim 
    537537! 
    538538    integer :: i,j,k,l,m,n,i0,j0,k0,l0,m0,n0 
    539539! 
    540     real,    dimension(nbdim) :: dim_gmin,   dim_gmax 
    541     real,    dimension(nbdim) :: dim_newmin, dim_newmax 
    542     real,    dimension(nbdim) :: dim_min 
     540    real(kind=8),    dimension(nbdim) :: dim_gmin,   dim_gmax 
     541    real(kind=8),    dimension(nbdim) :: dim_newmin, dim_newmax 
     542    real(kind=8),    dimension(nbdim) :: dim_min 
    543543    integer, dimension(nbdim) :: ind_gmin,ind_newmin, ind_newmax 
    544544! 
  • vendors/AGRIF/dev_r14608_AGRIF_domcfg/AGRIF_FILES/modtypes.F90

    r14107 r14622  
    380380    real                  :: Agrif_Efficiency = 0.7 
    381381    integer               :: MaxSearch = 5 
    382     real, dimension(3)    :: Agrif_mind 
     382    real(kind=8), dimension(3)    :: Agrif_mind 
    383383!> @} 
    384384!> \name parameters for the interpolation of the child grids 
     
    467467integer function Agrif_Ceiling ( x ) 
    468468!--------------------------------------------------------------------------------------------------- 
    469     real,   intent(in) :: x 
     469    real(kind=8),intent(in) :: x 
    470470! 
    471471    integer   :: i 
     
    487487    integer function Agrif_Int(x) 
    488488!--------------------------------------------------------------------------------------------------- 
    489     real,   intent(in) :: x 
     489    real(kind=8),intent(in) :: x 
    490490! 
    491491    integer :: i 
  • vendors/AGRIF/dev_r14608_AGRIF_domcfg/AGRIF_FILES/modupdate.F90

    r14107 r14622  
    5858    integer, dimension(6) :: ub_child 
    5959    integer, dimension(6) :: lb_parent 
    60     real   , dimension(6) ::  s_child           ! Child  grid position (s_root = 0) 
    61     real   , dimension(6) ::  s_parent          ! Parent grid position (s_root = 0) 
    62     real   , dimension(6) :: ds_child           ! Child  grid dx (ds_root = 1) 
    63     real   , dimension(6) :: ds_parent          ! Parent grid dx (ds_root = 1) 
     60    real(kind=8)   , dimension(6) ::  s_child           ! Child  grid position (s_root = 0) 
     61    real(kind=8)   , dimension(6) ::  s_parent          ! Parent grid position (s_root = 0) 
     62    real(kind=8)   , dimension(6) :: ds_child           ! Child  grid dx (ds_root = 1) 
     63    real(kind=8)   , dimension(6) :: ds_parent          ! Parent grid dx (ds_root = 1) 
    6464    logical, dimension(6) :: do_update          ! Indicates if we perform update for each dimension 
    6565    integer, dimension(6) :: posvar             ! Position of the variable on the cell (1 or 2) 
     
    160160    integer, dimension(nbdim), intent(in) :: posvar         !< Position of the variable on the cell (1 or 2) 
    161161    logical, dimension(nbdim), intent(in) :: do_update      !< Indicates if we update for each dimension 
    162     real,    dimension(nbdim), intent(in) :: s_child        !< Positions of the child grid 
    163     real,    dimension(nbdim), intent(in) :: s_parent       !< Positions of the parent grid 
    164     real,    dimension(nbdim), intent(in) :: ds_child       !< Space steps of the child grid 
    165     real,    dimension(nbdim), intent(in) :: ds_parent      !< Space steps of the parent grid 
     162    real(kind=8),    dimension(nbdim), intent(in) :: s_child        !< Positions of the child grid 
     163    real(kind=8),    dimension(nbdim), intent(in) :: s_parent       !< Positions of the parent grid 
     164    real(kind=8),    dimension(nbdim), intent(in) :: ds_child       !< Space steps of the child grid 
     165    real(kind=8),    dimension(nbdim), intent(in) :: ds_parent      !< Space steps of the parent grid 
    166166    procedure()                           :: procname       !< Data recovery procedure 
    167167! 
     
    274274    integer, dimension(nbdim), intent(in)   :: posvar       !< Position of the variable on the cell (1 or 2) 
    275275    logical, dimension(nbdim), intent(in)   :: do_update    !< Indicates if we update for each dimension 
    276     real,    dimension(nbdim), intent(in)   :: s_child      !< Positions of the child grid 
    277     real,    dimension(nbdim), intent(in)   :: s_parent     !< Positions of the parent grid 
    278     real,    dimension(nbdim), intent(in)   :: ds_child     !< Space steps of the child grid 
    279     real,    dimension(nbdim), intent(in)   :: ds_parent    !< Space steps of the parent grid 
     276    real(kind=8),    dimension(nbdim), intent(in)   :: s_child      !< Positions of the child grid 
     277    real(kind=8),    dimension(nbdim), intent(in)   :: s_parent     !< Positions of the parent grid 
     278    real(kind=8),    dimension(nbdim), intent(in)   :: ds_child     !< Space steps of the child grid 
     279    real(kind=8),    dimension(nbdim), intent(in)   :: ds_parent    !< Space steps of the parent grid 
    280280    procedure()                             :: procname     !< Data recovery procedure 
    281281    integer, dimension(6)  :: loctab_child      ! Indicates if the child grid has a common border 
     
    444444    integer, dimension(nbdim), intent(in)   :: lb_parent !< Index of the first point inside the domain for the parent 
    445445                                                            !!    grid variable 
    446     real,    dimension(nbdim), intent(in)   :: s_child      !< Positions of the child grid 
    447     real,    dimension(nbdim), intent(in)   :: s_parent     !< Positions of the parent grid 
    448     real,    dimension(nbdim), intent(in)   :: ds_child     !< Space steps of the child grid 
    449     real,    dimension(nbdim), intent(in)   :: ds_parent    !< Space steps of the parent grid 
     446    real(kind=8),    dimension(nbdim), intent(in)   :: s_child      !< Positions of the child grid 
     447    real(kind=8),    dimension(nbdim), intent(in)   :: s_parent     !< Positions of the parent grid 
     448    real(kind=8),    dimension(nbdim), intent(in)   :: ds_child     !< Space steps of the child grid 
     449    real(kind=8),    dimension(nbdim), intent(in)   :: ds_parent    !< Space steps of the parent grid 
    450450    procedure()                             :: procname     !< Data recovery procedure 
    451451    integer, optional,         intent(in)   :: nb, ndir 
     
    459459    integer, dimension(nbdim)       :: indmin, indmax 
    460460    integer, dimension(nbdim)       :: indminglob, indmaxglob 
    461     real   , dimension(nbdim)       :: s_Child_temp, s_Parent_temp 
     461    real(kind=8)   , dimension(nbdim)       :: s_Child_temp, s_Parent_temp 
    462462    integer, dimension(nbdim)       :: lowerbound,upperbound 
    463463    integer, dimension(nbdim)       :: pttruetabwhole, cetruetabwhole 
     
    17961796    integer,                   intent(in)   :: nbdim 
    17971797    integer, dimension(nbdim), intent(out)  :: indmin, indmax 
    1798     real,    dimension(nbdim), intent(out)  :: s_Parent_temp, s_Child_temp 
    1799     real,    dimension(nbdim), intent(in)   :: s_child,  ds_child 
    1800     real,    dimension(nbdim), intent(in)   :: s_parent, ds_parent 
     1798    real(kind=8),    dimension(nbdim), intent(out)  :: s_Parent_temp, s_Child_temp 
     1799    real(kind=8),    dimension(nbdim), intent(in)   :: s_child,  ds_child 
     1800    real(kind=8),    dimension(nbdim), intent(in)   :: s_parent, ds_parent 
    18011801    integer, dimension(nbdim), intent(in)   :: pttruetab, cetruetab 
    18021802    integer, dimension(nbdim), intent(in)   :: lb_child, lb_parent 
     
    18081808#endif 
    18091809! 
    1810     real,dimension(nbdim) :: dim_newmin,dim_newmax 
     1810    real(kind=8),dimension(nbdim) :: dim_newmin,dim_newmax 
    18111811    integer :: i 
    18121812#if defined AGRIF_MPI 
    1813     real    :: positionmin, positionmax 
     1813    real(kind=8)    :: positionmin, positionmax 
    18141814    integer :: imin, imax 
    18151815    integer :: coeffraf 
     
    19051905    integer,                            intent(in)  :: indmin, indmax 
    19061906    integer,                            intent(in)  :: lb_child, ub_child 
    1907     real,                               intent(in)  ::  s_child,  s_parent 
    1908     real,                               intent(in)  :: ds_child, ds_parent 
     1907    real(kind=8),                               intent(in)  ::  s_child,  s_parent 
     1908    real(kind=8),                               intent(in)  :: ds_child, ds_parent 
    19091909    real, dimension(indmin:indmax),     intent(out) :: tempP 
    19101910    real, dimension(lb_child:ub_child), intent(in)  :: tempC 
     
    19371937    integer, dimension(2),          intent(in)  :: indmin, indmax 
    19381938    integer, dimension(2),          intent(in)  :: lb_child, ub_child 
    1939     real,    dimension(2),          intent(in)  ::  s_child,  s_parent 
    1940     real,    dimension(2),          intent(in)  :: ds_child, ds_parent 
     1939    real(kind=8),    dimension(2),          intent(in)  ::  s_child,  s_parent 
     1940    real(kind=8),    dimension(2),          intent(in)  :: ds_child, ds_parent 
    19411941    real,    dimension(          & 
    19421942        indmin(1):indmax(1),     & 
  • vendors/AGRIF/dev_r14608_AGRIF_domcfg/AGRIF_FILES/modupdatebasic.F90

    r14107 r14622  
    4949    integer,             intent(in)     :: np           !< Length of parent array 
    5050    integer,             intent(in)     :: nc           !< Length of child  array 
    51     real,                intent(in)     :: s_parent     !< Parent grid position (s_root = 0) 
    52     real,                intent(in)     :: s_child      !< Child  grid position (s_root = 0) 
    53     real,                intent(in)     :: ds_parent    !< Parent grid dx (ds_root = 1) 
    54     real,                intent(in)     :: ds_child     !< Child  grid dx (ds_root = 1) 
     51    real(kind=8),        intent(in)     :: s_parent     !< Parent grid position (s_root = 0) 
     52    real(kind=8),        intent(in)     :: s_child      !< Child  grid position (s_root = 0) 
     53    real(kind=8),        intent(in)     :: ds_parent    !< Parent grid dx (ds_root = 1) 
     54    real(kind=8),        intent(in)     :: ds_child     !< Child  grid dx (ds_root = 1) 
    5555!--------------------------------------------------------------------------------------------------- 
    5656    integer :: i, locind_child_left, coeffraf 
     
    8484    integer,             intent(in)     :: np           !< Length of parent array 
    8585    integer,             intent(in)     :: nc           !< Length of child  array 
    86     real,                intent(in)     :: s_parent     !< Parent grid position (s_root = 0) 
    87     real,                intent(in)     :: s_child      !< Child  grid position (s_root = 0) 
    88     real,                intent(in)     :: ds_parent    !< Parent grid dx (ds_root = 1) 
    89     real,                intent(in)     :: ds_child     !< Child  grid dx (ds_root = 1) 
     86    real(kind=8),        intent(in)     :: s_parent     !< Parent grid position (s_root = 0) 
     87    real(kind=8),        intent(in)     :: s_child      !< Child  grid position (s_root = 0) 
     88    real(kind=8),        intent(in)     :: ds_parent    !< Parent grid dx (ds_root = 1) 
     89    real(kind=8),        intent(in)     :: ds_child     !< Child  grid dx (ds_root = 1) 
    9090    integer,             intent(in)     :: dir          !< Direction 
    9191!--------------------------------------------------------------------------------------------------- 
     
    157157    REAL, DIMENSION(nc), intent(in)     :: y 
    158158    INTEGER,             intent(in)     :: np,nc 
    159     REAL,                intent(in)     :: s_parent,  s_child 
    160     REAL,                intent(in)     :: ds_parent, ds_child 
     159    REAL(kind=8),        intent(in)     :: s_parent,  s_child 
     160    REAL(kind=8),        intent(in)     :: ds_parent, ds_child 
    161161! 
    162162    INTEGER :: i, ii, locind_child_left, coeffraf 
    163     REAL    :: xpos, invcoeffraf 
     163    REAL(kind=8)    :: xpos 
     164    REAL ::  invcoeffraf 
    164165    INTEGER :: nbnonnuls 
    165166    INTEGER :: diffmod 
     
    294295!--------------------------------------------------------------------------------------------------- 
    295296    INTEGER, intent(in) :: nc2, np, nc 
    296     REAL,    intent(in) :: s_parent,  s_child 
    297     REAL,    intent(in) :: ds_parent, ds_child 
     297    REAL(kind=8),    intent(in) :: s_parent,  s_child 
     298    REAL(kind=8),    intent(in) :: ds_parent, ds_child 
    298299    INTEGER, intent(in) :: dir 
    299300! 
    300301    INTEGER, DIMENSION(:,:), ALLOCATABLE :: indchildaverage_tmp 
    301302    INTEGER :: i, locind_child_left, coeffraf 
    302     REAL    :: xpos 
     303    REAL(kind=8)    :: xpos 
    303304    INTEGER :: diffmod 
    304305! 
     
    346347    REAL, DIMENSION(nc), intent(in)     :: y 
    347348    INTEGER,             intent(in)     :: np, nc 
    348     REAL,                intent(in)     :: s_parent,  s_child 
    349     REAL,                intent(in)     :: ds_parent, ds_child 
     349    REAL(kind=8),                intent(in)     :: s_parent,  s_child 
     350    REAL(kind=8),                intent(in)     :: ds_parent, ds_child 
    350351    INTEGER,             intent(in)     :: dir 
    351352! 
     
    403404    real, dimension(nc), intent(in)     :: y 
    404405    integer,             intent(in)     :: np, nc 
    405     real,                intent(in)     :: s_parent,  s_child 
    406     real,                intent(in)     :: ds_parent, ds_child 
    407 !--------------------------------------------------------------------------------------------------- 
    408     REAL    :: xpos, xposfin 
     406    real(kind=8),                intent(in)     :: s_parent,  s_child 
     407    real(kind=8),                intent(in)     :: ds_parent, ds_child 
     408!--------------------------------------------------------------------------------------------------- 
     409    REAL(kind=8)    :: xpos, xposfin 
    409410    INTEGER :: i, ii, diffmod 
    410411    INTEGER :: it1, it2 
Note: See TracChangeset for help on using the changeset viewer.