fix: Fixes for linter action and code style (#869)
[FMS.git] / test_fms / mpp / fill_halo.F90
blobbb8996ce380e3b2c0104f51b098a7ae7cda103ca
1 !***********************************************************************
2 !*                   GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
14 !* for more details.
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS.  If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @author Jessica Liptak
20 !> @brief This module contains routines to fill halos in different domain configurations
21 !! It is required by test_mpp_update_domains_real and test_mpp_update_domains_int.
22 module fill_halo
24 use :: platform_mod
26 implicit none
27 private
28 integer :: whalo = 2, ehalo = 2, shalo = 2, nhalo = 2
29 integer :: nx=64, ny=64, nz=10
31 public :: fill_halo_zero, fill_regular_refinement_halo, fill_regular_mosaic_halo
32 public :: fill_folded_north_halo, fill_folded_south_halo, fill_folded_east_halo, fill_folded_west_halo
34 !> Routines to fill halo regions of 64-bit and 32-bit real arrays on a regular grid
35 interface fill_regular_refinement_halo
36   module procedure fill_regular_refinement_halo_r8
37   module procedure fill_regular_refinement_halo_r4
38   module procedure fill_regular_refinement_halo_i8
39   module procedure fill_regular_refinement_halo_i4
40 end interface
42 !> Routines to fill halo regions of 64-bit and 32-bit real arrays with zeros
43 interface fill_halo_zero
44   module procedure fill_halo_zero_r8
45   module procedure fill_halo_zero_r4
46   module procedure fill_halo_zero_i8
47   module procedure fill_halo_zero_i4
48 end interface
50 !> Routines to fill halo regions of 64-bit and 32-bit real arrays on a mosaic grid
51 interface fill_regular_mosaic_halo
52   module procedure fill_regular_mosaic_halo_r8
53   module procedure fill_regular_mosaic_halo_r4
54   module procedure fill_regular_mosaic_halo_i8
55   module procedure fill_regular_mosaic_halo_i4
56 end interface fill_regular_mosaic_halo
58 !> Routines to fill halo regions of 64-bit and 32-bit real arrays on a domain with a folded north edge
59 interface fill_folded_north_halo
60   module procedure fill_folded_north_halo_r8
61   module procedure fill_folded_north_halo_r4
62   module procedure fill_folded_north_halo_i8
63   module procedure fill_folded_north_halo_i4
64 end interface fill_folded_north_halo
66 !> Routines to fill halo regions of 64-bit and 32-bit real arrays on a domain with a folded south edge
67 interface fill_folded_south_halo
68   module procedure fill_folded_south_halo_r8
69   module procedure fill_folded_south_halo_r4
70   module procedure fill_folded_south_halo_i8
71   module procedure fill_folded_south_halo_i4
72 end interface fill_folded_south_halo
74 !> Routines to fill halo regions of 64-bit and 32-bit real arrays on a domain with a folded east edge
75 interface fill_folded_east_halo
76   module procedure fill_folded_east_halo_r8
77   module procedure fill_folded_east_halo_r4
78   module procedure fill_folded_east_halo_i8
79   module procedure fill_folded_east_halo_i4
80 end interface fill_folded_east_halo
82 !> Routines to fill halo regions of 64-bit and 32-bit real arrays on a domain with a folded west edge
83 interface fill_folded_west_halo
84   module procedure fill_folded_west_halo_r8
85   module procedure fill_folded_west_halo_r4
86   module procedure fill_folded_west_halo_i8
87   module procedure fill_folded_west_halo_i4
88 end interface fill_folded_west_halo
90 contains
92   !> fill the halo region of a 64-bit real array with zeros
93   subroutine fill_halo_zero_r8(data, whalo, ehalo, shalo, nhalo, xshift, yshift, isc, iec, jsc, jec, isd, ied, &
94                               &  jsd, jed)
95     real(kind=r8_kind), dimension(isd:,jsd:,:), intent(inout) :: data
96     integer,                         intent(in) :: isc, iec, jsc, jec, isd, ied, jsd, jed
97     integer,                         intent(in) :: whalo, ehalo, shalo, nhalo, xshift, yshift
99     if(whalo >=0) then
100        data(iec+ehalo+1+xshift:ied+xshift,jsd:jed+yshift,:) = 0
101        data(isd:isc-whalo-1,jsd:jed+yshift,:) = 0
102     else
103        data(iec+1+xshift:iec-ehalo+xshift,jsc+shalo:jec-nhalo+yshift,:) = 0
104        data(isc+whalo:isc-1,jsc+shalo:jec-nhalo+yshift,:) = 0
105     end if
107     if(shalo>=0) then
108        data(isd:ied+xshift, jec+nhalo+1+yshift:jed+yshift,:) = 0
109        data(isd:ied+xshift, jsd:jsc-shalo-1,:) = 0
110     else
111        data(isc+whalo:iec-ehalo+xshift,jec+1+yshift:jec-nhalo+yshift,:) = 0
112        data(isc+whalo:iec-ehalo+xshift,jsc+shalo:jsc-1,:) = 0
113     end if
114   end subroutine fill_halo_zero_r8
116   !> fill the halo region of a 32-bit real array with zeros
117   subroutine fill_halo_zero_r4(data, whalo, ehalo, shalo, nhalo, xshift, yshift, isc, iec, jsc, jec, isd, ied, &
118                               &  jsd, jed)
119     real(kind=r4_kind), dimension(isd:,jsd:,:), intent(inout) :: data
120     integer,                         intent(in) :: isc, iec, jsc, jec, isd, ied, jsd, jed
121     integer,                         intent(in) :: whalo, ehalo, shalo, nhalo, xshift, yshift
123     if(whalo >=0) then
124        data(iec+ehalo+1+xshift:ied+xshift,jsd:jed+yshift,:) = 0
125        data(isd:isc-whalo-1,jsd:jed+yshift,:) = 0
126     else
127        data(iec+1+xshift:iec-ehalo+xshift,jsc+shalo:jec-nhalo+yshift,:) = 0
128        data(isc+whalo:isc-1,jsc+shalo:jec-nhalo+yshift,:) = 0
129     end if
131     if(shalo>=0) then
132        data(isd:ied+xshift, jec+nhalo+1+yshift:jed+yshift,:) = 0
133        data(isd:ied+xshift, jsd:jsc-shalo-1,:) = 0
134     else
135        data(isc+whalo:iec-ehalo+xshift,jec+1+yshift:jec-nhalo+yshift,:) = 0
136        data(isc+whalo:iec-ehalo+xshift,jsc+shalo:jsc-1,:) = 0
137     end if
138   end subroutine fill_halo_zero_r4
140 !> fill the halo region of a 64-bit integer array with zeros
141   subroutine fill_halo_zero_i8(data, whalo, ehalo, shalo, nhalo, xshift, yshift, isc, iec, jsc, jec, isd, ied, &
142                               &  jsd, jed)
143     integer(kind=i8_kind), dimension(isd:,jsd:,:), intent(inout) :: data
144     integer,                         intent(in) :: isc, iec, jsc, jec, isd, ied, jsd, jed
145     integer,                         intent(in) :: whalo, ehalo, shalo, nhalo, xshift, yshift
147     if(whalo >=0) then
148        data(iec+ehalo+1+xshift:ied+xshift,jsd:jed+yshift,:) = 0
149        data(isd:isc-whalo-1,jsd:jed+yshift,:) = 0
150     else
151        data(iec+1+xshift:iec-ehalo+xshift,jsc+shalo:jec-nhalo+yshift,:) = 0
152        data(isc+whalo:isc-1,jsc+shalo:jec-nhalo+yshift,:) = 0
153     end if
155     if(shalo>=0) then
156        data(isd:ied+xshift, jec+nhalo+1+yshift:jed+yshift,:) = 0
157        data(isd:ied+xshift, jsd:jsc-shalo-1,:) = 0
158     else
159        data(isc+whalo:iec-ehalo+xshift,jec+1+yshift:jec-nhalo+yshift,:) = 0
160        data(isc+whalo:iec-ehalo+xshift,jsc+shalo:jsc-1,:) = 0
161     end if
162   end subroutine fill_halo_zero_i8
164 !> fill the halo region of a 32-bit integer array with zeros
165   subroutine fill_halo_zero_i4(data, whalo, ehalo, shalo, nhalo, xshift, yshift, isc, iec, jsc, jec, isd, ied, &
166                               &  jsd, jed)
167     integer(kind=i4_kind), dimension(isd:,jsd:,:), intent(inout) :: data
168     integer,                         intent(in) :: isc, iec, jsc, jec, isd, ied, jsd, jed
169     integer,                         intent(in) :: whalo, ehalo, shalo, nhalo, xshift, yshift
171     if(whalo >=0) then
172        data(iec+ehalo+1+xshift:ied+xshift,jsd:jed+yshift,:) = 0
173        data(isd:isc-whalo-1,jsd:jed+yshift,:) = 0
174     else
175        data(iec+1+xshift:iec-ehalo+xshift,jsc+shalo:jec-nhalo+yshift,:) = 0
176        data(isc+whalo:isc-1,jsc+shalo:jec-nhalo+yshift,:) = 0
177     end if
179     if(shalo>=0) then
180        data(isd:ied+xshift, jec+nhalo+1+yshift:jed+yshift,:) = 0
181        data(isd:ied+xshift, jsd:jsc-shalo-1,:) = 0
182     else
183        data(isc+whalo:iec-ehalo+xshift,jec+1+yshift:jec-nhalo+yshift,:) = 0
184        data(isc+whalo:iec-ehalo+xshift,jsc+shalo:jsc-1,:) = 0
185     end if
186   end subroutine fill_halo_zero_i4
189   !> fill the halo region of 64-bit array on a regular grid
190   subroutine fill_regular_refinement_halo_r8( data, data_all, ni, nj, tm, te, tse, ts, &
191                                              tsw, tw, tnw, tn, tne, ioff, joff )
192     real(kind=r8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
193     real(kind=r8_kind), dimension(:,:,:,:),             intent(in)    :: data_all
194     integer, dimension(:),                intent(in)    :: ni, nj
195     integer,                              intent(in)    :: tm, te, tse, ts, tsw, tw, tnw, tn, tne
196     integer,                              intent(in)    :: ioff, joff
199     if(te>0) data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, 1:nj(tm)+joff,                   :) = &
200              data_all(1+ioff:ehalo+ioff,               1:nj(te)+joff,                   :,te)  ! east
201     if(ts>0) data    (1:ni(tm)+ioff,                   1-shalo:0,                       :) = &
202              data_all(1:ni(ts)+ioff,                   nj(ts)-shalo+1:nj(ts),           :,ts)  ! south
203     if(tw>0) data    (1-whalo:0,                       1:nj(tm)+joff,                   :) = &
204              data_all(ni(tw)-whalo+1:ni(tw),           1:nj(tw)+joff,                   :,tw)  ! west
205     if(tn>0) data    (1:ni(tm)+ioff,                   nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
206              data_all(1:ni(tn)+ioff,                   1+joff:nhalo+joff,               :,tn)  ! north
207     if(tse>0)data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, 1-shalo:0,                       :) = &
208              data_all(1+ioff:ehalo+ioff,               nj(tse)-shalo+1:nj(tse),         :,tse) ! southeast
209     if(tsw>0)data    (1-whalo:0,                       1-shalo:0,                       :) = &
210              data_all(ni(tsw)-whalo+1:ni(tsw),         nj(tsw)-shalo+1:nj(tsw),         :,tsw) ! southwest
211     if(tne>0)data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
212              data_all(1+ioff:ehalo+ioff,               1+joff:nhalo+joff,               :,tnw) ! northeast
213     if(tnw>0)data    (1-whalo:0,                       nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
214              data_all(ni(tnw)-whalo+1:ni(tnw),         1+joff:nhalo+joff,               :,tne) ! northwest
216   end subroutine fill_regular_refinement_halo_r8
218   !> fill the halo region of 32-bit array on a regular grid
219   subroutine fill_regular_refinement_halo_r4( data, data_all, ni, nj, tm, te, tse, ts, tsw, tw, tnw, tn, tne, &
220                                             &  ioff, joff )
221     real(kind=r4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
222     real(kind=r4_kind), dimension(:,:,:,:),             intent(in)    :: data_all
223     integer, dimension(:),                intent(in)    :: ni, nj
224     integer,                              intent(in)    :: tm, te, tse, ts, tsw, tw, tnw, tn, tne
225     integer,                              intent(in)    :: ioff, joff
228     if(te>0) data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, 1:nj(tm)+joff,                   :) = &
229              data_all(1+ioff:ehalo+ioff,               1:nj(te)+joff,                   :,te)  ! east
230     if(ts>0) data    (1:ni(tm)+ioff,                   1-shalo:0,                       :) = &
231              data_all(1:ni(ts)+ioff,                   nj(ts)-shalo+1:nj(ts),           :,ts)  ! south
232     if(tw>0) data    (1-whalo:0,                       1:nj(tm)+joff,                   :) = &
233              data_all(ni(tw)-whalo+1:ni(tw),           1:nj(tw)+joff,                   :,tw)  ! west
234     if(tn>0) data    (1:ni(tm)+ioff,                   nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
235              data_all(1:ni(tn)+ioff,                   1+joff:nhalo+joff,               :,tn)  ! north
236     if(tse>0)data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, 1-shalo:0,                       :) = &
237              data_all(1+ioff:ehalo+ioff,               nj(tse)-shalo+1:nj(tse),         :,tse) ! southeast
238     if(tsw>0)data    (1-whalo:0,                       1-shalo:0,                       :) = &
239              data_all(ni(tsw)-whalo+1:ni(tsw),         nj(tsw)-shalo+1:nj(tsw),         :,tsw) ! southwest
240     if(tne>0)data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
241              data_all(1+ioff:ehalo+ioff,               1+joff:nhalo+joff,               :,tnw) ! northeast
242     if(tnw>0)data    (1-whalo:0,                       nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
243              data_all(ni(tnw)-whalo+1:ni(tnw),         1+joff:nhalo+joff,               :,tne) ! northwest
245   end subroutine fill_regular_refinement_halo_r4
247 !> fill the halo region of 64-bit integer array on a regular grid
248   subroutine fill_regular_refinement_halo_i8( data, data_all, ni, nj, tm, te, tse, ts, tsw, &
249                                              tw, tnw, tn, tne, ioff, joff )
250     integer(kind=i8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
251     integer(kind=i8_kind), dimension(:,:,:,:),             intent(in)    :: data_all
252     integer, dimension(:),                intent(in)    :: ni, nj
253     integer,                              intent(in)    :: tm, te, tse, ts, tsw, tw, tnw, tn, tne
254     integer,                              intent(in)    :: ioff, joff
257     if(te>0) data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, 1:nj(tm)+joff,                   :) = &
258              data_all(1+ioff:ehalo+ioff,               1:nj(te)+joff,                   :,te)  ! east
259     if(ts>0) data    (1:ni(tm)+ioff,                   1-shalo:0,                       :) = &
260              data_all(1:ni(ts)+ioff,                   nj(ts)-shalo+1:nj(ts),           :,ts)  ! south
261     if(tw>0) data    (1-whalo:0,                       1:nj(tm)+joff,                   :) = &
262              data_all(ni(tw)-whalo+1:ni(tw),           1:nj(tw)+joff,                   :,tw)  ! west
263     if(tn>0) data    (1:ni(tm)+ioff,                   nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
264              data_all(1:ni(tn)+ioff,                   1+joff:nhalo+joff,               :,tn)  ! north
265     if(tse>0)data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, 1-shalo:0,                       :) = &
266              data_all(1+ioff:ehalo+ioff,               nj(tse)-shalo+1:nj(tse),         :,tse) ! southeast
267     if(tsw>0)data    (1-whalo:0,                       1-shalo:0,                       :) = &
268              data_all(ni(tsw)-whalo+1:ni(tsw),         nj(tsw)-shalo+1:nj(tsw),         :,tsw) ! southwest
269     if(tne>0)data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
270              data_all(1+ioff:ehalo+ioff,               1+joff:nhalo+joff,               :,tnw) ! northeast
271     if(tnw>0)data    (1-whalo:0,                       nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
272              data_all(ni(tnw)-whalo+1:ni(tnw),         1+joff:nhalo+joff,               :,tne) ! northwest
274   end subroutine fill_regular_refinement_halo_i8
276 !> fill the halo region of 32-bit integer array on a regular grid
277   subroutine fill_regular_refinement_halo_i4( data, data_all, ni, nj, tm, te, tse, ts, tsw, tw, tnw, tn, tne, &
278                                             &  ioff, joff )
279     integer(kind=i4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
280     integer(kind=i4_kind), dimension(:,:,:,:),             intent(in)    :: data_all
281     integer, dimension(:),                intent(in)    :: ni, nj
282     integer,                              intent(in)    :: tm, te, tse, ts, tsw, tw, tnw, tn, tne
283     integer,                              intent(in)    :: ioff, joff
286     if(te>0) data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, 1:nj(tm)+joff,                   :) = &
287              data_all(1+ioff:ehalo+ioff,               1:nj(te)+joff,                   :,te)  ! east
288     if(ts>0) data    (1:ni(tm)+ioff,                   1-shalo:0,                       :) = &
289              data_all(1:ni(ts)+ioff,                   nj(ts)-shalo+1:nj(ts),           :,ts)  ! south
290     if(tw>0) data    (1-whalo:0,                       1:nj(tm)+joff,                   :) = &
291              data_all(ni(tw)-whalo+1:ni(tw),           1:nj(tw)+joff,                   :,tw)  ! west
292     if(tn>0) data    (1:ni(tm)+ioff,                   nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
293              data_all(1:ni(tn)+ioff,                   1+joff:nhalo+joff,               :,tn)  ! north
294     if(tse>0)data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, 1-shalo:0,                       :) = &
295              data_all(1+ioff:ehalo+ioff,               nj(tse)-shalo+1:nj(tse),         :,tse) ! southeast
296     if(tsw>0)data    (1-whalo:0,                       1-shalo:0,                       :) = &
297              data_all(ni(tsw)-whalo+1:ni(tsw),         nj(tsw)-shalo+1:nj(tsw),         :,tsw) ! southwest
298     if(tne>0)data    (ni(tm)+1+ioff:ni(tm)+ehalo+ioff, nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
299              data_all(1+ioff:ehalo+ioff,               1+joff:nhalo+joff,               :,tnw) ! northeast
300     if(tnw>0)data    (1-whalo:0,                       nj(tm)+1+joff:nj(tm)+nhalo+joff, :) = &
301              data_all(ni(tnw)-whalo+1:ni(tnw),         1+joff:nhalo+joff,               :,tne) ! northwest
303   end subroutine fill_regular_refinement_halo_i4
305   ! Fill the halo points of a 64-bit real array on the regular mosaic grid
306   subroutine fill_regular_mosaic_halo_r8(data, data_all, te, tse, ts, tsw, tw, tnw, tn, tne)
307     real(kind=r8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
308     real(kind=r8_kind), dimension(:,:,:,:),             intent(in)    :: data_all
309     integer, intent(in)    :: te, tse, ts, tsw, tw, tnw, tn, tne
311     data(nx+1:nx+ehalo, 1:ny,          :) = data_all(1:ehalo,       1:ny,          :, te) ! east
312     data(1:nx,          1-shalo:0,     :) = data_all(1:nx,          ny-shalo+1:ny, :, ts) ! south
313     data(1-whalo:0,     1:ny,          :) = data_all(nx-whalo+1:nx, 1:ny,          :, tw) ! west
314     data(1:nx,          ny+1:ny+nhalo, :) = data_all(1:nx,          1:nhalo,       :, tn) ! north
315     data(nx+1:nx+ehalo, 1-shalo:0,     :) = data_all(1:ehalo,       ny-shalo+1:ny, :,tse) ! southeast
316     data(1-whalo:0,     1-shalo:0,     :) = data_all(nx-whalo+1:nx, ny-shalo+1:ny, :,tsw) ! southwest
317     data(nx+1:nx+ehalo, ny+1:ny+nhalo, :) = data_all(1:ehalo,       1:nhalo,       :,tnw) ! northeast
318     data(1-whalo:0,     ny+1:ny+nhalo, :) = data_all(nx-whalo+1:nx, 1:nhalo,       :,tne) ! northwest
319   end subroutine fill_regular_mosaic_halo_r8
321   !> Fill the halo points of a 32-bit real array on the regular mosaic grid
322   subroutine fill_regular_mosaic_halo_r4(data, data_all, te, tse, ts, tsw, tw, tnw, tn, tne)
323     real(kind=r4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
324     real(kind=r4_kind), dimension(:,:,:,:),             intent(in)    :: data_all
325     integer, intent(in)    :: te, tse, ts, tsw, tw, tnw, tn, tne
327     data(nx+1:nx+ehalo, 1:ny,          :) = data_all(1:ehalo,       1:ny,          :, te) ! east
328     data(1:nx,          1-shalo:0,     :) = data_all(1:nx,          ny-shalo+1:ny, :, ts) ! south
329     data(1-whalo:0,     1:ny,          :) = data_all(nx-whalo+1:nx, 1:ny,          :, tw) ! west
330     data(1:nx,          ny+1:ny+nhalo, :) = data_all(1:nx,          1:nhalo,       :, tn) ! north
331     data(nx+1:nx+ehalo, 1-shalo:0,     :) = data_all(1:ehalo,       ny-shalo+1:ny, :,tse) ! southeast
332     data(1-whalo:0,     1-shalo:0,     :) = data_all(nx-whalo+1:nx, ny-shalo+1:ny, :,tsw) ! southwest
333     data(nx+1:nx+ehalo, ny+1:ny+nhalo, :) = data_all(1:ehalo,       1:nhalo,       :,tnw) ! northeast
334     data(1-whalo:0,     ny+1:ny+nhalo, :) = data_all(nx-whalo+1:nx, 1:nhalo,       :,tne) ! northwest
335   end subroutine fill_regular_mosaic_halo_r4
337   ! Fill the halo points of a 64-bit integer array on the regular mosaic grid
338   subroutine fill_regular_mosaic_halo_i8(data, data_all, te, tse, ts, tsw, tw, tnw, tn, tne)
339     integer(kind=i8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
340     integer(kind=i8_kind), dimension(:,:,:,:),             intent(in)    :: data_all
341     integer, intent(in)    :: te, tse, ts, tsw, tw, tnw, tn, tne
343     data(nx+1:nx+ehalo, 1:ny,          :) = data_all(1:ehalo,       1:ny,          :, te) ! east
344     data(1:nx,          1-shalo:0,     :) = data_all(1:nx,          ny-shalo+1:ny, :, ts) ! south
345     data(1-whalo:0,     1:ny,          :) = data_all(nx-whalo+1:nx, 1:ny,          :, tw) ! west
346     data(1:nx,          ny+1:ny+nhalo, :) = data_all(1:nx,          1:nhalo,       :, tn) ! north
347     data(nx+1:nx+ehalo, 1-shalo:0,     :) = data_all(1:ehalo,       ny-shalo+1:ny, :,tse) ! southeast
348     data(1-whalo:0,     1-shalo:0,     :) = data_all(nx-whalo+1:nx, ny-shalo+1:ny, :,tsw) ! southwest
349     data(nx+1:nx+ehalo, ny+1:ny+nhalo, :) = data_all(1:ehalo,       1:nhalo,       :,tnw) ! northeast
350     data(1-whalo:0,     ny+1:ny+nhalo, :) = data_all(nx-whalo+1:nx, 1:nhalo,       :,tne) ! northwest
351   end subroutine fill_regular_mosaic_halo_i8
353   !> Fill the halo points of a 64-bit integer array on the regular mosaic grid
354   subroutine fill_regular_mosaic_halo_i4(data, data_all, te, tse, ts, tsw, tw, tnw, tn, tne)
355     integer(kind=i4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
356     integer(kind=i4_kind), dimension(:,:,:,:),             intent(in)    :: data_all
357     integer, intent(in)    :: te, tse, ts, tsw, tw, tnw, tn, tne
359     data(nx+1:nx+ehalo, 1:ny,          :) = data_all(1:ehalo,       1:ny,          :, te) ! east
360     data(1:nx,          1-shalo:0,     :) = data_all(1:nx,          ny-shalo+1:ny, :, ts) ! south
361     data(1-whalo:0,     1:ny,          :) = data_all(nx-whalo+1:nx, 1:ny,          :, tw) ! west
362     data(1:nx,          ny+1:ny+nhalo, :) = data_all(1:nx,          1:nhalo,       :, tn) ! north
363     data(nx+1:nx+ehalo, 1-shalo:0,     :) = data_all(1:ehalo,       ny-shalo+1:ny, :,tse) ! southeast
364     data(1-whalo:0,     1-shalo:0,     :) = data_all(nx-whalo+1:nx, ny-shalo+1:ny, :,tsw) ! southwest
365     data(nx+1:nx+ehalo, ny+1:ny+nhalo, :) = data_all(1:ehalo,       1:nhalo,       :,tnw) ! northeast
366     data(1-whalo:0,     ny+1:ny+nhalo, :) = data_all(nx-whalo+1:nx, 1:nhalo,       :,tne) ! northwest
367   end subroutine fill_regular_mosaic_halo_i4
369   !> Fill the halo region of a 64-bit array real on a domain with a folded north edge
370   subroutine fill_folded_north_halo_r8(data, ioff, joff, ishift, jshift, sign)
371     real(kind=r8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
372     integer, intent(in) :: ioff, joff, ishift, jshift, sign
373     ! local
374     integer :: nxp, nyp, m1, m2
376     nxp = nx+ishift
377     nyp = ny+jshift
378     m1 = ishift - ioff
379     m2 = 2*ishift - ioff
381     data(1-whalo:0,1:nyp,:) = data(nx-whalo+1:nx,1:ny+jshift,:) ! west
382     data(nx+1:nx+ehalo+ishift,1:nyp,:) = data(1:ehalo+ishift,1:ny+jshift,:) ! east
383     if(m1 .GE. 1-whalo) &
384       data(1-whalo:m1,nyp+1:nyp+nhalo,:) = sign*data(whalo+m2:1+ishift:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
385     data(m1+1:nx+m2,nyp+1:nyp+nhalo,:) = sign*data(nx+ishift:1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
386     data(nx+m2+1:nxp+ehalo,nyp+1:nyp+nhalo,:) = sign*data(nx:nx-ehalo+m1+1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
388   end subroutine fill_folded_north_halo_r8
390   !> Fill the halo region of a 32-bit real array on a domain with a folded north edge
391   subroutine fill_folded_north_halo_r4(data, ioff, joff, ishift, jshift, sign)
392     real(kind=r4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
393     integer, intent(in) :: ioff, joff, ishift, jshift, sign
394     ! local
395     integer :: nxp, nyp, m1, m2
397     nxp = nx+ishift
398     nyp = ny+jshift
399     m1 = ishift - ioff
400     m2 = 2*ishift - ioff
402     data(1-whalo:0,1:nyp,:) = data(nx-whalo+1:nx,1:ny+jshift,:) ! west
403     data(nx+1:nx+ehalo+ishift,1:nyp,:) = data(1:ehalo+ishift,1:ny+jshift,:) ! east
405     if(m1 .GE. 1-whalo) &
406       data(1-whalo:m1,nyp+1:nyp+nhalo,:) = sign*data(whalo+m2:1+ishift:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
407     data(m1+1:nx+m2,nyp+1:nyp+nhalo,:) = sign*data(nx+ishift:1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
408     data(nx+m2+1:nxp+ehalo,nyp+1:nyp+nhalo,:) = sign*data(nx:nx-ehalo+m1+1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
410   end subroutine fill_folded_north_halo_r4
412   !> Fill the halo region of a 64-bit integer array on a domain with a folded north edge
413   subroutine fill_folded_north_halo_i8(data, ioff, joff, ishift, jshift, sign)
414     integer(kind=i8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
415     integer, intent(in) :: ioff, joff, ishift, jshift, sign
416     ! local
417     integer :: nxp, nyp, m1, m2
419     nxp = nx+ishift
420     nyp = ny+jshift
421     m1 = ishift - ioff
422     m2 = 2*ishift - ioff
424     data(1-whalo:0,1:nyp,:) = data(nx-whalo+1:nx,1:ny+jshift,:) ! west
425     data(nx+1:nx+ehalo+ishift,1:nyp,:) = data(1:ehalo+ishift,1:ny+jshift,:) ! east
426     if(m1 .GE. 1-whalo) &
427       data(1-whalo:m1,nyp+1:nyp+nhalo,:) = sign*data(whalo+m2:1+ishift:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
428     data(m1+1:nx+m2,nyp+1:nyp+nhalo,:) = sign*data(nx+ishift:1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
429     data(nx+m2+1:nxp+ehalo,nyp+1:nyp+nhalo,:) = sign*data(nx:nx-ehalo+m1+1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
431   end subroutine fill_folded_north_halo_i8
433   !> Fill the halo region of a 32-bit integer array on a domain with a folded north edge
434   subroutine fill_folded_north_halo_i4(data, ioff, joff, ishift, jshift, sign)
435     integer(kind=i4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
436     integer, intent(in) :: ioff, joff, ishift, jshift, sign
437     ! local
438     integer :: nxp, nyp, m1, m2
440     nxp = nx+ishift
441     nyp = ny+jshift
442     m1 = ishift - ioff
443     m2 = 2*ishift - ioff
445     data(1-whalo:0,1:nyp,:) = data(nx-whalo+1:nx,1:ny+jshift,:) ! west
446     data(nx+1:nx+ehalo+ishift,1:nyp,:) = data(1:ehalo+ishift,1:ny+jshift,:) ! east
448     if(m1 .GE. 1-whalo) &
449       data(1-whalo:m1,nyp+1:nyp+nhalo,:) = sign*data(whalo+m2:1+ishift:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
450     data(m1+1:nx+m2,nyp+1:nyp+nhalo,:) = sign*data(nx+ishift:1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
451     data(nx+m2+1:nxp+ehalo,nyp+1:nyp+nhalo,:) = sign*data(nx:nx-ehalo+m1+1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
453   end subroutine fill_folded_north_halo_i4
455   !> Fill the halo region of a 64-bit real array on a domain with a folded south edge
456   subroutine fill_folded_south_halo_r8(data, ioff, joff, ishift, jshift, sign)
457     real(kind=r8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
458     integer, intent(in) :: ioff, joff, ishift, jshift, sign
459     ! local
460     integer  :: nxp, nyp, m1, m2
462     nxp = nx+ishift
463     nyp = ny+jshift
464     m1 = ishift - ioff
465     m2 = 2*ishift - ioff
467     data(1-whalo:0,1:nyp,:) = data(nx-whalo+1:nx,1:nyp,:) ! west
468     data(nx+1:nx+ehalo+ishift,1:nyp,:) = data(1:ehalo+ishift, 1:nyp,:) ! east
469     if(m1 .GE. 1-whalo) &
470       data(1-whalo:m1,1-shalo:0,:) = sign*data(whalo+m2:1+ishift:-1, shalo+jshift:1+jshift:-1,:)
472     data(m1+1:nx+m2,1-shalo:0,:) = sign*data(nxp:1:-1,shalo+jshift:1+jshift:-1,:)
473     data(nx+m2+1:nxp+ehalo,1-shalo:0,:) = sign*data(nx:nx-ehalo+m1+1:-1,shalo+jshift:1+jshift:-1,:)
475   end subroutine fill_folded_south_halo_r8
477   !> Fill the halo region of a 32-bit real array on a domain with a folded south edge
478   subroutine fill_folded_south_halo_r4(data, ioff, joff, ishift, jshift, sign)
479     real(kind=r4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
480     integer, intent(in) :: ioff, joff, ishift, jshift, sign
481     ! local
482     integer  :: nxp, nyp, m1, m2
484     nxp = nx+ishift
485     nyp = ny+jshift
486     m1 = ishift - ioff
487     m2 = 2*ishift - ioff
489     data(1-whalo:0,1:nyp,:) = data(nx-whalo+1:nx,1:nyp,:) ! west
490     data(nx+1:nx+ehalo+ishift,1:nyp,:) = data(1:ehalo+ishift, 1:nyp,:) ! east
491     if(m1 .GE. 1-whalo) &
492       data(1-whalo:m1,1-shalo:0,:) = sign*data(whalo+m2:1+ishift:-1, shalo+jshift:1+jshift:-1,:)
494     data(m1+1:nx+m2,1-shalo:0,:) = sign*data(nxp:1:-1,shalo+jshift:1+jshift:-1,:)
495     data(nx+m2+1:nxp+ehalo,1-shalo:0,:) = sign*data(nx:nx-ehalo+m1+1:-1,shalo+jshift:1+jshift:-1,:)
497   end subroutine fill_folded_south_halo_r4
499   !> Fill the halo region of a 64-bit intger array on a domain with a folded south edge
500   subroutine fill_folded_south_halo_i8(data, ioff, joff, ishift, jshift, sign)
501     integer(kind=i8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
502     integer, intent(in) :: ioff, joff, ishift, jshift, sign
503     ! local
504     integer  :: nxp, nyp, m1, m2
506     nxp = nx+ishift
507     nyp = ny+jshift
508     m1 = ishift - ioff
509     m2 = 2*ishift - ioff
511     data(1-whalo:0,1:nyp,:) = data(nx-whalo+1:nx,1:nyp,:) ! west
512     data(nx+1:nx+ehalo+ishift,1:nyp,:) = data(1:ehalo+ishift, 1:nyp,:) ! east
513     if(m1 .GE. 1-whalo) &
514       data(1-whalo:m1,1-shalo:0,:) = sign*data(whalo+m2:1+ishift:-1, shalo+jshift:1+jshift:-1,:)
516     data(m1+1:nx+m2,1-shalo:0,:) = sign*data(nxp:1:-1,shalo+jshift:1+jshift:-1,:)
517     data(nx+m2+1:nxp+ehalo,1-shalo:0,:) = sign*data(nx:nx-ehalo+m1+1:-1,shalo+jshift:1+jshift:-1,:)
519   end subroutine fill_folded_south_halo_i8
521   !> Fill the halo region of a 32-bit integer array on a domain with a folded south edge
522   subroutine fill_folded_south_halo_i4(data, ioff, joff, ishift, jshift, sign)
523     integer(kind=i4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
524     integer, intent(in) :: ioff, joff, ishift, jshift, sign
525     ! local
526     integer  :: nxp, nyp, m1, m2
528     nxp = nx+ishift
529     nyp = ny+jshift
530     m1 = ishift - ioff
531     m2 = 2*ishift - ioff
533     data(1-whalo:0,1:nyp,:) = data(nx-whalo+1:nx,1:nyp,:) ! west
534     data(nx+1:nx+ehalo+ishift,1:nyp,:) = data(1:ehalo+ishift, 1:nyp,:) ! east
535     if(m1 .GE. 1-whalo) &
536       data(1-whalo:m1,1-shalo:0,:) = sign*data(whalo+m2:1+ishift:-1, shalo+jshift:1+jshift:-1,:)
538     data(m1+1:nx+m2,1-shalo:0,:) = sign*data(nxp:1:-1,shalo+jshift:1+jshift:-1,:)
539     data(nx+m2+1:nxp+ehalo,1-shalo:0,:) = sign*data(nx:nx-ehalo+m1+1:-1,shalo+jshift:1+jshift:-1,:)
541   end subroutine fill_folded_south_halo_i4
543   !> Fill the halo region of a 64-bit real array on a domain with a folded west edge
544   subroutine fill_folded_west_halo_r8(data, ioff, joff, ishift, jshift, sign)
545     real(kind=r8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
546     integer, intent(in) :: ioff, joff, ishift, jshift, sign
547     ! local
548     integer :: nxp, nyp, m1, m2
550     nxp = nx+ishift
551     nyp = ny+jshift
552     m1 = jshift - joff
553     m2 = 2*jshift - joff
555     data(1:nxp, 1-shalo:0,:) = data(1:nxp, ny-shalo+1:ny, :) ! south
556     data(1:nxp, ny+1:nyp+nhalo,:) = data(1:nxp, 1:nhalo+jshift,:) ! north
557     if(m1 .GE. 1-shalo) &
558       data(1-whalo:0, 1-shalo:m1,:) = sign*data(whalo+ishift:1+ishift:-1,shalo+m2:1+jshift:-1,:)
559     data(1-whalo:0, m1+1:ny+m2,:) = sign*data(whalo+ishift:1+ishift:-1, nyp:1:-1, :)
560     data(1-whalo:0, ny+m2+1:nyp+nhalo,:) = sign*data(whalo+ishift:1+ishift:-1, ny:ny-nhalo+m1+1:-1,:)
562   end subroutine fill_folded_west_halo_r8
564   !> Fill the halo region of a 32-bit real array on a domain with a folded west edge
565   subroutine fill_folded_west_halo_r4(data, ioff, joff, ishift, jshift, sign)
566     real(kind=r4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
567     integer, intent(in) :: ioff, joff, ishift, jshift, sign
568     ! local
569     integer :: nxp, nyp, m1, m2
571     nxp = nx+ishift
572     nyp = ny+jshift
573     m1 = jshift - joff
574     m2 = 2*jshift - joff
576     data(1:nxp, 1-shalo:0,:) = data(1:nxp, ny-shalo+1:ny, :) ! south
577     data(1:nxp, ny+1:nyp+nhalo,:) = data(1:nxp, 1:nhalo+jshift,:) ! north
578     if(m1 .GE. 1-shalo) &
579       data(1-whalo:0, 1-shalo:m1,:) = sign*data(whalo+ishift:1+ishift:-1,shalo+m2:1+jshift:-1,:)
580     data(1-whalo:0, m1+1:ny+m2,:) = sign*data(whalo+ishift:1+ishift:-1, nyp:1:-1, :)
581     data(1-whalo:0, ny+m2+1:nyp+nhalo,:) = sign*data(whalo+ishift:1+ishift:-1, ny:ny-nhalo+m1+1:-1,:)
583   end subroutine fill_folded_west_halo_r4
585   !> Fill the halo region of a 64-bit integer array on a domain with a folded west edge
586   subroutine fill_folded_west_halo_i8(data, ioff, joff, ishift, jshift, sign)
587     integer(kind=i8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
588     integer, intent(in) :: ioff, joff, ishift, jshift, sign
589     ! local
590     integer :: nxp, nyp, m1, m2
592     nxp = nx+ishift
593     nyp = ny+jshift
594     m1 = jshift - joff
595     m2 = 2*jshift - joff
597     data(1:nxp, 1-shalo:0,:) = data(1:nxp, ny-shalo+1:ny, :) ! south
598     data(1:nxp, ny+1:nyp+nhalo,:) = data(1:nxp, 1:nhalo+jshift,:) ! north
599     if(m1 .GE. 1-shalo) &
600       data(1-whalo:0, 1-shalo:m1,:) = sign*data(whalo+ishift:1+ishift:-1,shalo+m2:1+jshift:-1,:)
601     data(1-whalo:0, m1+1:ny+m2,:) = sign*data(whalo+ishift:1+ishift:-1, nyp:1:-1, :)
602     data(1-whalo:0, ny+m2+1:nyp+nhalo,:) = sign*data(whalo+ishift:1+ishift:-1, ny:ny-nhalo+m1+1:-1,:)
604   end subroutine fill_folded_west_halo_i8
606   !> Fill the halo region of a 32-bit integer array on a domain with a folded west edge
607   subroutine fill_folded_west_halo_i4(data, ioff, joff, ishift, jshift, sign)
608     integer(kind=i4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
609     integer, intent(in) :: ioff, joff, ishift, jshift, sign
610     ! local
611     integer :: nxp, nyp, m1, m2
613     nxp = nx+ishift
614     nyp = ny+jshift
615     m1 = jshift - joff
616     m2 = 2*jshift - joff
618     data(1:nxp, 1-shalo:0,:) = data(1:nxp, ny-shalo+1:ny, :) ! south
619     data(1:nxp, ny+1:nyp+nhalo,:) = data(1:nxp, 1:nhalo+jshift,:) ! north
620     if(m1 .GE. 1-shalo) &
621       data(1-whalo:0, 1-shalo:m1,:) = sign*data(whalo+ishift:1+ishift:-1,shalo+m2:1+jshift:-1,:)
622     data(1-whalo:0, m1+1:ny+m2,:) = sign*data(whalo+ishift:1+ishift:-1, nyp:1:-1, :)
623     data(1-whalo:0, ny+m2+1:nyp+nhalo,:) = sign*data(whalo+ishift:1+ishift:-1, ny:ny-nhalo+m1+1:-1,:)
625   end subroutine fill_folded_west_halo_i4
627   !> Fill the halo region of a 64-bit real array on a domain with a folded east edge
628   subroutine fill_folded_east_halo_r8(data, ioff, joff, ishift, jshift, sign)
629     real(kind=r8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
630     integer, intent(in) :: ioff, joff, ishift, jshift, sign
631     ! local
632     integer :: nxp, nyp, m1, m2
634     nxp = nx+ishift
635     nyp = ny+jshift
636     m1 = jshift - joff
637     m2 = 2*jshift - joff
639     data(1:nxp, 1-shalo:0, :)  = data(1:nxp, ny-shalo+1:ny, :) ! south
640     data(1:nxp, ny+1:nyp+nhalo, :) = data(1:nxp, 1:nhalo+jshift,:) ! north
641     if(m1 .GE. 1-shalo) &
642       data(nxp+1:nxp+ehalo, 1-shalo:m1, :) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, shalo+m2:1+jshift:-1,:)
644     data(nxp+1:nxp+ehalo, m1+1:ny+m2, :) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, nyp:1:-1, :)
645     data(nxp+1:nxp+ehalo, ny+m2+1:nyp+nhalo,:) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, ny:ny-nhalo+m1+1:-1,:)
647   end subroutine fill_folded_east_halo_r8
649   !> Fill the halo region of a 32-bit real array on a domain with a folded east edge
650   subroutine fill_folded_east_halo_r4(data, ioff, joff, ishift, jshift, sign)
651     real(kind=r4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
652     integer, intent(in) :: ioff, joff, ishift, jshift, sign
653     ! local
654     integer :: nxp, nyp, m1, m2
656     nxp = nx+ishift
657     nyp = ny+jshift
658     m1 = jshift - joff
659     m2 = 2*jshift - joff
661     data(1:nxp, 1-shalo:0, :)  = data(1:nxp, ny-shalo+1:ny, :) ! south
662     data(1:nxp, ny+1:nyp+nhalo, :) = data(1:nxp, 1:nhalo+jshift,:) ! north
663     if(m1 .GE. 1-shalo) &
664       data(nxp+1:nxp+ehalo, 1-shalo:m1, :) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, shalo+m2:1+jshift:-1,:)
666     data(nxp+1:nxp+ehalo, m1+1:ny+m2, :) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, nyp:1:-1, :)
667     data(nxp+1:nxp+ehalo, ny+m2+1:nyp+nhalo,:) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, ny:ny-nhalo+m1+1:-1,:)
669   end subroutine fill_folded_east_halo_r4
671   !> Fill the halo region of a 64-bit integer array on a domain with a folded east edge
672   subroutine fill_folded_east_halo_i8(data, ioff, joff, ishift, jshift, sign)
673     integer(kind=i8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
674     integer, intent(in) :: ioff, joff, ishift, jshift, sign
675     ! local
676     integer :: nxp, nyp, m1, m2
678     nxp = nx+ishift
679     nyp = ny+jshift
680     m1 = jshift - joff
681     m2 = 2*jshift - joff
683     data(1:nxp, 1-shalo:0, :)  = data(1:nxp, ny-shalo+1:ny, :) ! south
684     data(1:nxp, ny+1:nyp+nhalo, :) = data(1:nxp, 1:nhalo+jshift,:) ! north
685     if(m1 .GE. 1-shalo) &
686       data(nxp+1:nxp+ehalo, 1-shalo:m1, :) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, shalo+m2:1+jshift:-1,:)
688     data(nxp+1:nxp+ehalo, m1+1:ny+m2, :) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, nyp:1:-1, :)
689     data(nxp+1:nxp+ehalo, ny+m2+1:nyp+nhalo,:) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, ny:ny-nhalo+m1+1:-1,:)
691   end subroutine fill_folded_east_halo_i8
693   !> Fill the halo region of a 32-bit integer array on a domain with a folded east edge
694   subroutine fill_folded_east_halo_i4(data, ioff, joff, ishift, jshift, sign)
695     integer(kind=i4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: data
696     integer, intent(in) :: ioff, joff, ishift, jshift, sign
697     ! local
698     integer :: nxp, nyp, m1, m2
700     nxp = nx+ishift
701     nyp = ny+jshift
702     m1 = jshift - joff
703     m2 = 2*jshift - joff
705     data(1:nxp, 1-shalo:0, :)  = data(1:nxp, ny-shalo+1:ny, :) ! south
706     data(1:nxp, ny+1:nyp+nhalo, :) = data(1:nxp, 1:nhalo+jshift,:) ! north
707     if(m1 .GE. 1-shalo) &
708       data(nxp+1:nxp+ehalo, 1-shalo:m1, :) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, shalo+m2:1+jshift:-1,:)
710     data(nxp+1:nxp+ehalo, m1+1:ny+m2, :) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, nyp:1:-1, :)
711     data(nxp+1:nxp+ehalo, ny+m2+1:nyp+nhalo,:) = sign*data(nxp-ioff:nxp-ehalo-ioff+1:-1, ny:ny-nhalo+m1+1:-1,:)
713   end subroutine fill_folded_east_halo_i4
715 end module fill_halo