Merge remote-tracking branch 'origin/release-v4.5'
[WRF.git] / frame / module_clear_halos.F
blob88de325ecf6de84a3d5eeebabd2be314214554e7
1 module module_clear_halos
2   implicit none
3 contains
4   ! --------------------------------------------------------------------
5   subroutine clear_ij_full_domain(grid,how)
6     ! Convenience function - wrapper around clear_ij_halos.  Clears
7     ! full domain with badval.  See clear_ij_halos for details.
8     use module_domain, only: domain,get_ijk_from_grid,fieldlist
9     type(domain), intent(inout) :: grid
10     integer, intent(in) :: how
11     !
12     call clear_ij_halos(grid,how,full_domain=.true.)
13   end subroutine clear_ij_full_domain
14   ! --------------------------------------------------------------------
15   subroutine clear_ij_halos(grid,how,full_domain)
16     ! Clears halo regions OR full domain with badval.  Select full
17     ! domain with full_domain=.true.  Select badval type with "how"
18     ! parameter:
20     ! how=1 -- badval=0
21     ! how=2 -- badval=quiet NaN or -maxint
22     ! how=3 -- badval=signaling NaN or -maxint
24     ! Fills outside domain with 0 UNLESS fill_domain=.true.  If
25     ! fill_domain=true., entire array is filled with badval.
27     use module_domain, only: domain,get_ijk_from_grid,fieldlist
28     use module_configure, only: PARAM_FIRST_SCALAR
29 #ifndef NO_IEEE_MODULE
30     use,intrinsic :: ieee_arithmetic
31 #endif
32     implicit none
34     logical, intent(in), optional :: full_domain
35     integer, intent(in) :: how
36     type(domain), intent(inout) :: grid
38     type( fieldlist ), pointer :: p
39     integer :: itrace, i,j, &
40          ids, ide, jds, jde, kds, kde,    &
41          ims, ime, jms, jme, kms, kme,    &
42          ips, ipe, jps, jpe, kps, kpe
43     logical :: fulldom
44     real :: badR, badR_N,badR_NE,badR_NW,badR_S,badR_SW,badR_SE,badR_E,badR_W
45 #if (RWORDSIZE==4)
46     double precision :: badD, badD_N,badD_NE,badD_NW,badD_S,badD_SW,badD_SE,badD_E,badD_W
47 #else
48     real             :: badD, badD_N,badD_NE,badD_NW,badD_S,badD_SW,badD_SE,badD_E,badD_W
49 #endif
50     integer :: badI, badI_N,badI_NE,badI_NW,badI_S,badI_SW,badI_SE,badI_E,badI_W
52     select case(how)
53     case(0)
54        return
55     case(1)
56        call wrf_message('Fill I and J halos with 0.')
57        badR = 0
58        badD = 0
59        badI = 0
60     case(2)
61        call wrf_message('Fill I and J halos with -maxint or quiet NaN.')
62 #ifndef NO_IEEE_MODULE
63        badR = ieee_value(badR,ieee_quiet_nan)
64        badD = ieee_value(badD,ieee_quiet_nan)
65        badI = -huge(badI)
66 #else
67        badR = -huge(badR) 
68        badD = -huge(badD)
69        badI = -huge(badI)
70 #endif
71     case(3)
72        call wrf_message('Fill I and J halos with -maxint or signalling NaN.')
73 #ifndef NO_IEEE_MODULE
74        badR = ieee_value(badR,ieee_signaling_nan)
75        badD = ieee_value(badD,ieee_signaling_nan)
76        badI = -huge(badI)
77 #else
78        badR = -huge(badR) 
79        badD = -huge(badD)
80        badI = -huge(badI)
81 #endif
82     case default
83        if(fulldom) then
84           call wrf_message('Invalid value for clear_ij_full_domain/clear_ij_halos "how" parameter.  Will not clear domain.')
85        else
86           call wrf_message('Invalid value for clear_ij_halos "how" parameter.  Will not clear halos.')
87        endif
88        return
89     end select
91     fulldom=.false.
92     if(present(full_domain)) fulldom=full_domain
93     if(fulldom) then
94        call wrf_message('Filling entire memory area, not just halos.')
95     endif
97     badR_N =badR ; badD_N =badD ; badI_N =badI
98     badR_NE=badR ; badD_NE=badD ; badI_NE=badI
99     badR_NW=badR ; badD_NW=badD ; badI_NW=badI
100     badR_S =badR ; badD_S =badD ; badI_S =badI
101     badR_SE=badR ; badD_SE=badD ; badI_SE=badI
102     badR_SW=badR ; badD_SW=badD ; badI_SW=badI
103     badR_E =badR ; badD_E =badD ; badI_E =badI
104     badR_W =badR ; badD_W =badD ; badI_W =badI
106     CALL get_ijk_from_grid (  grid ,      &
107          ids, ide, jds, jde, kds, kde,    &
108          ims, ime, jms, jme, kms, kme,    &
109          ips, ipe, jps, jpe, kps, kpe     )
111     if(ips==ids) then
112        badR_S =0 ; badD_S =0 ; badI_S =0
113        badR_SE=0 ; badD_SE=0 ; badI_SE=0
114        badR_SW=0 ; badD_SW=0 ; badI_SW=0
115     endif
116     if(ipe==ide) then
117        badR_N =0 ; badD_N =0 ; badI_N =0
118        badR_NE=0 ; badD_NE=0 ; badI_NE=0
119        badR_NW=0 ; badD_NW=0 ; badI_NW=0
120     endif
121     if(jps==jds) then
122        badR_NW=0 ; badD_NW=0 ; badI_NW=0
123        badR_SW=0 ; badD_SW=0 ; badI_SW=0
124        badR_W =0 ; badD_W =0 ; badI_W =0
125     endif
126     if(jpe==jde) then
127        badR_NE=0 ; badD_NE=0 ; badI_NE=0
128        badR_SE=0 ; badD_SE=0 ; badI_SE=0
129        badR_E =0 ; badD_E =0 ; badI_E =0
130     endif
132     if(.not.associated(grid%head_statevars)) then
133        call wrf_message('grid%head_statevars is not associated')
134        return
135     elseif(.not.associated(grid%head_statevars%next)) then
136        call wrf_message('grid%head_statevars%next is not associated')
137        return
138     endif
139     p => grid%head_statevars%next
140     DO WHILE ( ASSOCIATED( p ) ) 
141        IF ( p%ProcOrient .NE. 'X' .AND. p%ProcOrient .NE. 'Y' ) THEN
142           IF ( p%Ndim .EQ. 2 ) THEN
143              IF (      p%MemoryOrder(1:1) .EQ. 'X' .AND.  p%MemoryOrder(2:2) .EQ.  'Y' ) THEN
144                 IF      ( p%Type .EQ. 'r' ) THEN
145                    IF ( SIZE(p%rfield_2d,1)*SIZE(p%rfield_2d,2) .GT. 1 ) THEN
146                       if(fulldom) then
147                          p%rfield_2d=badR
148                       else
149                          p%rfield_2d(ims:ips-1,jps:jpe) = badR_S
150                          p%rfield_2d(ims:ips-1,jms:jps-1) = badR_SW
151                          p%rfield_2d(ims:ips-1,jpe+1:jme) = badR_SE
152                          p%rfield_2d(ipe+1:ime,jps:jpe) = badR_N
153                          p%rfield_2d(ipe+1:ime,jms:jps-1) = badR_NW
154                          p%rfield_2d(ipe+1:ime,jpe+1:jme) = badR_NE
155                          p%rfield_2d(ips:ipe,jms:jps-1) = badR_W
156                          p%rfield_2d(ips:ipe,jpe+1:jme) = badR_E
157                       endif
158                    ENDIF
159                 ELSE IF ( p%Type .EQ. 'd' ) THEN
160                    IF ( SIZE(p%dfield_2d,1)*SIZE(p%dfield_2d,2) .GT. 1 ) THEN
161                       if(fulldom) then
162                          p%dfield_2d=badD
163                       else
164                          p%dfield_2d(ims:ips-1,jps:jpe) = badD_S
165                          p%dfield_2d(ims:ips-1,jms:jps-1) = badD_SW
166                          p%dfield_2d(ims:ips-1,jpe+1:jme) = badD_SE
167                          p%dfield_2d(ipe+1:ime,jps:jpe) = badD_N
168                          p%dfield_2d(ipe+1:ime,jms:jps-1) = badD_NW
169                          p%dfield_2d(ipe+1:ime,jpe+1:jme) = badD_NE
170                          p%dfield_2d(ips:ipe,jms:jps-1) = badD_W
171                          p%dfield_2d(ips:ipe,jpe+1:jme) = badD_E
172                       endif
173                    ENDIF
174                 ELSE IF ( p%Type .EQ. 'i' ) THEN
175                    IF ( SIZE(p%ifield_2d,1)*SIZE(p%ifield_2d,2) .GT. 1 ) THEN
176                       if(fulldom) then
177                          p%ifield_2d=badI
178                       else
179                          p%ifield_2d(ims:ips-1,jps:jpe) = badI_S
180                          p%ifield_2d(ims:ips-1,jms:jps-1) = badI_SW
181                          p%ifield_2d(ims:ips-1,jpe+1:jme) = badI_SE
182                          p%ifield_2d(ipe+1:ime,jps:jpe) = badI_N
183                          p%ifield_2d(ipe+1:ime,jms:jps-1) = badI_NW
184                          p%ifield_2d(ipe+1:ime,jpe+1:jme) = badI_NE
185                          p%ifield_2d(ips:ipe,jms:jps-1) = badI_W
186                          p%ifield_2d(ips:ipe,jpe+1:jme) = badI_E
187                       endif
188                    ENDIF
189                 ENDIF
190              ENDIF
191           ELSE IF ( p%Ndim .EQ. 3 ) THEN
192              IF (      p%MemoryOrder(1:1) .EQ. 'X' .AND.  p%MemoryOrder(3:3) .EQ.  'Y' ) THEN
193                 IF      ( p%Type .EQ. 'r' ) THEN
194                    IF ( SIZE(p%rfield_3d,1)*SIZE(p%rfield_3d,3) .GT. 1 ) THEN
195                       if(fulldom) then
196                          p%rfield_3d=badR
197                       else
198                          p%rfield_3d(ims:ips-1,:,jps:jpe) = badR_S
199                          p%rfield_3d(ims:ips-1,:,jms:jps-1) = badR_SW
200                          p%rfield_3d(ims:ips-1,:,jpe+1:jme) = badR_SE
201                          p%rfield_3d(ipe+1:ime,:,jps:jpe) = badR_N
202                          p%rfield_3d(ipe+1:ime,:,jms:jps-1) = badR_NW
203                          p%rfield_3d(ipe+1:ime,:,jpe+1:jme) = badR_NE
204                          p%rfield_3d(ips:ipe,:,jms:jps-1) = badR_W
205                          p%rfield_3d(ips:ipe,:,jpe+1:jme) = badR_E
206                       endif
207                    ENDIF
208                 ELSE IF ( p%Type .EQ. 'd' ) THEN
209                    IF ( SIZE(p%dfield_3d,1)*SIZE(p%dfield_3d,3) .GT. 1 ) THEN
210                       if(fulldom) then
211                          p%dfield_3d=badD
212                       else
213                          p%dfield_3d(ims:ips-1,:,jps:jpe) = badD_S
214                          p%dfield_3d(ims:ips-1,:,jms:jps-1) = badD_SW
215                          p%dfield_3d(ims:ips-1,:,jpe+1:jme) = badD_SE
216                          p%dfield_3d(ipe+1:ime,:,jps:jpe) = badD_N
217                          p%dfield_3d(ipe+1:ime,:,jms:jps-1) = badD_NW
218                          p%dfield_3d(ipe+1:ime,:,jpe+1:jme) = badD_NE
219                          p%dfield_3d(ips:ipe,:,jms:jps-1) = badD_W
220                          p%dfield_3d(ips:ipe,:,jpe+1:jme) = badD_E
221                       endif
222                    ENDIF
223                 ELSE IF ( p%Type .EQ. 'i' ) THEN
224                    IF ( SIZE(p%ifield_3d,1)*SIZE(p%ifield_3d,3) .GT. 1 ) THEN
225                       if(fulldom) then
226                          p%ifield_3d=badI
227                       else
228                          p%ifield_3d(ims:ips-1,:,jps:jpe) = badI_S
229                          p%ifield_3d(ims:ips-1,:,jms:jps-1) = badI_SW
230                          p%ifield_3d(ims:ips-1,:,jpe+1:jme) = badI_SE
231                          p%ifield_3d(ipe+1:ime,:,jps:jpe) = badI_N
232                          p%ifield_3d(ipe+1:ime,:,jms:jps-1) = badI_NW
233                          p%ifield_3d(ipe+1:ime,:,jpe+1:jme) = badI_NE
234                          p%ifield_3d(ips:ipe,:,jms:jps-1) = badI_W
235                          p%ifield_3d(ips:ipe,:,jpe+1:jme) = badI_E
236                       endif
237                    ENDIF
238                 ENDIF
239              ELSE IF (  p%MemoryOrder(1:2) .EQ. 'XY' ) THEN
240                 IF      ( p%Type .EQ. 'r' ) THEN
241                    IF ( SIZE(p%rfield_3d,1)*SIZE(p%rfield_3d,2) .GT. 1 ) THEN
242                       if(fulldom) then
243                          p%rfield_3d=badR
244                       else
245                          p%rfield_3d(ims:ips-1,jps:jpe,:) = badR_S
246                          p%rfield_3d(ims:ips-1,jms:jps-1,:) = badR_SW
247                          p%rfield_3d(ims:ips-1,jpe+1:jme,:) = badR_SE
248                          p%rfield_3d(ipe+1:ime,jps:jpe,:) = badR_N
249                          p%rfield_3d(ipe+1:ime,jms:jps-1,:) = badR_NW
250                          p%rfield_3d(ipe+1:ime,jpe+1:jme,:) = badR_NE
251                          p%rfield_3d(ips:ipe,jms:jps-1,:) = badR_W
252                          p%rfield_3d(ips:ipe,jpe+1:jme,:) = badR_E
253                       endif
254                    ENDIF
255                 ELSE IF ( p%Type .EQ. 'd' ) THEN
256                    IF ( SIZE(p%dfield_3d,1)*SIZE(p%dfield_3d,2) .GT. 1 ) THEN
257                       if(fulldom) then
258                          p%dfield_3d=badD
259                       else
260                          p%dfield_3d(ims:ips-1,jps:jpe,:) = badD_S
261                          p%dfield_3d(ims:ips-1,jms:jps-1,:) = badD_SW
262                          p%dfield_3d(ims:ips-1,jpe+1:jme,:) = badD_SE
263                          p%dfield_3d(ipe+1:ime,jps:jpe,:) = badD_N
264                          p%dfield_3d(ipe+1:ime,jms:jps-1,:) = badD_NW
265                          p%dfield_3d(ipe+1:ime,jpe+1:jme,:) = badD_NE
266                          p%dfield_3d(ips:ipe,jms:jps-1,:) = badD_W
267                          p%dfield_3d(ips:ipe,jpe+1:jme,:) = badD_E
268                       endif
269                    ENDIF
270                 ELSE IF ( p%Type .EQ. 'i' ) THEN
271                    IF ( SIZE(p%ifield_3d,1)*SIZE(p%ifield_3d,2) .GT. 1 ) THEN
272                       if(fulldom) then
273                          p%ifield_3d=badI
274                       else
275                          p%ifield_3d(ims:ips-1,jps:jpe,:) = badI_S
276                          p%ifield_3d(ims:ips-1,jms:jps-1,:) = badI_SW
277                          p%ifield_3d(ims:ips-1,jpe+1:jme,:) = badI_SE
278                          p%ifield_3d(ipe+1:ime,jps:jpe,:) = badI_N
279                          p%ifield_3d(ipe+1:ime,jms:jps-1,:) = badI_NW
280                          p%ifield_3d(ipe+1:ime,jpe+1:jme,:) = badI_NE
281                          p%ifield_3d(ips:ipe,jms:jps-1,:) = badI_W
282                          p%ifield_3d(ips:ipe,jpe+1:jme,:) = badI_E
283                       endif
284                    ENDIF
285                 ENDIF
286              ENDIF
287           ELSE IF ( p%Ndim .EQ. 4 ) THEN
288              IF (      p%MemoryOrder(1:1) .EQ. 'X' .AND.  p%MemoryOrder(3:3) .EQ.  'Y' ) THEN
289                 IF      ( p%Type .EQ. 'r' ) THEN
290                    IF ( SIZE(p%rfield_4d,1)*SIZE(p%rfield_4d,3) .GT. 1 ) THEN
291                       DO itrace = PARAM_FIRST_SCALAR , p%num_table(grid%id)
292                          if(fulldom) then
293                             p%rfield_4d(:,:,:,itrace)=badR
294                          else
295                             p%rfield_4d(ims:ips-1,:,jps:jpe,itrace) = badR_S
296                             p%rfield_4d(ims:ips-1,:,jms:jps-1,itrace) = badR_SW
297                             p%rfield_4d(ims:ips-1,:,jpe+1:jme,itrace) = badR_SE
298                             p%rfield_4d(ipe+1:ime,:,jps:jpe,itrace) = badR_N
299                             p%rfield_4d(ipe+1:ime,:,jms:jps-1,itrace) = badR_NW
300                             p%rfield_4d(ipe+1:ime,:,jpe+1:jme,itrace) = badR_NE
301                             p%rfield_4d(ips:ipe,:,jms:jps-1,itrace) = badR_W
302                             p%rfield_4d(ips:ipe,:,jpe+1:jme,itrace) = badR_E
303                          endif
304                       ENDDO
305                    ENDIF
306                 ELSE IF ( p%Type .EQ. 'd' ) THEN
307                    IF ( SIZE(p%dfield_4d,1)*SIZE(p%dfield_4d,3) .GT. 1 ) THEN
308                       DO itrace = PARAM_FIRST_SCALAR , p%num_table(grid%id)
309                          if(fulldom) then
310                             p%dfield_4d(:,:,:,itrace)=badD
311                          else
312                             p%dfield_4d(ims:ips-1,:,jps:jpe,itrace) = badD_S
313                             p%dfield_4d(ims:ips-1,:,jms:jps-1,itrace) = badD_SW
314                             p%dfield_4d(ims:ips-1,:,jpe+1:jme,itrace) = badD_SE
315                             p%dfield_4d(ipe+1:ime,:,jps:jpe,itrace) = badD_N
316                             p%dfield_4d(ipe+1:ime,:,jms:jps-1,itrace) = badD_NW
317                             p%dfield_4d(ipe+1:ime,:,jpe+1:jme,itrace) = badD_NE
318                             p%dfield_4d(ips:ipe,:,jms:jps-1,itrace) = badD_W
319                             p%dfield_4d(ips:ipe,:,jpe+1:jme,itrace) = badD_E
320                          endif
321                       ENDDO
322                    ENDIF
323                 ELSE IF ( p%Type .EQ. 'i' ) THEN
324                    IF ( SIZE(p%ifield_4d,1)*SIZE(p%ifield_4d,3) .GT. 1 ) THEN
325                       DO itrace = PARAM_FIRST_SCALAR , p%num_table(grid%id)
326                          if(fulldom) then
327                             p%ifield_4d(:,:,:,itrace)=badI
328                          else
329                             p%ifield_4d(ims:ips-1,:,jps:jpe,itrace) = badI_S
330                             p%ifield_4d(ims:ips-1,:,jms:jps-1,itrace) = badI_SW
331                             p%ifield_4d(ims:ips-1,:,jpe+1:jme,itrace) = badI_SE
332                             p%ifield_4d(ipe+1:ime,:,jps:jpe,itrace) = badI_N
333                             p%ifield_4d(ipe+1:ime,:,jms:jps-1,itrace) = badI_NW
334                             p%ifield_4d(ipe+1:ime,:,jpe+1:jme,itrace) = badI_NE
335                             p%ifield_4d(ips:ipe,:,jms:jps-1,itrace) = badI_W
336                             p%ifield_4d(ips:ipe,:,jpe+1:jme,itrace) = badI_E
337                          endif
338                       ENDDO
339                    ENDIF
340                 ENDIF
341              ELSE IF (  p%MemoryOrder(1:2) .EQ. 'XY' ) THEN
342                 IF      ( p%Type .EQ. 'r' ) THEN
343                    IF ( SIZE(p%rfield_4d,1)*SIZE(p%rfield_4d,2) .GT. 1 ) THEN
344                       DO itrace = PARAM_FIRST_SCALAR , p%num_table(grid%id)
345                          if(fulldom) then
346                             p%rfield_4d(:,:,:,itrace)=badR
347                          else
348                             p%rfield_4d(ims:ips-1,jps:jpe,:,itrace) = badR_S
349                             p%rfield_4d(ims:ips-1,jms:jps-1,:,itrace) = badR_SW
350                             p%rfield_4d(ims:ips-1,jpe+1:jme,:,itrace) = badR_SE
351                             p%rfield_4d(ipe+1:ime,jps:jpe,:,itrace) = badR_N
352                             p%rfield_4d(ipe+1:ime,jms:jps-1,:,itrace) = badR_NW
353                             p%rfield_4d(ipe+1:ime,jpe+1:jme,:,itrace) = badR_NE
354                             p%rfield_4d(ips:ipe,jms:jps-1,:,itrace) = badR_W
355                             p%rfield_4d(ips:ipe,jpe+1:jme,:,itrace) = badR_E
356                          endif
357                       ENDDO
358                    ENDIF
359                 ELSE IF ( p%Type .EQ. 'd' ) THEN
360                    IF ( SIZE(p%dfield_4d,1)*SIZE(p%dfield_4d,2) .GT. 1 ) THEN
361                       DO itrace = PARAM_FIRST_SCALAR , p%num_table(grid%id)
362                          if(fulldom) then
363                             p%dfield_4d(:,:,:,itrace)=badD
364                          else
365                             p%dfield_4d(ims:ips-1,jps:jpe,:,itrace) = badD_S
366                             p%dfield_4d(ims:ips-1,jms:jps-1,:,itrace) = badD_SW
367                             p%dfield_4d(ims:ips-1,jpe+1:jme,:,itrace) = badD_SE
368                             p%dfield_4d(ipe+1:ime,jps:jpe,:,itrace) = badD_N
369                             p%dfield_4d(ipe+1:ime,jms:jps-1,:,itrace) = badD_NW
370                             p%dfield_4d(ipe+1:ime,jpe+1:jme,:,itrace) = badD_NE
371                             p%dfield_4d(ips:ipe,jms:jps-1,:,itrace) = badD_W
372                             p%dfield_4d(ips:ipe,jpe+1:jme,:,itrace) = badD_E
373                          endif
374                       ENDDO
375                    ENDIF
376                 ELSE IF ( p%Type .EQ. 'i' ) THEN
377                    IF ( SIZE(p%ifield_4d,1)*SIZE(p%ifield_4d,2) .GT. 1 ) THEN
378                       DO itrace = PARAM_FIRST_SCALAR , p%num_table(grid%id)
379                          if(fulldom) then
380                             p%ifield_4d(:,:,:,itrace)=badI
381                          else
382                             p%ifield_4d(ims:ips-1,jps:jpe,:,itrace) = badI_S
383                             p%ifield_4d(ims:ips-1,jms:jps-1,:,itrace) = badI_SW
384                             p%ifield_4d(ims:ips-1,jpe+1:jme,:,itrace) = badI_SE
385                             p%ifield_4d(ipe+1:ime,jps:jpe,:,itrace) = badI_N
386                             p%ifield_4d(ipe+1:ime,jms:jps-1,:,itrace) = badI_NW
387                             p%ifield_4d(ipe+1:ime,jpe+1:jme,:,itrace) = badI_NE
388                             p%ifield_4d(ips:ipe,jms:jps-1,:,itrace) = badI_W
389                             p%ifield_4d(ips:ipe,jpe+1:jme,:,itrace) = badI_E
390                          endif
391                       ENDDO
392                    ENDIF
393                 ENDIF
394              ENDIF
395           ENDIF
396        ENDIF
397        p => p%next
398     ENDDO
399   end subroutine clear_ij_halos
400 end module module_clear_halos