MOM6
mom_domains::pass_var Interface Reference

Detailed Description

Do a halo update on an array.

Definition at line 49 of file MOM_domains.F90.

Private functions

subroutine pass_var_3d (array, MOM_dom, sideflag, complete, position, halo, clock)
 pass_var_3d does a halo update for a three-dimensional array. More...
 
subroutine pass_var_2d (array, MOM_dom, sideflag, complete, position, halo, inner_halo, clock)
 pass_var_2d does a halo update for a two-dimensional array. More...
 

Functions and subroutines

◆ pass_var_2d()

subroutine mom_domains::pass_var::pass_var_2d ( real, dimension(:,:), intent(inout)  array,
type(mom_domain_type), intent(inout)  MOM_dom,
integer, intent(in), optional  sideflag,
logical, intent(in), optional  complete,
integer, intent(in), optional  position,
integer, intent(in), optional  halo,
integer, intent(in), optional  inner_halo,
integer, intent(in), optional  clock 
)
private

pass_var_2d does a halo update for a two-dimensional array.

Parameters
[in,out]arrayThe array which is having its halos points exchanged.
[in,out]mom_domThe MOM_domain_type containing the mpp_domain needed to determine where data should be sent.
[in]sideflagAn optional integer indicating which directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST, TO_NORTH, and TO_SOUTH. For example, TO_EAST sends the data to the processor to the east, so the halos on the western side are filled. TO_ALL is the default if sideflag is omitted.
[in]completeAn optional argument indicating whether the halo updates should be completed before progress resumes. Omitting complete is the same as setting complete to .true.
[in]positionAn optional argument indicating the position. This is usally CORNER, but is CENTER by default.
[in]haloThe size of the halo to update - the full halo by default.
[in]inner_haloThe size of an inner halo to avoid updating, or 0 to avoid updating symmetric memory computational domain points. Setting this >=0 also enforces that complete=.true.
[in]clockThe handle for a cpu time clock that should be started then stopped to time this routine.

Definition at line 185 of file MOM_domains.F90.

