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 10087 for vendors/AGRIF/CMEMS_2020/AGRIF_FILES/modarrays.F90 – NEMO

Ignore:
Timestamp:
2018-09-05T15:33:44+02:00 (6 years ago)
Author:
rblod
Message:

update AGRIF library

File:
1 edited

Legend:

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

    r5656 r10087  
    5555                               proc_id,         & 
    5656                               coords,          & 
    57                                lb_tab_true, ub_tab_true, memberin ) 
     57                               lb_tab_true, ub_tab_true, memberin,  & 
     58                               indminglob3,indmaxglob3) 
    5859!--------------------------------------------------------------------------------------------------- 
    5960    integer,                   intent(in)  :: nbdim         !< Number of dimensions 
     
    6162    integer, dimension(nbdim), intent(in)  :: ub_var        !< Local upper boundary on the current processor 
    6263    integer, dimension(nbdim), intent(in)  :: lb_tab        !< Global lower boundary of the variable 
     64    integer, dimension(nbdim),OPTIONAL     :: indminglob3,indmaxglob3 !< True bounds for MPI USE 
    6365    integer, dimension(nbdim), intent(in)  :: ub_tab        !< Global upper boundary of the variable 
    6466    integer,                   intent(in)  :: proc_id       !< Current processor 
     
    7880        call Agrif_InvLoc( lb_var(i), proc_id, coord_i, lb_glob_index ) 
    7981        call Agrif_InvLoc( ub_var(i), proc_id, coord_i, ub_glob_index ) 
     82        if (present(indminglob3)) then 
     83          indminglob3(i)=lb_glob_index 
     84          indmaxglob3(i)=ub_glob_index 
     85        endif 
    8086#else 
    8187        lb_glob_index = lb_var(i) 
    8288        ub_glob_index = ub_var(i) 
    8389#endif 
     90 
    8491        lb_tab_true(i) = max(lb_tab(i), lb_glob_index) 
    8592        ub_tab_true(i) = min(ub_tab(i), ub_glob_index) 
     93 
    8694    enddo 
    8795! 
     
    123131! 
    124132    iminmaxg(1:nbdim,2) = - iminmaxg(1:nbdim,2) 
    125     call MPI_ALLREDUCE(iminmaxg, lubglob, 2*nbdim, MPI_INTEGER, MPI_MIN, Agrif_mpi_comm, code) 
     133    call MPI_ALLREDUCE(iminmaxg, lubglob, 2*nbdim, MPI_INTEGER, MPI_MIN, & 
     134                       Agrif_mpi_comm, code) 
    126135    lubglob(1:nbdim,2)  = - lubglob(1:nbdim,2) 
    127136#endif 
     
    225234    case (1) ; call Agrif_set_array_tozero_1D(variable%array1) 
    226235    case (2) ; call Agrif_set_array_tozero_2D(variable%array2) 
    227     case (3) ; call Agrif_set_array_tozero_3D(variable%array3) 
     236    case (3) ; call Agrif_set_array_tozero_reshape(variable%array3,size(variable%array3)) 
     237!case (3) ; call Agrif_set_array_tozero_3D(variable%array3) 
    228238    case (4) ; call Agrif_set_array_tozero_4D(variable%array4) 
    229239    case (5) ; call Agrif_set_array_tozero_5D(variable%array5) 
     
    266276!=================================================================================================== 
    267277! 
     278!=================================================================================================== 
     279! 
     280!=================================================================================================== 
     281!  subroutine agrif_set_array_cond 
     282! 
     283!> Compute the masking of \b variablein, according to the required dimension. 
     284!--------------------------------------------------------------------------------------------------- 
     285subroutine agrif_set_array_cond ( variablein, variableout, value, nbdim ) 
     286!--------------------------------------------------------------------------------------------------- 
     287    type(Agrif_Variable), intent(in)    :: variablein     !< Variablein 
     288    type(Agrif_Variable), intent(inout) :: variableout    !< Variableout 
     289    real,intent(in) :: value                            !< special value 
     290    integer, intent(in)                 :: nbdim        !< Dimension of the array 
     291 
     292! 
     293    select case (nbdim) 
     294    case (1) ; call agrif_set_array_cond_1D(variablein%array1,variableout%array1,value) 
     295    case (2) ; call agrif_set_array_cond_2D(variablein%array2,variableout%array2,value) 
     296    case (3) ; call agrif_set_array_cond_reshape(variablein%array3,variableout%array3,value,size(variablein%array3)) 
     297!    case (3) ; call agrif_set_array_cond_3D(variablein%array3,variableout%array3,value) 
     298    case (4) ; call agrif_set_array_cond_4D(variablein%array4,variableout%array4,value) 
     299    case (5) ; call agrif_set_array_cond_5D(variablein%array5,variableout%array5,value) 
     300    case (6) ; call agrif_set_array_cond_6D(variablein%array6,variableout%array6,value) 
     301    end select 
     302!--------------------------------------------------------------------------------------------------- 
     303contains 
     304!--------------------------------------------------------------------------------------------------- 
     305subroutine agrif_set_array_cond_1D(arrayin,arrayout,value) 
     306real,dimension(:),intent(in) :: arrayin 
     307real,dimension(:),intent(out) :: arrayout 
     308real :: value 
     309 
     310where (arrayin == value) 
     311  arrayout = 0. 
     312elsewhere 
     313  arrayout = 1. 
     314end where 
     315 
     316end subroutine agrif_set_array_cond_1D 
     317! 
     318subroutine agrif_set_array_cond_2D(arrayin,arrayout,value) 
     319real,dimension(:,:),intent(in) :: arrayin 
     320real,dimension(:,:),intent(out) :: arrayout 
     321real :: value 
     322 
     323where (arrayin == value) 
     324  arrayout = 0. 
     325elsewhere 
     326  arrayout = 1. 
     327end where 
     328 
     329end subroutine agrif_set_array_cond_2D 
     330! 
     331subroutine agrif_set_array_cond_3D(arrayin,arrayout,value) 
     332real,dimension(:,:,:),intent(in) :: arrayin 
     333real,dimension(:,:,:),intent(out) :: arrayout 
     334real :: value 
     335 
     336where (arrayin == value) 
     337  arrayout = 0. 
     338elsewhere 
     339  arrayout = 1. 
     340end where 
     341 
     342end subroutine agrif_set_array_cond_3D 
     343! 
     344subroutine agrif_set_array_cond_4D(arrayin,arrayout,value) 
     345real,dimension(:,:,:,:),intent(in) :: arrayin 
     346real,dimension(:,:,:,:),intent(out) :: arrayout 
     347real :: value 
     348 
     349where (arrayin == value) 
     350  arrayout = 0. 
     351elsewhere 
     352  arrayout = 1. 
     353end where 
     354 
     355end subroutine agrif_set_array_cond_4D 
     356! 
     357subroutine agrif_set_array_cond_5D(arrayin,arrayout,value) 
     358real,dimension(:,:,:,:,:),intent(in) :: arrayin 
     359real,dimension(:,:,:,:,:),intent(out) :: arrayout 
     360real :: value 
     361 
     362where (arrayin == value) 
     363  arrayout = 0. 
     364elsewhere 
     365  arrayout = 1. 
     366end where 
     367 
     368end subroutine agrif_set_array_cond_5D 
     369! 
     370subroutine agrif_set_array_cond_6D(arrayin,arrayout,value) 
     371real,dimension(:,:,:,:,:,:),intent(in) :: arrayin 
     372real,dimension(:,:,:,:,:,:),intent(out) :: arrayout 
     373real :: value 
     374 
     375where (arrayin == value) 
     376  arrayout = 0. 
     377elsewhere 
     378  arrayout = 1. 
     379end where 
     380 
     381end subroutine agrif_set_array_cond_6D 
     382!--------------------------------------------------------------------------------------------------- 
     383end subroutine agrif_set_array_cond 
    268384!=================================================================================================== 
    269385!  subroutine Agrif_var_copy_array 
     
    330446        real, dimension(l(1):,l(2):,l(3):), intent(out) :: tabout 
    331447        real, dimension(m(1):,m(2):,m(3):), intent(in)  :: tabin 
    332         tabout(inf1(1):sup1(1), & 
    333                inf1(2):sup1(2), & 
    334                inf1(3):sup1(3)) = tabin(inf2(1):sup2(1), & 
    335                                         inf2(2):sup2(2), & 
    336                                         inf2(3):sup2(3)) 
     448        integer :: i,j,k 
     449 
     450 
     451!$OMP PARALLEL DO DEFAULT(none) PRIVATE(i,j,k) & 
     452!$OMP SHARED(inf1,inf2,sup1,sup2,tabin,tabout) & 
     453!$OMP SCHEDULE(RUNTIME) 
     454        do k=inf1(3),sup1(3) 
     455        do j=inf1(2),sup1(2) 
     456        do i=inf1(1),sup1(1) 
     457!          tabout(i,j,k) = tabin(i+inf2(1)-inf1(1),j+inf2(2)-inf1(2),k+inf2(3)-inf1(3)) 
     458          tabout(i,j,k) = tabin(i,j,k) 
     459        enddo 
     460        enddo 
     461        enddo 
     462!$OMP END PARALLEL DO 
     463 
     464 
     465!        tabout(inf1(1):sup1(1), & 
     466!               inf1(2):sup1(2), & 
     467!               inf1(3):sup1(3)) = tabin(inf2(1):sup2(1), & 
     468!                                        inf2(2):sup2(2), & 
     469!                                        inf2(3):sup2(3)) 
    337470    end subroutine Agrif_copy_array_3d 
    338471! 
     
    631764    integer, dimension(6), intent(out)          :: lb_child     !< Lower bound on the child grid 
    632765    integer, dimension(6), intent(out)          :: lb_parent    !< Lower bound on the parent grid 
    633     real, dimension(6),    intent(out)          :: s_child      !< Child  grid position (s_root = 0) 
    634     real, dimension(6),    intent(out)          :: s_parent     !< Parent grid position (s_root = 0) 
    635     real, dimension(6),    intent(out)          :: ds_child     !< Child  grid dx (ds_root = 1) 
    636     real, dimension(6),    intent(out)          :: ds_parent    !< Parent grid dx (ds_root = 1) 
     766    real(kind=8), dimension(6),    intent(out)  :: s_child      !< Child  grid position (s_root = 0) 
     767    real(kind=8), dimension(6),    intent(out)  :: s_parent     !< Parent grid position (s_root = 0) 
     768    real(kind=8), dimension(6),    intent(out)  :: ds_child     !< Child  grid dx (ds_root = 1) 
     769    real(kind=8), dimension(6),    intent(out)  :: ds_parent    !< Parent grid dx (ds_root = 1) 
    637770    integer,               intent(out)          :: nbdim        !< Number of dimensions 
    638771    logical,               intent(in)           :: interp       !< .true. if preprocess for interpolation, \n 
     
    671804            else 
    672805                ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(1) - 1 
    673                 s_child(n)  = s_child(n)  + 0.5*ds_child(n) 
    674                 s_parent(n) = s_parent(n) + 0.5*ds_parent(n) 
     806                s_child(n)  = s_child(n)  + 0.5d0*ds_child(n) 
     807                s_parent(n) = s_parent(n) + 0.5d0*ds_parent(n) 
    675808            endif 
    676809! 
     
    689822            else 
    690823                ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(2) - 1 
    691                 s_child(n)  = s_child(n)  + 0.5*ds_child(n) 
    692                 s_parent(n) = s_parent(n) + 0.5*ds_parent(n) 
     824                s_child(n)  = s_child(n)  + 0.5d0*ds_child(n) 
     825                s_parent(n) = s_parent(n) + 0.5d0*ds_parent(n) 
    693826            endif 
    694827! 
     
    727860!           No interpolation but only a copy of the values of the grid variable 
    728861            lb_parent(n) = lb_child(n) 
    729             s_child(n)   = 0. 
    730             s_parent(n)  = 0. 
    731             ds_child(n)  = 1. 
    732             ds_parent(n) = 1. 
     862            s_child(n)   = 0.d0 
     863            s_parent(n)  = 0.d0 
     864            ds_child(n)  = 1.d0 
     865            ds_parent(n) = 1.d0 
    733866! 
    734867        end select 
     
    803936    do i = 1,nbdim 
    804937! 
     938     if (coords(i) == 0) then 
     939       nbloc(i) = 1 
     940       locbounds(i,1,1) = lb_glob(i) 
     941       locbounds(i,2,1) = ub_glob(i) 
     942       locbounds(i,1,2) = lb_glob(i) 
     943       locbounds(i,2,2) = ub_glob(i) 
     944     else 
    805945        call Agrif_InvLoc(lb_var(i), rank, coords(i), i1) 
    806946! 
     
    816956            endif 
    817957        enddo 
     958     endif 
    818959    enddo 
    819960 
     
    825966! 
    826967end module Agrif_Arrays 
     968 
     969 
     970subroutine agrif_set_array_cond_reshape(arrayin,arrayout,value,n) 
     971integer :: n 
     972real,dimension(n) :: arrayin,arrayout 
     973real :: value 
     974 
     975integer :: i 
     976 
     977do i=1,n 
     978  if (arrayin(i) == value) then 
     979    arrayout(i) = 0. 
     980  else 
     981    arrayout(i) = 1. 
     982  endif 
     983enddo 
     984 
     985end subroutine agrif_set_array_cond_reshape 
     986 
     987subroutine agrif_set_array_tozero_reshape(array,n) 
     988integer :: n 
     989real,dimension(n) :: array 
     990 
     991integer :: i 
     992 
     993!$OMP PARALLEL DO DEFAULT(none) PRIVATE(i) & 
     994!$OMP SHARED(array,n) 
     995do i=1,n 
     996    array(i) = 0. 
     997enddo 
     998!$OMP END PARALLEL DO 
     999 
     1000end subroutine agrif_set_array_tozero_reshape 
Note: See TracChangeset for help on using the changeset viewer.