Changeset 6258 for branches/2015/dev_r5803_UKMO_AGRIF_Vert_interp/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modarrays.F90
- Timestamp:
- 2016-01-15T13:11:56+01:00 (8 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/2015/dev_r5803_UKMO_AGRIF_Vert_interp/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modarrays.F90
r5819 r6258 1 1 ! 2 ! $Id $2 ! $Id: modarrays.F 662 2007-05-25 15:58:52Z opalod $ 3 3 ! 4 4 ! AGRIF (Adaptive Grid Refinement In Fortran) … … 55 55 proc_id, & 56 56 coords, & 57 lb_tab_true, ub_tab_true, memberin ) 57 lb_tab_true, ub_tab_true, memberin, & 58 indminglob3,indmaxglob3) 58 59 !--------------------------------------------------------------------------------------------------- 59 60 integer, intent(in) :: nbdim !< Number of dimensions … … 61 62 integer, dimension(nbdim), intent(in) :: ub_var !< Local upper boundary on the current processor 62 63 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 63 65 integer, dimension(nbdim), intent(in) :: ub_tab !< Global upper boundary of the variable 64 66 integer, intent(in) :: proc_id !< Current processor … … 78 80 call Agrif_InvLoc( lb_var(i), proc_id, coord_i, lb_glob_index ) 79 81 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 80 86 #else 81 87 lb_glob_index = lb_var(i) … … 84 90 lb_tab_true(i) = max(lb_tab(i), lb_glob_index) 85 91 ub_tab_true(i) = min(ub_tab(i), ub_glob_index) 92 86 93 enddo 87 94 ! … … 123 130 ! 124 131 iminmaxg(1:nbdim,2) = - iminmaxg(1:nbdim,2) 125 call MPI_ALLREDUCE(iminmaxg, lubglob, 2*nbdim, MPI_INTEGER, MPI_MIN, Agrif_mpi_comm, code) 132 call MPI_ALLREDUCE(iminmaxg, lubglob, 2*nbdim, MPI_INTEGER, MPI_MIN, & 133 Agrif_mpi_comm, code) 126 134 lubglob(1:nbdim,2) = - lubglob(1:nbdim,2) 127 135 #endif … … 803 811 do i = 1,nbdim 804 812 ! 813 if (coords(i) == 0) then 814 nbloc(i) = 1 815 locbounds(i,1,1) = lb_glob(i) 816 locbounds(i,2,1) = ub_glob(i) 817 locbounds(i,1,2) = lb_glob(i) 818 locbounds(i,2,2) = ub_glob(i) 819 else 805 820 call Agrif_InvLoc(lb_var(i), rank, coords(i), i1) 806 821 ! … … 816 831 endif 817 832 enddo 833 endif 818 834 enddo 819 835
Note: See TracChangeset
for help on using the changeset viewer.