MOM6
mom_domains::pass_var Interface Reference

Detailed Description

Do a halo update on an array.

Definition at line 54 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...
 

Detailed Description

Do a halo update on an array.

Definition at line 54 of file MOM_domains.F90.

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 CENTER by default and is often CORNER, but could also be EAST_FACE or NORTH_FACE.
[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 190 of file MOM_domains.F90.

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

◆ 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 CENTER by default and is often CORNER, but could also be EAST_FACE or NORTH_FACE.
[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 144 of file MOM_domains.F90.

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

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