185  real, dimension(:,:), intent(inout) :: array !< The array which is having its halos points
186  !! exchanged.
187  type(MOM_domain_type), intent(inout) :: MOM_dom !< The MOM_domain_type containing the mpp_domain
188  !! needed to determine where data should be sent.
189  integer, optional, intent(in) :: sideflag !< An optional integer indicating which
190  !! directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST,
191  !! TO_NORTH, and TO_SOUTH. For example, TO_EAST sends the data to the processor to the east,
192  !! so the halos on the western side are filled. TO_ALL is the default if sideflag is omitted.
193  logical, optional, intent(in) :: complete !< An optional argument indicating whether the
194  !! halo updates should be completed before
195  !! progress resumes. Omitting complete is the
196  !! same as setting complete to .true.
197  integer, optional, intent(in) :: position !< An optional argument indicating the position.
198  !! This is usally CORNER, but is CENTER
199  !! by default.
200  integer, optional, intent(in) :: halo !< The size of the halo to update - the full halo
201  !! by default.
202  integer, optional, intent(in) :: inner_halo !< The size of an inner halo to avoid updating,
203  !! or 0 to avoid updating symmetric memory
204  !! computational domain points. Setting this >=0
205  !! also enforces that complete=.true.
206  integer, optional, intent(in) :: clock !< The handle for a cpu time clock that should be
207  !! started then stopped to time this routine.
208 
209  ! Local variables
210  real, allocatable, dimension(:,:) :: tmp
211  integer :: pos, i_halo, j_halo
212  integer :: isc, iec, jsc, jec, isd, ied, jsd, jed, IscB, IecB, JscB, JecB
213  integer :: inner, i, j, isfw, iefw, isfe, iefe, jsfs, jefs, jsfn, jefn
214  integer :: dirflag
215  logical :: block_til_complete
216 
217  if (present(clock)) then ; if (clock>0) call cpu_clock_begin(clock) ; endif
218 
219  dirflag = to_all ! 60
220  if (present(sideflag)) then ; if (sideflag > 0) dirflag = sideflag ; endif
221  block_til_complete = .true. ; if (present(complete)) block_til_complete = complete
222  pos = center ; if (present(position)) pos = position
223 
224  if (present(inner_halo)) then ; if (inner_halo >= 0) then
225  ! Store the original values.
226  allocate(tmp(size(array,1), size(array,2)))
227  tmp(:,:) = array(:,:)
228  block_til_complete = .true.
229  endif ; endif
230 
231  if (present(halo) .and. mom_dom%thin_halo_updates) then
232  call mpp_update_domains(array, mom_dom%mpp_domain, flags=dirflag, &
233  complete=block_til_complete, position=position, &
234  whalo=halo, ehalo=halo, shalo=halo, nhalo=halo)
235  else
236  call mpp_update_domains(array, mom_dom%mpp_domain, flags=dirflag, &
237  complete=block_til_complete, position=position)
238  endif
239 
240  if (present(inner_halo)) then ; if (inner_halo >= 0) then
241  call mpp_get_compute_domain(mom_dom%mpp_domain, isc, iec, jsc, jec)
242  call mpp_get_data_domain(mom_dom%mpp_domain, isd, ied, jsd, jed)
243  ! Convert to local indices for arrays starting at 1.
244  isc = isc - (isd-1) ; iec = iec - (isd-1) ; ied = ied - (isd-1) ; isd = 1
245  jsc = jsc - (jsd-1) ; jec = jec - (jsd-1) ; jed = jed - (jsd-1) ; jsd = 1
246  i_halo = min(inner_halo, isc-1) ; j_halo = min(inner_halo, jsc-1)
247 
248  ! Figure out the array index extents of the eastern, western, northern and southern regions to copy.
249  if (pos == center) then
250  if (size(array,1) == ied) then
251  isfw = isc - i_halo ; iefw = isc ; isfe = iec ; iefe = iec + i_halo
252  else ; call mom_error(fatal, "pass_var_2d: wrong i-size for CENTER array.") ; endif
253  if (size(array,2) == jed) then
254  isfw = isc - i_halo ; iefw = isc ; isfe = iec ; iefe = iec + i_halo
255  else ; call mom_error(fatal, "pass_var_2d: wrong j-size for CENTER array.") ; endif
256  elseif (pos == corner) then
257  if (size(array,1) == ied) then
258  isfw = max(isc - (i_halo+1), 1) ; iefw = isc ; isfe = iec ; iefe = iec + i_halo
259  elseif (size(array,1) == ied+1) then
260  isfw = isc - i_halo ; iefw = isc+1 ; isfe = iec+1 ; iefe = min(iec + 1 + i_halo, ied+1)
261  else ; call mom_error(fatal, "pass_var_2d: wrong i-size for CORNER array.") ; endif
262  if (size(array,2) == jed) then
263  jsfs = max(jsc - (j_halo+1), 1) ; jefs = jsc ; jsfn = jec ; jefn = jec + j_halo
264  elseif (size(array,2) == jed+1) then
265  jsfs = jsc - j_halo ; jefs = jsc+1 ; jsfn = jec+1 ; jefn = min(jec + 1 + j_halo, jed+1)
266  else ; call mom_error(fatal, "pass_var_2d: wrong j-size for CORNER array.") ; endif
267  else
268  call mom_error(fatal, "pass_var_2d: Unrecognized position")
269  endif
270 
271  ! Copy back the stored inner halo points
272  do j=jsfs,jefn ; do i=isfw,iefw ; array(i,j) = tmp(i,j) ; enddo ; enddo
273  do j=jsfs,jefn ; do i=isfe,iefe ; array(i,j) = tmp(i,j) ; enddo ; enddo
274  do j=jsfs,jefs ; do i=isfw,iefe ; array(i,j) = tmp(i,j) ; enddo ; enddo
275  do j=jsfn,jefn ; do i=isfw,iefe ; array(i,j) = tmp(i,j) ; enddo ; enddo
276 
277  deallocate(tmp)
278  endif ; endif
279 
280  if (present(clock)) then ; if (clock>0) call cpu_clock_end(clock) ; endif
281 

