fix: mosaic2 tests with strict debug flags (#1597)
[FMS.git] / test_fms / mpp / fill_halo.F90
blob63013ddbf807a9520a7d8a2c9d5816336253e792
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(halo_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) :: halo_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        halo_data(iec+ehalo+1+xshift:ied+xshift,jsd:jed+yshift,:) = 0
101        halo_data(isd:isc-whalo-1,jsd:jed+yshift,:) = 0
102     else
103        halo_data(iec+1+xshift:iec-ehalo+xshift,jsc+shalo:jec-nhalo+yshift,:) = 0
104        halo_data(isc+whalo:isc-1,jsc+shalo:jec-nhalo+yshift,:) = 0
105     end if
107     if(shalo>=0) then
108        halo_data(isd:ied+xshift, jec+nhalo+1+yshift:jed+yshift,:) = 0
109        halo_data(isd:ied+xshift, jsd:jsc-shalo-1,:) = 0
110     else
111        halo_data(isc+whalo:iec-ehalo+xshift,jec+1+yshift:jec-nhalo+yshift,:) = 0
112        halo_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(halo_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) :: halo_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       halo_data(iec+ehalo+1+xshift:ied+xshift,jsd:jed+yshift,:) = 0
125       halo_data(isd:isc-whalo-1,jsd:jed+yshift,:) = 0
126     else
127       halo_data(iec+1+xshift:iec-ehalo+xshift,jsc+shalo:jec-nhalo+yshift,:) = 0
128       halo_data(isc+whalo:isc-1,jsc+shalo:jec-nhalo+yshift,:) = 0
129     end if
131     if(shalo>=0) then
132       halo_data(isd:ied+xshift, jec+nhalo+1+yshift:jed+yshift,:) = 0
133       halo_data(isd:ied+xshift, jsd:jsc-shalo-1,:) = 0
134     else
135       halo_data(isc+whalo:iec-ehalo+xshift,jec+1+yshift:jec-nhalo+yshift,:) = 0
136       halo_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(halo_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) :: halo_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       halo_data(iec+ehalo+1+xshift:ied+xshift,jsd:jed+yshift,:) = 0
149       halo_data(isd:isc-whalo-1,jsd:jed+yshift,:) = 0
150     else
151       halo_data(iec+1+xshift:iec-ehalo+xshift,jsc+shalo:jec-nhalo+yshift,:) = 0
152       halo_data(isc+whalo:isc-1,jsc+shalo:jec-nhalo+yshift,:) = 0
153     end if
155     if(shalo>=0) then
156       halo_data(isd:ied+xshift, jec+nhalo+1+yshift:jed+yshift,:) = 0
157       halo_data(isd:ied+xshift, jsd:jsc-shalo-1,:) = 0
158     else
159       halo_data(isc+whalo:iec-ehalo+xshift,jec+1+yshift:jec-nhalo+yshift,:) = 0
160       halo_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(halo_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) :: halo_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       halo_data(iec+ehalo+1+xshift:ied+xshift,jsd:jed+yshift,:) = 0
173       halo_data(isd:isc-whalo-1,jsd:jed+yshift,:) = 0
174     else
175       halo_data(iec+1+xshift:iec-ehalo+xshift,jsc+shalo:jec-nhalo+yshift,:) = 0
176       halo_data(isc+whalo:isc-1,jsc+shalo:jec-nhalo+yshift,:) = 0
177     end if
179     if(shalo>=0) then
180       halo_data(isd:ied+xshift, jec+nhalo+1+yshift:jed+yshift,:) = 0
181       halo_data(isd:ied+xshift, jsd:jsc-shalo-1,:) = 0
182     else
183       halo_data(isc+whalo:iec-ehalo+xshift,jec+1+yshift:jec-nhalo+yshift,:) = 0
184       halo_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( halo_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) :: halo_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) halo_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) halo_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) halo_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) halo_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)halo_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)halo_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)halo_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)halo_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( halo_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) :: halo_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) halo_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) halo_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) halo_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) halo_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)halo_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)halo_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)halo_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)halo_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( halo_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) :: halo_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) halo_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) halo_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) halo_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) halo_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)halo_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)halo_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)halo_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)halo_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( halo_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) :: halo_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) halo_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) halo_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) halo_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) halo_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)halo_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)halo_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)halo_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)halo_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(halo_data, data_all, te, tse, ts, tsw, tw, tnw, tn, tne)
307     real(kind=r8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
308     real(kind=r8_kind), dimension(:,:,:,:),             intent(in)    :: data_all
309     integer, intent(in)    :: te, tse, ts, tsw, tw, tnw, tn, tne
311     halo_data(nx+1:nx+ehalo, 1:ny,          :) = data_all(1:ehalo,       1:ny,          :, te) ! east
312     halo_data(1:nx,          1-shalo:0,     :) = data_all(1:nx,          ny-shalo+1:ny, :, ts) ! south
313     halo_data(1-whalo:0,     1:ny,          :) = data_all(nx-whalo+1:nx, 1:ny,          :, tw) ! west
314     halo_data(1:nx,          ny+1:ny+nhalo, :) = data_all(1:nx,          1:nhalo,       :, tn) ! north
315     halo_data(nx+1:nx+ehalo, 1-shalo:0,     :) = data_all(1:ehalo,       ny-shalo+1:ny, :,tse) ! southeast
316     halo_data(1-whalo:0,     1-shalo:0,     :) = data_all(nx-whalo+1:nx, ny-shalo+1:ny, :,tsw) ! southwest
317     halo_data(nx+1:nx+ehalo, ny+1:ny+nhalo, :) = data_all(1:ehalo,       1:nhalo,       :,tnw) ! northeast
318     halo_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(halo_data, data_all, te, tse, ts, tsw, tw, tnw, tn, tne)
323     real(kind=r4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
324     real(kind=r4_kind), dimension(:,:,:,:),             intent(in)    :: data_all
325     integer, intent(in)    :: te, tse, ts, tsw, tw, tnw, tn, tne
327     halo_data(nx+1:nx+ehalo, 1:ny,          :) = data_all(1:ehalo,       1:ny,          :, te) ! east
328     halo_data(1:nx,          1-shalo:0,     :) = data_all(1:nx,          ny-shalo+1:ny, :, ts) ! south
329     halo_data(1-whalo:0,     1:ny,          :) = data_all(nx-whalo+1:nx, 1:ny,          :, tw) ! west
330     halo_data(1:nx,          ny+1:ny+nhalo, :) = data_all(1:nx,          1:nhalo,       :, tn) ! north
331     halo_data(nx+1:nx+ehalo, 1-shalo:0,     :) = data_all(1:ehalo,       ny-shalo+1:ny, :,tse) ! southeast
332     halo_data(1-whalo:0,     1-shalo:0,     :) = data_all(nx-whalo+1:nx, ny-shalo+1:ny, :,tsw) ! southwest
333     halo_data(nx+1:nx+ehalo, ny+1:ny+nhalo, :) = data_all(1:ehalo,       1:nhalo,       :,tnw) ! northeast
334     halo_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(halo_data, data_all, te, tse, ts, tsw, tw, tnw, tn, tne)
339     integer(kind=i8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
340     integer(kind=i8_kind), dimension(:,:,:,:),             intent(in)    :: data_all
341     integer, intent(in)    :: te, tse, ts, tsw, tw, tnw, tn, tne
343     halo_data(nx+1:nx+ehalo, 1:ny,          :) = data_all(1:ehalo,       1:ny,          :, te) ! east
344     halo_data(1:nx,          1-shalo:0,     :) = data_all(1:nx,          ny-shalo+1:ny, :, ts) ! south
345     halo_data(1-whalo:0,     1:ny,          :) = data_all(nx-whalo+1:nx, 1:ny,          :, tw) ! west
346     halo_data(1:nx,          ny+1:ny+nhalo, :) = data_all(1:nx,          1:nhalo,       :, tn) ! north
347     halo_data(nx+1:nx+ehalo, 1-shalo:0,     :) = data_all(1:ehalo,       ny-shalo+1:ny, :,tse) ! southeast
348     halo_data(1-whalo:0,     1-shalo:0,     :) = data_all(nx-whalo+1:nx, ny-shalo+1:ny, :,tsw) ! southwest
349     halo_data(nx+1:nx+ehalo, ny+1:ny+nhalo, :) = data_all(1:ehalo,       1:nhalo,       :,tnw) ! northeast
350     halo_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(halo_data, data_all, te, tse, ts, tsw, tw, tnw, tn, tne)
355     integer(kind=i4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
356     integer(kind=i4_kind), dimension(:,:,:,:),             intent(in)    :: data_all
357     integer, intent(in)    :: te, tse, ts, tsw, tw, tnw, tn, tne
359     halo_data(nx+1:nx+ehalo, 1:ny,          :) = data_all(1:ehalo,       1:ny,          :, te) ! east
360     halo_data(1:nx,          1-shalo:0,     :) = data_all(1:nx,          ny-shalo+1:ny, :, ts) ! south
361     halo_data(1-whalo:0,     1:ny,          :) = data_all(nx-whalo+1:nx, 1:ny,          :, tw) ! west
362     halo_data(1:nx,          ny+1:ny+nhalo, :) = data_all(1:nx,          1:nhalo,       :, tn) ! north
363     halo_data(nx+1:nx+ehalo, 1-shalo:0,     :) = data_all(1:ehalo,       ny-shalo+1:ny, :,tse) ! southeast
364     halo_data(1-whalo:0,     1-shalo:0,     :) = data_all(nx-whalo+1:nx, ny-shalo+1:ny, :,tsw) ! southwest
365     halo_data(nx+1:nx+ehalo, ny+1:ny+nhalo, :) = data_all(1:ehalo,       1:nhalo,       :,tnw) ! northeast
366     halo_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(halo_data, ioff, joff, ishift, jshift, sign)
371     real(kind=r8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_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     halo_data(1-whalo:0,1:nyp,:) = halo_data(nx-whalo+1:nx,1:ny+jshift,:) ! west
382     halo_data(nx+1:nx+ehalo+ishift,1:nyp,:) = halo_data(1:ehalo+ishift,1:ny+jshift,:) ! east
383     if(m1 .GE. 1-whalo) &
384       halo_data(1-whalo:m1,nyp+1:nyp+nhalo,:) = sign*halo_data(whalo+m2:1+ishift:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
385     halo_data(m1+1:nx+m2,nyp+1:nyp+nhalo,:) = sign*halo_data(nx+ishift:1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
386     halo_data(nx+m2+1:nxp+ehalo,nyp+1:nyp+nhalo,:) = &
387         sign*halo_data(nx:nx-ehalo+m1+1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
389   end subroutine fill_folded_north_halo_r8
391   !> Fill the halo region of a 32-bit real array on a domain with a folded north edge
392   subroutine fill_folded_north_halo_r4(halo_data, ioff, joff, ishift, jshift, sign)
393     real(kind=r4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
394     integer, intent(in) :: ioff, joff, ishift, jshift, sign
395     ! local
396     integer :: nxp, nyp, m1, m2
398     nxp = nx+ishift
399     nyp = ny+jshift
400     m1 = ishift - ioff
401     m2 = 2*ishift - ioff
403     halo_data(1-whalo:0,1:nyp,:) = halo_data(nx-whalo+1:nx,1:ny+jshift,:) ! west
404     halo_data(nx+1:nx+ehalo+ishift,1:nyp,:) = halo_data(1:ehalo+ishift,1:ny+jshift,:) ! east
406     if(m1 .GE. 1-whalo) &
407       halo_data(1-whalo:m1,nyp+1:nyp+nhalo,:) = sign*halo_data(whalo+m2:1+ishift:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
408     halo_data(m1+1:nx+m2,nyp+1:nyp+nhalo,:) = sign*halo_data(nx+ishift:1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
409     halo_data(nx+m2+1:nxp+ehalo,nyp+1:nyp+nhalo,:) = &
410         sign*halo_data(nx:nx-ehalo+m1+1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
412   end subroutine fill_folded_north_halo_r4
414   !> Fill the halo region of a 64-bit integer array on a domain with a folded north edge
415   subroutine fill_folded_north_halo_i8(halo_data, ioff, joff, ishift, jshift, sign)
416     integer(kind=i8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
417     integer, intent(in) :: ioff, joff, ishift, jshift, sign
418     ! local
419     integer :: nxp, nyp, m1, m2
421     nxp = nx+ishift
422     nyp = ny+jshift
423     m1 = ishift - ioff
424     m2 = 2*ishift - ioff
426     halo_data(1-whalo:0,1:nyp,:) = halo_data(nx-whalo+1:nx,1:ny+jshift,:) ! west
427     halo_data(nx+1:nx+ehalo+ishift,1:nyp,:) = halo_data(1:ehalo+ishift,1:ny+jshift,:) ! east
428     if(m1 .GE. 1-whalo) &
429       halo_data(1-whalo:m1,nyp+1:nyp+nhalo,:) = sign*halo_data(whalo+m2:1+ishift:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
430     halo_data(m1+1:nx+m2,nyp+1:nyp+nhalo,:) = sign*halo_data(nx+ishift:1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
431     halo_data(nx+m2+1:nxp+ehalo,nyp+1:nyp+nhalo,:) = &
432         sign*halo_data(nx:nx-ehalo+m1+1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
434   end subroutine fill_folded_north_halo_i8
436   !> Fill the halo region of a 32-bit integer array on a domain with a folded north edge
437   subroutine fill_folded_north_halo_i4(halo_data, ioff, joff, ishift, jshift, sign)
438     integer(kind=i4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
439     integer, intent(in) :: ioff, joff, ishift, jshift, sign
440     ! local
441     integer :: nxp, nyp, m1, m2
443     nxp = nx+ishift
444     nyp = ny+jshift
445     m1 = ishift - ioff
446     m2 = 2*ishift - ioff
448     halo_data(1-whalo:0,1:nyp,:) = halo_data(nx-whalo+1:nx,1:ny+jshift,:) ! west
449     halo_data(nx+1:nx+ehalo+ishift,1:nyp,:) = halo_data(1:ehalo+ishift,1:ny+jshift,:) ! east
451     if(m1 .GE. 1-whalo) &
452       halo_data(1-whalo:m1,nyp+1:nyp+nhalo,:) = sign*halo_data(whalo+m2:1+ishift:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
453     halo_data(m1+1:nx+m2,nyp+1:nyp+nhalo,:) = sign*halo_data(nx+ishift:1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
454     halo_data(nx+m2+1:nxp+ehalo,nyp+1:nyp+nhalo,:) = &
455         sign*halo_data(nx:nx-ehalo+m1+1:-1,nyp-joff:nyp-nhalo-joff+1:-1,:)
457   end subroutine fill_folded_north_halo_i4
459   !> Fill the halo region of a 64-bit real array on a domain with a folded south edge
460   subroutine fill_folded_south_halo_r8(halo_data, ioff, joff, ishift, jshift, sign)
461     real(kind=r8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
462     integer, intent(in) :: ioff, joff, ishift, jshift, sign
463     ! local
464     integer  :: nxp, nyp, m1, m2
466     nxp = nx+ishift
467     nyp = ny+jshift
468     m1 = ishift - ioff
469     m2 = 2*ishift - ioff
471     halo_data(1-whalo:0,1:nyp,:) = halo_data(nx-whalo+1:nx,1:nyp,:) ! west
472     halo_data(nx+1:nx+ehalo+ishift,1:nyp,:) = halo_data(1:ehalo+ishift, 1:nyp,:) ! east
473     if(m1 .GE. 1-whalo) &
474       halo_data(1-whalo:m1,1-shalo:0,:) = sign*halo_data(whalo+m2:1+ishift:-1, shalo+jshift:1+jshift:-1,:)
476     halo_data(m1+1:nx+m2,1-shalo:0,:) = sign*halo_data(nxp:1:-1,shalo+jshift:1+jshift:-1,:)
477     halo_data(nx+m2+1:nxp+ehalo,1-shalo:0,:) = sign*halo_data(nx:nx-ehalo+m1+1:-1,shalo+jshift:1+jshift:-1,:)
479   end subroutine fill_folded_south_halo_r8
481   !> Fill the halo region of a 32-bit real array on a domain with a folded south edge
482   subroutine fill_folded_south_halo_r4(halo_data, ioff, joff, ishift, jshift, sign)
483     real(kind=r4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
484     integer, intent(in) :: ioff, joff, ishift, jshift, sign
485     ! local
486     integer  :: nxp, nyp, m1, m2
488     nxp = nx+ishift
489     nyp = ny+jshift
490     m1 = ishift - ioff
491     m2 = 2*ishift - ioff
493     halo_data(1-whalo:0,1:nyp,:) = halo_data(nx-whalo+1:nx,1:nyp,:) ! west
494     halo_data(nx+1:nx+ehalo+ishift,1:nyp,:) = halo_data(1:ehalo+ishift, 1:nyp,:) ! east
495     if(m1 .GE. 1-whalo) &
496       halo_data(1-whalo:m1,1-shalo:0,:) = sign*halo_data(whalo+m2:1+ishift:-1, shalo+jshift:1+jshift:-1,:)
498     halo_data(m1+1:nx+m2,1-shalo:0,:) = sign*halo_data(nxp:1:-1,shalo+jshift:1+jshift:-1,:)
499     halo_data(nx+m2+1:nxp+ehalo,1-shalo:0,:) = sign*halo_data(nx:nx-ehalo+m1+1:-1,shalo+jshift:1+jshift:-1,:)
501   end subroutine fill_folded_south_halo_r4
503   !> Fill the halo region of a 64-bit intger array on a domain with a folded south edge
504   subroutine fill_folded_south_halo_i8(halo_data, ioff, joff, ishift, jshift, sign)
505     integer(kind=i8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
506     integer, intent(in) :: ioff, joff, ishift, jshift, sign
507     ! local
508     integer  :: nxp, nyp, m1, m2
510     nxp = nx+ishift
511     nyp = ny+jshift
512     m1 = ishift - ioff
513     m2 = 2*ishift - ioff
515     halo_data(1-whalo:0,1:nyp,:) = halo_data(nx-whalo+1:nx,1:nyp,:) ! west
516     halo_data(nx+1:nx+ehalo+ishift,1:nyp,:) = halo_data(1:ehalo+ishift, 1:nyp,:) ! east
517     if(m1 .GE. 1-whalo) &
518       halo_data(1-whalo:m1,1-shalo:0,:) = sign*halo_data(whalo+m2:1+ishift:-1, shalo+jshift:1+jshift:-1,:)
520     halo_data(m1+1:nx+m2,1-shalo:0,:) = sign*halo_data(nxp:1:-1,shalo+jshift:1+jshift:-1,:)
521     halo_data(nx+m2+1:nxp+ehalo,1-shalo:0,:) = sign*halo_data(nx:nx-ehalo+m1+1:-1,shalo+jshift:1+jshift:-1,:)
523   end subroutine fill_folded_south_halo_i8
525   !> Fill the halo region of a 32-bit integer array on a domain with a folded south edge
526   subroutine fill_folded_south_halo_i4(halo_data, ioff, joff, ishift, jshift, sign)
527     integer(kind=i4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
528     integer, intent(in) :: ioff, joff, ishift, jshift, sign
529     ! local
530     integer  :: nxp, nyp, m1, m2
532     nxp = nx+ishift
533     nyp = ny+jshift
534     m1 = ishift - ioff
535     m2 = 2*ishift - ioff
537     halo_data(1-whalo:0,1:nyp,:) = halo_data(nx-whalo+1:nx,1:nyp,:) ! west
538     halo_data(nx+1:nx+ehalo+ishift,1:nyp,:) = halo_data(1:ehalo+ishift, 1:nyp,:) ! east
539     if(m1 .GE. 1-whalo) &
540       halo_data(1-whalo:m1,1-shalo:0,:) = sign*halo_data(whalo+m2:1+ishift:-1, shalo+jshift:1+jshift:-1,:)
542     halo_data(m1+1:nx+m2,1-shalo:0,:) = sign*halo_data(nxp:1:-1,shalo+jshift:1+jshift:-1,:)
543     halo_data(nx+m2+1:nxp+ehalo,1-shalo:0,:) = sign*halo_data(nx:nx-ehalo+m1+1:-1,shalo+jshift:1+jshift:-1,:)
545   end subroutine fill_folded_south_halo_i4
547   !> Fill the halo region of a 64-bit real array on a domain with a folded west edge
548   subroutine fill_folded_west_halo_r8(halo_data, ioff, joff, ishift, jshift, sign)
549     real(kind=r8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
550     integer, intent(in) :: ioff, joff, ishift, jshift, sign
551     ! local
552     integer :: nxp, nyp, m1, m2
554     nxp = nx+ishift
555     nyp = ny+jshift
556     m1 = jshift - joff
557     m2 = 2*jshift - joff
559     halo_data(1:nxp, 1-shalo:0,:) = halo_data(1:nxp, ny-shalo+1:ny, :) ! south
560     halo_data(1:nxp, ny+1:nyp+nhalo,:) = halo_data(1:nxp, 1:nhalo+jshift,:) ! north
561     if(m1 .GE. 1-shalo) &
562       halo_data(1-whalo:0, 1-shalo:m1,:) = sign*halo_data(whalo+ishift:1+ishift:-1,shalo+m2:1+jshift:-1,:)
563     halo_data(1-whalo:0, m1+1:ny+m2,:) = sign*halo_data(whalo+ishift:1+ishift:-1, nyp:1:-1, :)
564     halo_data(1-whalo:0, ny+m2+1:nyp+nhalo,:) = sign*halo_data(whalo+ishift:1+ishift:-1, ny:ny-nhalo+m1+1:-1,:)
566   end subroutine fill_folded_west_halo_r8
568   !> Fill the halo region of a 32-bit real array on a domain with a folded west edge
569   subroutine fill_folded_west_halo_r4(halo_data, ioff, joff, ishift, jshift, sign)
570     real(kind=r4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
571     integer, intent(in) :: ioff, joff, ishift, jshift, sign
572     ! local
573     integer :: nxp, nyp, m1, m2
575     nxp = nx+ishift
576     nyp = ny+jshift
577     m1 = jshift - joff
578     m2 = 2*jshift - joff
580     halo_data(1:nxp, 1-shalo:0,:) = halo_data(1:nxp, ny-shalo+1:ny, :) ! south
581     halo_data(1:nxp, ny+1:nyp+nhalo,:) = halo_data(1:nxp, 1:nhalo+jshift,:) ! north
582     if(m1 .GE. 1-shalo) &
583       halo_data(1-whalo:0, 1-shalo:m1,:) = sign*halo_data(whalo+ishift:1+ishift:-1,shalo+m2:1+jshift:-1,:)
584     halo_data(1-whalo:0, m1+1:ny+m2,:) = sign*halo_data(whalo+ishift:1+ishift:-1, nyp:1:-1, :)
585     halo_data(1-whalo:0, ny+m2+1:nyp+nhalo,:) = sign*halo_data(whalo+ishift:1+ishift:-1, ny:ny-nhalo+m1+1:-1,:)
587   end subroutine fill_folded_west_halo_r4
589   !> Fill the halo region of a 64-bit integer array on a domain with a folded west edge
590   subroutine fill_folded_west_halo_i8(halo_data, ioff, joff, ishift, jshift, sign)
591     integer(kind=i8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
592     integer, intent(in) :: ioff, joff, ishift, jshift, sign
593     ! local
594     integer :: nxp, nyp, m1, m2
596     nxp = nx+ishift
597     nyp = ny+jshift
598     m1 = jshift - joff
599     m2 = 2*jshift - joff
601     halo_data(1:nxp, 1-shalo:0,:) = halo_data(1:nxp, ny-shalo+1:ny, :) ! south
602     halo_data(1:nxp, ny+1:nyp+nhalo,:) = halo_data(1:nxp, 1:nhalo+jshift,:) ! north
603     if(m1 .GE. 1-shalo) &
604       halo_data(1-whalo:0, 1-shalo:m1,:) = sign*halo_data(whalo+ishift:1+ishift:-1,shalo+m2:1+jshift:-1,:)
605     halo_data(1-whalo:0, m1+1:ny+m2,:) = sign*halo_data(whalo+ishift:1+ishift:-1, nyp:1:-1, :)
606     halo_data(1-whalo:0, ny+m2+1:nyp+nhalo,:) = sign*halo_data(whalo+ishift:1+ishift:-1, ny:ny-nhalo+m1+1:-1,:)
608   end subroutine fill_folded_west_halo_i8
610   !> Fill the halo region of a 32-bit integer array on a domain with a folded west edge
611   subroutine fill_folded_west_halo_i4(halo_data, ioff, joff, ishift, jshift, sign)
612     integer(kind=i4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
613     integer, intent(in) :: ioff, joff, ishift, jshift, sign
614     ! local
615     integer :: nxp, nyp, m1, m2
617     nxp = nx+ishift
618     nyp = ny+jshift
619     m1 = jshift - joff
620     m2 = 2*jshift - joff
622     halo_data(1:nxp, 1-shalo:0,:) = halo_data(1:nxp, ny-shalo+1:ny, :) ! south
623     halo_data(1:nxp, ny+1:nyp+nhalo,:) = halo_data(1:nxp, 1:nhalo+jshift,:) ! north
624     if(m1 .GE. 1-shalo) &
625       halo_data(1-whalo:0, 1-shalo:m1,:) = sign*halo_data(whalo+ishift:1+ishift:-1,shalo+m2:1+jshift:-1,:)
626     halo_data(1-whalo:0, m1+1:ny+m2,:) = sign*halo_data(whalo+ishift:1+ishift:-1, nyp:1:-1, :)
627     halo_data(1-whalo:0, ny+m2+1:nyp+nhalo,:) = sign*halo_data(whalo+ishift:1+ishift:-1, ny:ny-nhalo+m1+1:-1,:)
629   end subroutine fill_folded_west_halo_i4
631   !> Fill the halo region of a 64-bit real array on a domain with a folded east edge
632   subroutine fill_folded_east_halo_r8(halo_data, ioff, joff, ishift, jshift, sign)
633     real(kind=r8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
634     integer, intent(in) :: ioff, joff, ishift, jshift, sign
635     ! local
636     integer :: nxp, nyp, m1, m2
638     nxp = nx+ishift
639     nyp = ny+jshift
640     m1 = jshift - joff
641     m2 = 2*jshift - joff
643     halo_data(1:nxp, 1-shalo:0, :)  = halo_data(1:nxp, ny-shalo+1:ny, :) ! south
644     halo_data(1:nxp, ny+1:nyp+nhalo, :) = halo_data(1:nxp, 1:nhalo+jshift,:) ! north
645     if(m1 .GE. 1-shalo) &
646       halo_data(nxp+1:nxp+ehalo, 1-shalo:m1, :) = sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, shalo+m2:1+jshift:-1,:)
648     halo_data(nxp+1:nxp+ehalo, m1+1:ny+m2, :) = sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, nyp:1:-1, :)
649     halo_data(nxp+1:nxp+ehalo, ny+m2+1:nyp+nhalo,:) = &
650       sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, ny:ny-nhalo+m1+1:-1,:)
652   end subroutine fill_folded_east_halo_r8
654   !> Fill the halo region of a 32-bit real array on a domain with a folded east edge
655   subroutine fill_folded_east_halo_r4(halo_data, ioff, joff, ishift, jshift, sign)
656     real(kind=r4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
657     integer, intent(in) :: ioff, joff, ishift, jshift, sign
658     ! local
659     integer :: nxp, nyp, m1, m2
661     nxp = nx+ishift
662     nyp = ny+jshift
663     m1 = jshift - joff
664     m2 = 2*jshift - joff
666     halo_data(1:nxp, 1-shalo:0, :)  = halo_data(1:nxp, ny-shalo+1:ny, :) ! south
667     halo_data(1:nxp, ny+1:nyp+nhalo, :) = halo_data(1:nxp, 1:nhalo+jshift,:) ! north
668     if(m1 .GE. 1-shalo) &
669       halo_data(nxp+1:nxp+ehalo, 1-shalo:m1, :) = sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, shalo+m2:1+jshift:-1,:)
671     halo_data(nxp+1:nxp+ehalo, m1+1:ny+m2, :) = sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, nyp:1:-1, :)
672     halo_data(nxp+1:nxp+ehalo, ny+m2+1:nyp+nhalo,:) = &
673       sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, ny:ny-nhalo+m1+1:-1,:)
675   end subroutine fill_folded_east_halo_r4
677   !> Fill the halo region of a 64-bit integer array on a domain with a folded east edge
678   subroutine fill_folded_east_halo_i8(halo_data, ioff, joff, ishift, jshift, sign)
679     integer(kind=i8_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
680     integer, intent(in) :: ioff, joff, ishift, jshift, sign
681     ! local
682     integer :: nxp, nyp, m1, m2
684     nxp = nx+ishift
685     nyp = ny+jshift
686     m1 = jshift - joff
687     m2 = 2*jshift - joff
689     halo_data(1:nxp, 1-shalo:0, :)  = halo_data(1:nxp, ny-shalo+1:ny, :) ! south
690     halo_data(1:nxp, ny+1:nyp+nhalo, :) = halo_data(1:nxp, 1:nhalo+jshift,:) ! north
691     if(m1 .GE. 1-shalo) &
692       halo_data(nxp+1:nxp+ehalo, 1-shalo:m1, :) = sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, shalo+m2:1+jshift:-1,:)
694     halo_data(nxp+1:nxp+ehalo, m1+1:ny+m2, :) = sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, nyp:1:-1, :)
695     halo_data(nxp+1:nxp+ehalo, ny+m2+1:nyp+nhalo,:) = &
696       sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, ny:ny-nhalo+m1+1:-1,:)
698   end subroutine fill_folded_east_halo_i8
700   !> Fill the halo region of a 32-bit integer array on a domain with a folded east edge
701   subroutine fill_folded_east_halo_i4(halo_data, ioff, joff, ishift, jshift, sign)
702     integer(kind=i4_kind), dimension(1-whalo:,1-shalo:,:), intent(inout) :: halo_data
703     integer, intent(in) :: ioff, joff, ishift, jshift, sign
704     ! local
705     integer :: nxp, nyp, m1, m2
707     nxp = nx+ishift
708     nyp = ny+jshift
709     m1 = jshift - joff
710     m2 = 2*jshift - joff
712     halo_data(1:nxp, 1-shalo:0, :)  = halo_data(1:nxp, ny-shalo+1:ny, :) ! south
713     halo_data(1:nxp, ny+1:nyp+nhalo, :) = halo_data(1:nxp, 1:nhalo+jshift,:) ! north
714     if(m1 .GE. 1-shalo) &
715       halo_data(nxp+1:nxp+ehalo, 1-shalo:m1, :) = sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, shalo+m2:1+jshift:-1,:)
717     halo_data(nxp+1:nxp+ehalo, m1+1:ny+m2, :) = sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, nyp:1:-1, :)
718     halo_data(nxp+1:nxp+ehalo, ny+m2+1:nyp+nhalo,:) = &
719       sign*halo_data(nxp-ioff:nxp-ehalo-ioff+1:-1, ny:ny-nhalo+m1+1:-1,:)
721   end subroutine fill_folded_east_halo_i4
723 end module fill_halo