◆ pass_var_3d()

subroutine mom_domains::pass_var::pass_var_3d ( real, dimension(:,:,:), intent(inout)  array,
type(mom_domain_type), intent(inout)  MOM_dom,
integer, intent(in), optional  sideflag,
logical, intent(in), optional  complete,
integer, intent(in), optional  position,
integer, intent(in), optional  halo,
integer, intent(in), optional  clock 
)
private

pass_var_3d does a halo update for a three-dimensional array.

Parameters
[in,out]arrayThe array which is having its halos points exchanged.
[in,out]mom_domThe MOM_domain_type containing the mpp_domain needed to determine where data should be sent.
[in]sideflagAn optional integer indicating which directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST, TO_NORTH, and TO_SOUTH. For example, TO_EAST sends the data to the processor to the east, sothe halos on the western side are filled. TO_ALL is the default if sideflag is omitted.
[in]completeAn optional argument indicating whether the halo updates should be completed before progress resumes. Omitting complete is the same as setting complete to .true.
[in]positionAn optional argument indicating the position. This is usally CORNER, but is CENTER by default.
[in]haloThe size of the halo to update - the full halo by default.
[in]clockThe handle for a cpu time clock that should be started then stopped to time this routine.

Definition at line 139 of file MOM_domains.F90.

139  real, dimension(:,:,:), intent(inout) :: array !< The array which is having its halos points
140  !! exchanged.
141  type(MOM_domain_type), intent(inout) :: MOM_dom !< The MOM_domain_type containing the mpp_domain
142  !! needed to determine where data should be
143  !! sent.
144  integer, optional, intent(in) :: sideflag !< An optional integer indicating which
145  !! directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST,
146  !! TO_NORTH, and TO_SOUTH. For example, TO_EAST sends the data to the processor to the east,
147  !! sothe halos on the western side are filled. TO_ALL is the default if sideflag is omitted.
148  logical, optional, intent(in) :: complete !< An optional argument indicating whether the
149  !! halo updates should be completed before
150  !! progress resumes. Omitting complete is the
151  !! same as setting complete to .true.
152  integer, optional, intent(in) :: position !< An optional argument indicating the position.
153  !! This is usally CORNER, but is CENTER by
154  !! default.
155  integer, optional, intent(in) :: halo !< The size of the halo to update - the full
156  !! halo by default.
157  integer, optional, intent(in) :: clock !< The handle for a cpu time clock that should be
158  !! started then stopped to time this routine.
159 
160  integer :: dirflag
161  logical :: block_til_complete
162 
163  if (present(clock)) then ; if (clock>0) call cpu_clock_begin(clock) ; endif
164 
165  dirflag = to_all ! 60
166  if (present(sideflag)) then ; if (sideflag > 0) dirflag = sideflag ; endif
167  block_til_complete = .true.
168  if (present(complete)) block_til_complete = complete
169 
170  if (present(halo) .and. mom_dom%thin_halo_updates) then
171  call mpp_update_domains(array, mom_dom%mpp_domain, flags=dirflag, &
172  complete=block_til_complete, position=position, &
173  whalo=halo, ehalo=halo, shalo=halo, nhalo=halo)
174  else
175  call mpp_update_domains(array, mom_dom%mpp_domain, flags=dirflag, &
176  complete=block_til_complete, position=position)
177  endif
178 
179  if (present(clock)) then ; if (clock>0) call cpu_clock_end(clock) ; endif
180 

The documentation for this interface was generated from the following file: