New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
modarrays.F90 in branches/UKMO/dev_r5518_GO6_package/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES – NEMO

source: branches/UKMO/dev_r5518_GO6_package/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modarrays.F90 @ 7993

Last change on this file since 7993 was 7993, checked in by frrh, 7 years ago

Merge in missing revisions 6428:2477 inclusive and 6482 from nemo_v3_6_STABLE
branch. In ptic, this includes the fix for restartability of runoff fields in coupled
models. Evolution of coupled models will therefor be affected.

These changes donot affect evolution of the current stand-alone NEMO-CICE GO6
standard configuration.

Work and testing documented in Met Office GMED ticket 320.

File size: 38.9 KB
Line 
1!
2! $Id$
3!
4!     AGRIF (Adaptive Grid Refinement In Fortran)
5!
6!     Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
7!                Christophe Vouland (Christophe.Vouland@imag.fr)
8!
9!     This program is free software; you can redistribute it and/or modify
10!     it under the terms of the GNU General Public License as published by
11!     the Free Software Foundation; either version 2 of the License, or
12!     (at your option) any later version.
13!
14!     This program is distributed in the hope that it will be useful,
15!     but WITHOUT ANY WARRANTY; without even the implied warranty of
16!     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
17!     GNU General Public License for more details.
18!
19!     You should have received a copy of the GNU General Public License
20!     along with this program; if not, write to the Free Software
21!     Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA.
22!
23!
24!> Module Agrif_Arrays
25!
26module Agrif_Arrays
27!
28    use Agrif_Types
29    use Agrif_Grids
30!
31    implicit none
32!
33#if defined AGRIF_MPI
34    interface
35        subroutine Agrif_InvLoc ( indloc, proc_id, dir, indglob )
36            integer, intent(in)     :: indloc   !< local index
37            integer, intent(in)     :: proc_id  !< rank of the proc calling this function
38            integer, intent(in)     :: dir      !< direction of the index
39            integer, intent(out)    :: indglob  !< global index
40        end subroutine Agrif_InvLoc
41    end interface
42    private :: Agrif_InvLoc
43#endif
44!
45contains
46!
47!===================================================================================================
48! subroutine Agrif_Childbounds
49!
50!> Computes the global indices of the child grid
51!---------------------------------------------------------------------------------------------------
52subroutine Agrif_Childbounds ( nbdim,           &
53                               lb_var, ub_var,  &
54                               lb_tab, ub_tab,  &
55                               proc_id,         &
56                               coords,          &
57                               lb_tab_true, ub_tab_true, memberin )
58!---------------------------------------------------------------------------------------------------
59    integer,                   intent(in)  :: nbdim         !< Number of dimensions
60    integer, dimension(nbdim), intent(in)  :: lb_var        !< Local lower boundary on the current processor
61    integer, dimension(nbdim), intent(in)  :: ub_var        !< Local upper boundary on the current processor
62    integer, dimension(nbdim), intent(in)  :: lb_tab        !< Global lower boundary of the variable
63    integer, dimension(nbdim), intent(in)  :: ub_tab        !< Global upper boundary of the variable
64    integer,                   intent(in)  :: proc_id       !< Current processor
65    integer, dimension(nbdim), intent(in)  :: coords
66    integer, dimension(nbdim), intent(out) :: lb_tab_true   !< Global value of lb_var on the current processor
67    integer, dimension(nbdim), intent(out) :: ub_tab_true   !< Global value of ub_var on the current processor
68    logical,                   intent(out) :: memberin
69!
70    integer :: i, coord_i
71    integer :: lb_glob_index, ub_glob_index ! Lower and upper global indices
72!
73    do i = 1, nbdim
74!
75        coord_i = coords(i)
76!
77#if defined AGRIF_MPI
78        call Agrif_InvLoc( lb_var(i), proc_id, coord_i, lb_glob_index )
79        call Agrif_InvLoc( ub_var(i), proc_id, coord_i, ub_glob_index )
80#else
81        lb_glob_index = lb_var(i)
82        ub_glob_index = ub_var(i)
83#endif
84        lb_tab_true(i) = max(lb_tab(i), lb_glob_index)
85        ub_tab_true(i) = min(ub_tab(i), ub_glob_index)
86    enddo
87!
88    memberin = .true.
89    do i = 1,nbdim
90        if (ub_tab_true(i) < lb_tab_true(i)) then
91            memberin = .false.
92            exit
93        endif
94    enddo
95!---------------------------------------------------------------------------------------------------
96end subroutine Agrif_Childbounds
97!===================================================================================================
98!
99!===================================================================================================
100subroutine Agrif_get_var_global_bounds( var, lubglob, nbdim )
101!---------------------------------------------------------------------------------------------------
102    type(Agrif_Variable),        intent(in)  :: var
103    integer, dimension(nbdim,2), intent(out) :: lubglob
104    integer,                     intent(in)  :: nbdim
105!
106#if defined AGRIF_MPI
107    include 'mpif.h'
108    integer, dimension(nbdim)   :: lb, ub
109    integer, dimension(nbdim,2) :: iminmaxg
110    integer :: i, code, coord_i
111#endif
112!
113#if !defined AGRIF_MPI
114    call Agrif_get_var_bounds_array(var, lubglob(:,1), lubglob(:,2), nbdim)
115#else
116    call Agrif_get_var_bounds_array(var, lb, ub, nbdim)
117
118    do i = 1,nbdim
119        coord_i = var % root_var % coords(i)
120        call Agrif_InvLoc( lb(i), Agrif_Procrank, coord_i, iminmaxg(i,1) )
121        call Agrif_InvLoc( ub(i), Agrif_Procrank, coord_i, iminmaxg(i,2) )
122    enddo
123!
124    iminmaxg(1:nbdim,2) = - iminmaxg(1:nbdim,2)
125    call MPI_ALLREDUCE(iminmaxg, lubglob, 2*nbdim, MPI_INTEGER, MPI_MIN, Agrif_mpi_comm, code)
126    lubglob(1:nbdim,2)  = - lubglob(1:nbdim,2)
127#endif
128!---------------------------------------------------------------------------------------------------
129end subroutine Agrif_get_var_global_bounds
130!===================================================================================================
131!
132!===================================================================================================
133!  subroutine Agrif_get_var_bounds
134!
135!> Gets the lower and the upper boundaries of a variable, for one particular direction.
136!---------------------------------------------------------------------------------------------------
137subroutine Agrif_get_var_bounds ( variable, lower, upper, index )
138!---------------------------------------------------------------------------------------------------
139    type(Agrif_Variable), intent(in)    :: variable   !< Variable for which we want to extract boundaries
140    integer,              intent(out)   :: lower      !< Lower bound
141    integer,              intent(out)   :: upper      !< Upper bound
142    integer,              intent(in)    :: index      !< Direction for wich we want to know the boundaries
143!
144    lower = variable % lb(index)
145    upper = variable % ub(index)
146!---------------------------------------------------------------------------------------------------
147end subroutine Agrif_get_var_bounds
148!===================================================================================================
149!
150!===================================================================================================
151!  subroutine Agrif_get_var_bounds_array
152!
153!> Gets the lower and the upper boundaries of a table.
154!---------------------------------------------------------------------------------------------------
155subroutine Agrif_get_var_bounds_array ( variable, lower, upper, nbdim )
156!---------------------------------------------------------------------------------------------------
157    type(Agrif_Variable),      intent(in)   :: variable   !< Variable for which we want to extract boundaries
158    integer, dimension(nbdim), intent(out)  :: lower      !< Lower bounds array
159    integer, dimension(nbdim), intent(out)  :: upper      !< Upper bounds array
160    integer, intent(in)                     :: nbdim      !< Numer of dimensions of the variable
161!
162    lower = variable % lb(1:nbdim)
163    upper = variable % ub(1:nbdim)
164!---------------------------------------------------------------------------------------------------
165end subroutine Agrif_get_var_bounds_array
166!===================================================================================================
167!
168!===================================================================================================
169!  subroutine Agrif_array_allocate
170!
171!> Allocates data array in \b variable, according to the required dimension.
172!---------------------------------------------------------------------------------------------------
173subroutine Agrif_array_allocate ( variable, lb, ub, nbdim )
174!---------------------------------------------------------------------------------------------------
175    type(Agrif_Variable),      intent(inout) :: variable    !< Variable struct for allocation
176    integer, dimension(nbdim), intent(in)    :: lb          !< Lower bound
177    integer, dimension(nbdim), intent(in)    :: ub          !< Upper bound
178    integer,                   intent(in)    :: nbdim       !< Dimension of the array
179!
180    select case (nbdim)
181    case (1) ; allocate(variable%array1(lb(1):ub(1)))
182    case (2) ; allocate(variable%array2(lb(1):ub(1),lb(2):ub(2)))
183    case (3) ; allocate(variable%array3(lb(1):ub(1),lb(2):ub(2),lb(3):ub(3)))
184    case (4) ; allocate(variable%array4(lb(1):ub(1),lb(2):ub(2),lb(3):ub(3),lb(4):ub(4)))
185    case (5) ; allocate(variable%array5(lb(1):ub(1),lb(2):ub(2),lb(3):ub(3),lb(4):ub(4),lb(5):ub(5)))
186    case (6) ; allocate(variable%array6(lb(1):ub(1),lb(2):ub(2),lb(3):ub(3),lb(4):ub(4),lb(5):ub(5),lb(6):ub(6)))
187    end select
188!---------------------------------------------------------------------------------------------------
189end subroutine Agrif_array_allocate
190!===================================================================================================
191!
192!===================================================================================================
193!  subroutine Agrif_array_deallocate
194!
195!> Dellocates data array in \b variable, according to the required dimension.
196!---------------------------------------------------------------------------------------------------
197subroutine Agrif_array_deallocate ( variable, nbdim )
198!---------------------------------------------------------------------------------------------------
199    type(Agrif_Variable), intent(inout) :: variable     !< Variable struct for deallocation
200    integer,              intent(in)    :: nbdim        !< Dimension of the array
201!
202    select case (nbdim)
203    case (1) ; deallocate(variable%array1)
204    case (2) ; deallocate(variable%array2)
205    case (3) ; deallocate(variable%array3)
206    case (4) ; deallocate(variable%array4)
207    case (5) ; deallocate(variable%array5)
208    case (6) ; deallocate(variable%array6)
209    end select
210!---------------------------------------------------------------------------------------------------
211end subroutine Agrif_array_deallocate
212!===================================================================================================
213!
214!===================================================================================================
215!  subroutine Agrif_var_set_array_tozero
216!
217!> Reset the value of the data array in \b variable, according to the required dimension.
218!---------------------------------------------------------------------------------------------------
219subroutine Agrif_var_set_array_tozero ( variable, nbdim )
220!---------------------------------------------------------------------------------------------------
221    type(Agrif_Variable), intent(inout) :: variable     !< Variable
222    integer, intent(in)                 :: nbdim        !< Dimension of the array you want to reset
223!
224    select case (nbdim)
225    case (1) ; call Agrif_set_array_tozero_1D(variable%array1)
226    case (2) ; call Agrif_set_array_tozero_2D(variable%array2)
227    case (3) ; call Agrif_set_array_tozero_3D(variable%array3)
228    case (4) ; call Agrif_set_array_tozero_4D(variable%array4)
229    case (5) ; call Agrif_set_array_tozero_5D(variable%array5)
230    case (6) ; call Agrif_set_array_tozero_6D(variable%array6)
231    end select
232!---------------------------------------------------------------------------------------------------
233contains
234!---------------------------------------------------------------------------------------------------
235    subroutine Agrif_set_array_tozero_1D ( array )
236        real, dimension(:),           intent(out) :: array
237        array = 0.
238    end subroutine Agrif_set_array_tozero_1D
239!
240    subroutine Agrif_set_array_tozero_2D ( array )
241        real, dimension(:,:),         intent(out) :: array
242        array = 0.
243    end subroutine Agrif_set_array_tozero_2D
244!
245    subroutine Agrif_set_array_tozero_3D ( array )
246        real, dimension(:,:,:),       intent(out) :: array
247        array = 0.
248    end subroutine Agrif_set_array_tozero_3D
249!
250    subroutine Agrif_set_array_tozero_4D ( array )
251        real, dimension(:,:,:,:),     intent(out) :: array
252        array = 0.
253    end subroutine Agrif_set_array_tozero_4D
254!
255    subroutine Agrif_set_array_tozero_5D ( array )
256        real, dimension(:,:,:,:,:),   intent(out) :: array
257        array = 0.
258    end subroutine Agrif_set_array_tozero_5D
259!
260    subroutine Agrif_set_array_tozero_6D ( array )
261        real, dimension(:,:,:,:,:,:), intent(out) :: array
262        array = 0.
263    end subroutine Agrif_set_array_tozero_6D
264!---------------------------------------------------------------------------------------------------
265end subroutine Agrif_var_set_array_tozero
266!===================================================================================================
267!
268!===================================================================================================
269!  subroutine Agrif_var_copy_array
270!
271!> Copy a part of data array from var2 to var1
272!---------------------------------------------------------------------------------------------------
273subroutine Agrif_var_copy_array ( var1, inf1, sup1, var2, inf2, sup2, nbdim )
274!---------------------------------------------------------------------------------------------------
275    type(Agrif_Variable),      intent(inout) :: var1    !< Modified variable
276    integer, dimension(nbdim), intent(in)    :: inf1    !< Lower boundary for var1
277    integer, dimension(nbdim), intent(in)    :: sup1    !< Upper boundary for var1
278    type(Agrif_Variable),      intent(in)    :: var2    !< Input variable
279    integer, dimension(nbdim), intent(in)    :: inf2    !< Lower boundary for var2
280    integer, dimension(nbdim), intent(in)    :: sup2    !< Upper boundary for var2
281    integer,                   intent(in)    :: nbdim   !< Dimension of the array
282!
283    select case (nbdim)
284    case (1) ; var1%array1(inf1(1):sup1(1)) = var2%array1(inf2(1):sup2(1))
285    case (2) ; call Agrif_copy_array_2d( var1%array2, var2%array2, &
286                                         lbound(var1%array2), lbound(var2%array2), inf1, sup1, inf2, sup2 )
287    case (3) ; call Agrif_copy_array_3d( var1%array3, var2%array3, &
288                                         lbound(var1%array3), lbound(var2%array3), inf1, sup1, inf2, sup2 )
289    case (4) ; call Agrif_copy_array_4d( var1%array4, var2%array4, &
290                                         lbound(var1%array4), lbound(var2%array4), inf1, sup1, inf2, sup2 )
291    case (5) ; var1%array5(inf1(1):sup1(1),   &
292                           inf1(2):sup1(2),   &
293                           inf1(3):sup1(3),   &
294                           inf1(4):sup1(4),   &
295                           inf1(5):sup1(5)) = var2%array5(inf2(1):sup2(1), &
296                                                          inf2(2):sup2(2), &
297                                                          inf2(3):sup2(3), &
298                                                          inf2(4):sup2(4), &
299                                                          inf2(5):sup2(5))
300    case (6) ; var1%array6(inf1(1):sup1(1),   &
301                           inf1(2):sup1(2),   &
302                           inf1(3):sup1(3),   &
303                           inf1(4):sup1(4),   &
304                           inf1(5):sup1(5),   &
305                           inf1(6):sup1(6)) = var2%array6(inf2(1):sup2(1), &
306                                                          inf2(2):sup2(2), &
307                                                          inf2(3):sup2(3), &
308                                                          inf2(4):sup2(4), &
309                                                          inf2(5):sup2(5), &
310                                                          inf2(6):sup2(6))
311    end select
312!---------------------------------------------------------------------------------------------------
313contains
314!---------------------------------------------------------------------------------------------------
315    subroutine Agrif_copy_array_2d ( tabout, tabin, l, m, inf1, sup1, inf2, sup2 )
316        integer, dimension(2), intent(in)   :: l, m
317        integer, dimension(2), intent(in)   :: inf1, sup1
318        integer, dimension(2), intent(in)   :: inf2, sup2
319        real, dimension(l(1):,l(2):), intent(out) :: tabout
320        real, dimension(m(1):,m(2):), intent(in)  :: tabin
321        tabout(inf1(1):sup1(1), &
322               inf1(2):sup1(2)) = tabin(inf2(1):sup2(1), &
323                                        inf2(2):sup2(2))
324    end subroutine Agrif_copy_array_2d
325!
326    subroutine Agrif_copy_array_3d ( tabout, tabin, l, m, inf1, sup1, inf2, sup2 )
327        integer, dimension(3), intent(in)   :: l, m
328        integer, dimension(3), intent(in)   :: inf1, sup1
329        integer, dimension(3), intent(in)   :: inf2,sup2
330        real, dimension(l(1):,l(2):,l(3):), intent(out) :: tabout
331        real, dimension(m(1):,m(2):,m(3):), intent(in)  :: tabin
332        tabout(inf1(1):sup1(1), &
333               inf1(2):sup1(2), &
334               inf1(3):sup1(3)) = tabin(inf2(1):sup2(1), &
335                                        inf2(2):sup2(2), &
336                                        inf2(3):sup2(3))
337    end subroutine Agrif_copy_array_3d
338!
339    subroutine Agrif_copy_array_4d ( tabout, tabin, l, m, inf1, sup1, inf2, sup2 )
340        integer, dimension(4), intent(in)   :: l, m
341        integer, dimension(4), intent(in)   :: inf1, sup1
342        integer, dimension(4), intent(in)   :: inf2, sup2
343        real, dimension(l(1):,l(2):,l(3):,l(4):), intent(out) :: tabout
344        real, dimension(m(1):,m(2):,m(3):,m(4):), intent(in)  :: tabin
345        tabout(inf1(1):sup1(1), &
346               inf1(2):sup1(2), &
347               inf1(3):sup1(3), &
348               inf1(4):sup1(4)) = tabin(inf2(1):sup2(1), &
349                                        inf2(2):sup2(2), &
350                                        inf2(3):sup2(3), &
351                                        inf2(4):sup2(4))
352    end subroutine Agrif_copy_array_4d
353!---------------------------------------------------------------------------------------------------
354end subroutine Agrif_var_copy_array
355!===================================================================================================
356!
357!===================================================================================================
358!  subroutine Agrif_var_full_copy_array
359!
360!> Copy the full data array from var2 to var1
361!---------------------------------------------------------------------------------------------------
362subroutine Agrif_var_full_copy_array ( var1, var2, nbdim )
363!---------------------------------------------------------------------------------------------------
364    type(Agrif_Variable), intent(inout) :: var1    !< Modified variable
365    type(Agrif_Variable), intent(in)    :: var2    !< Input variable
366    integer,              intent(in)    :: nbdim   !< Dimension of the array
367!
368    select case (nbdim)
369    case (1) ; var1 % array1 = var2 % array1
370    case (2) ; var1 % array2 = var2 % array2
371    case (3) ; var1 % array3 = var2 % array3
372    case (4) ; var1 % array4 = var2 % array4
373    case (5) ; var1 % array5 = var2 % array5
374    case (6) ; var1 % array6 = var2 % array6
375    end select
376!---------------------------------------------------------------------------------------------------
377end subroutine Agrif_var_full_copy_array
378!===================================================================================================
379!
380!===================================================================================================
381!  subroutine GiveAgrif_SpecialValueToTab_mpi
382!
383!> Copy \b value in data array \b var2 where it is present in \b var1.
384!---------------------------------------------------------------------------------------------------
385subroutine GiveAgrif_SpecialValueToTab_mpi ( var1, var2, bounds, value, nbdim )
386!---------------------------------------------------------------------------------------------------
387    type(Agrif_Variable),      intent(in)    :: var1    !< Modified variable
388    type(Agrif_Variable),      intent(inout) :: var2    !< Input variable
389    integer, dimension(:,:,:), intent(in)    :: bounds  !< Bound for both arrays
390    real,                      intent(in)    :: value   !< Special value
391    integer,                   intent(in)    :: nbdim   !< Dimension of the array
392!
393    select case (nbdim)
394    case (1)
395        where (var1 % array1(bounds(1,1,2):bounds(1,2,2)) == value )
396            var2 % array1(bounds(1,1,1):bounds(1,2,1)) = value
397        end where
398    case (2)
399        where (var1 % array2(bounds(1,1,2):bounds(1,2,2), &
400                             bounds(2,1,2):bounds(2,2,2)) == value)
401            var2 % array2(bounds(1,1,1):bounds(1,2,1), &
402                          bounds(2,1,1):bounds(2,2,1)) = value
403        end where
404    case (3)
405        where (var1 % array3(bounds(1,1,2):bounds(1,2,2), &
406                             bounds(2,1,2):bounds(2,2,2), &
407                             bounds(3,1,2):bounds(3,2,2)) == value)
408            var2 % array3(bounds(1,1,1):bounds(1,2,1), &
409                          bounds(2,1,1):bounds(2,2,1), &
410                          bounds(3,1,1):bounds(3,2,1)) = value
411        end where
412    case (4)
413        where (var1 % array4(bounds(1,1,2):bounds(1,2,2), &
414                             bounds(2,1,2):bounds(2,2,2), &
415                             bounds(3,1,2):bounds(3,2,2), &
416                             bounds(4,1,2):bounds(4,2,2)) == value)
417            var2 % array4(bounds(1,1,1):bounds(1,2,1), &
418                          bounds(2,1,1):bounds(2,2,1), &
419                          bounds(3,1,1):bounds(3,2,1), &
420                          bounds(4,1,1):bounds(4,2,1)) = value
421        end where
422    case (5)
423        where (var1 % array5(bounds(1,1,2):bounds(1,2,2), &
424                             bounds(2,1,2):bounds(2,2,2), &
425                             bounds(3,1,2):bounds(3,2,2), &
426                             bounds(4,1,2):bounds(4,2,2), &
427                             bounds(5,1,2):bounds(5,2,2)) == value)
428            var2 % array5(bounds(1,1,1):bounds(1,2,1), &
429                          bounds(2,1,1):bounds(2,2,1), &
430                          bounds(3,1,1):bounds(3,2,1), &
431                          bounds(4,1,1):bounds(4,2,1), &
432                          bounds(5,1,1):bounds(5,2,1)) = value
433        end where
434    case (6)
435        where (var1 % array6(bounds(1,1,2):bounds(1,2,2), &
436                             bounds(2,1,2):bounds(2,2,2), &
437                             bounds(3,1,2):bounds(3,2,2), &
438                             bounds(4,1,2):bounds(4,2,2), &
439                             bounds(5,1,2):bounds(5,2,2), &
440                             bounds(6,1,2):bounds(6,2,2)) == value)
441            var2 % array6(bounds(1,1,1):bounds(1,2,1), &
442                          bounds(2,1,1):bounds(2,2,1), &
443                          bounds(3,1,1):bounds(3,2,1), &
444                          bounds(4,1,1):bounds(4,2,1), &
445                          bounds(5,1,1):bounds(5,2,1), &
446                          bounds(6,1,1):bounds(6,2,1)) = value
447        end where
448    end select
449!---------------------------------------------------------------------------------------------------
450end subroutine GiveAgrif_SpecialValueToTab_mpi
451!===================================================================================================
452!
453! no more used ???
454#if 0
455!===================================================================================================
456!  subroutine GiveAgrif_SpecialValueToTab
457!---------------------------------------------------------------------------------------------------
458subroutine GiveAgrif_SpecialValueToTab ( var1, var2, &
459                                         lower, upper, Value, nbdim )
460!---------------------------------------------------------------------------------------------------
461    TYPE(Agrif_Variable),      pointer      :: var1
462    TYPE(Agrif_Variable),      pointer      :: var2
463    INTEGER,                   intent(in)   :: nbdim
464    INTEGER, DIMENSION(nbdim), intent(in)   :: lower, upper
465    REAL,                      intent(in)   :: Value
466!
467    select case (nbdim)
468    case (1)
469        where (var1 % array1( lower(1):upper(1)) == Value)
470            var2 % array1(lower(1):upper(1)) = Value
471        end where
472    case (2)
473        where (var1 % array2( lower(1):upper(1), &
474                                   lower(2):upper(2)) == Value)
475            var2 % array2(lower(1):upper(1), &
476                               lower(2):upper(2)) = Value
477        end where
478    case (3)
479        where (var1 % array3( lower(1):upper(1), &
480                                   lower(2):upper(2), &
481                                   lower(3):upper(3)) == Value)
482            var2 % array3(lower(1):upper(1), &
483                               lower(2):upper(2), &
484                               lower(3):upper(3)) = Value
485        end where
486    case (4)
487        where (var1 % array4( lower(1):upper(1), &
488                                   lower(2):upper(2), &
489                                   lower(3):upper(3), &
490                                   lower(4):upper(4)) == Value)
491            var2 % array4(lower(1):upper(1), &
492                               lower(2):upper(2), &
493                               lower(3):upper(3), &
494                               lower(4):upper(4)) = Value
495        end where
496    case (5)
497        where (var1 % array5( lower(1):upper(1), &
498                                   lower(2):upper(2), &
499                                   lower(3):upper(3), &
500                                   lower(4):upper(4), &
501                                   lower(5):upper(5)) == Value)
502            var2 % array5(lower(1):upper(1), &
503                               lower(2):upper(2), &
504                               lower(3):upper(3), &
505                               lower(4):upper(4), &
506                               lower(5):upper(5)) = Value
507        end where
508    case (6)
509        where (var1 % array6( lower(1):upper(1), &
510                                   lower(2):upper(2), &
511                                   lower(3):upper(3), &
512                                   lower(4):upper(4), &
513                                   lower(5):upper(5), &
514                                   lower(6):upper(6)) == Value)
515            var2 % array6(lower(1):upper(1), &
516                               lower(2):upper(2), &
517                               lower(3):upper(3), &
518                               lower(4):upper(4), &
519                               lower(5):upper(5), &
520                               lower(6):upper(6)) = Value
521        end where
522    end select
523!---------------------------------------------------------------------------------------------------
524end subroutine GiveAgrif_SpecialValueToTab
525!===================================================================================================
526#endif
527!
528#if defined AGRIF_MPI
529!===================================================================================================
530!  subroutine Agrif_var_replace_value
531!
532!> Replace \b value by \var2 content in \var1 data array.
533!---------------------------------------------------------------------------------------------------
534subroutine Agrif_var_replace_value ( var1, var2, lower, upper, value, nbdim )
535!---------------------------------------------------------------------------------------------------
536    type(Agrif_Variable),      intent(inout) :: var1    !< Modified variable
537    type(Agrif_Variable),      intent(in)    :: var2    !< Input variable
538    integer, dimension(nbdim), intent(in)    :: lower   !< Lower bound
539    integer, dimension(nbdim), intent(in)    :: upper   !< Upper bound
540    real,                      intent(in)    :: value   !< Special value
541    integer,                   intent(in)    :: nbdim   !< Dimension of the array
542!
543    integer :: i,j,k,l,m,n
544!
545    select case (nbdim)
546    case (1)
547        do i = lower(1),upper(1)
548            if (var1%array1(i) == value) then
549                var1%array1(i) = var2%array1(i)
550            endif
551        enddo
552    case (2)
553        do j = lower(2),upper(2)
554        do i = lower(1),upper(1)
555            if (var1%array2(i,j) == value) then
556                var1%array2(i,j) = var2%array2(i,j)
557            endif
558        enddo
559        enddo
560    case (3)
561        do k = lower(3),upper(3)
562        do j = lower(2),upper(2)
563        do i = lower(1),upper(1)
564            if (var1%array3(i,j,k) == value) then
565                var1%array3(i,j,k) = var2%array3(i,j,k)
566            endif
567        enddo
568        enddo
569        enddo
570    case (4)
571        do l = lower(4),upper(4)
572        do k = lower(3),upper(3)
573        do j = lower(2),upper(2)
574        do i = lower(1),upper(1)
575            if (var1%array4(i,j,k,l) == value) then
576                var1%array4(i,j,k,l) = var2%array4(i,j,k,l)
577            endif
578        enddo
579        enddo
580        enddo
581        enddo
582    case (5)
583        do m = lower(5),upper(5)
584        do l = lower(4),upper(4)
585        do k = lower(3),upper(3)
586        do j = lower(2),upper(2)
587        do i = lower(1),upper(1)
588            if (var1%array5(i,j,k,l,m) == value) then
589                var1%array5(i,j,k,l,m) = var2%array5(i,j,k,l,m)
590            endif
591        enddo
592        enddo
593        enddo
594        enddo
595        enddo
596    case (6)
597        do n = lower(6),upper(6)
598        do m = lower(5),upper(5)
599        do l = lower(4),upper(4)
600        do k = lower(3),upper(3)
601        do j = lower(2),upper(2)
602        do i = lower(1),upper(1)
603            if (var1%array6(i,j,k,l,m,n) == value) then
604                var1%array6(i,j,k,l,m,n) = var2%array6(i,j,k,l,m,n)
605            endif
606        enddo
607        enddo
608        enddo
609        enddo
610        enddo
611        enddo
612    end select
613!---------------------------------------------------------------------------------------------------
614end subroutine Agrif_var_replace_value
615!===================================================================================================
616#endif
617!
618!===================================================================================================
619!  subroutine PreProcessToInterpOrUpdate
620!---------------------------------------------------------------------------------------------------
621subroutine PreProcessToInterpOrUpdate ( parent, child,              &
622                                        nb_child, ub_child,         &
623                                        lb_child, lb_parent,        &
624                                        s_child,  s_parent,         &
625                                        ds_child, ds_parent, nbdim, interp )
626!---------------------------------------------------------------------------------------------------
627    type(Agrif_Variable), pointer, intent(in)   :: parent       !< Variable on the parent grid
628    type(Agrif_Variable), pointer, intent(in)   :: child        !< Variable on the child grid
629    integer, dimension(6), intent(out)          :: nb_child     !< Number of cells on the child grid
630    integer, dimension(6), intent(out)          :: ub_child     !< Upper bound on the child grid
631    integer, dimension(6), intent(out)          :: lb_child     !< Lower bound on the child grid
632    integer, dimension(6), intent(out)          :: lb_parent    !< Lower bound on the parent grid
633    real, dimension(6),    intent(out)          :: s_child      !< Child  grid position (s_root = 0)
634    real, dimension(6),    intent(out)          :: s_parent     !< Parent grid position (s_root = 0)
635    real, dimension(6),    intent(out)          :: ds_child     !< Child  grid dx (ds_root = 1)
636    real, dimension(6),    intent(out)          :: ds_parent    !< Parent grid dx (ds_root = 1)
637    integer,               intent(out)          :: nbdim        !< Number of dimensions
638    logical,               intent(in)           :: interp       !< .true. if preprocess for interpolation, \n
639                                                                !! .false. if preprocess for update
640!
641    type(Agrif_Variable), pointer :: root_var
642    type(Agrif_Grid),     pointer :: Agrif_Child_Gr
643    type(Agrif_Grid),     pointer :: Agrif_Parent_Gr
644    integer :: n
645!
646    Agrif_Child_Gr  => Agrif_Curgrid
647    Agrif_Parent_Gr => Agrif_Curgrid % parent
648!
649    root_var => child % root_var
650!
651!   Number of dimensions of the current grid
652    nbdim = root_var % nbdim
653!
654    do n = 1,nbdim
655!
656!       Value of interptab(n) can be either x,y,z or N for a no space dimension
657        select case(root_var % interptab(n))
658!
659        case('x')
660!
661            lb_child(n)  = root_var % point(n)
662            lb_parent(n) = root_var % point(n)
663            nb_child(n)  = Agrif_Child_Gr % nb(1)
664            s_child(n)   = Agrif_Child_Gr  % Agrif_x(1)
665            s_parent(n)  = Agrif_Parent_Gr % Agrif_x(1)
666            ds_child(n)  = Agrif_Child_Gr  % Agrif_dx(1)
667            ds_parent(n) = Agrif_Parent_Gr % Agrif_dx(1)
668!
669            if ( root_var % posvar(n) == 1 ) then
670                ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(1)
671            else
672                ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(1) - 1
673                s_child(n)  = s_child(n)  + 0.5*ds_child(n)
674                s_parent(n) = s_parent(n) + 0.5*ds_parent(n)
675            endif
676!
677        case('y')
678!
679            lb_child(n)  = root_var % point(n)
680            lb_parent(n) = root_var % point(n)
681            nb_child(n)  = Agrif_Child_Gr % nb(2)
682            s_child(n)   = Agrif_Child_Gr  % Agrif_x(2)
683            s_parent(n)  = Agrif_Parent_Gr % Agrif_x(2)
684            ds_child(n)  = Agrif_Child_Gr  % Agrif_dx(2)
685            ds_parent(n) = Agrif_Parent_Gr % Agrif_dx(2)
686!
687            if (root_var % posvar(n)==1) then
688                ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(2)
689            else
690                ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(2) - 1
691                s_child(n)  = s_child(n)  + 0.5*ds_child(n)
692                s_parent(n) = s_parent(n) + 0.5*ds_parent(n)
693            endif
694!
695        case('z')
696!
697            lb_child(n)  = root_var % point(n)
698            lb_parent(n) = root_var % point(n)
699            nb_child(n)  = Agrif_Child_Gr % nb(3)
700            s_child(n)   = Agrif_Child_Gr  % Agrif_x(3)
701            s_parent(n)  = Agrif_Parent_Gr % Agrif_x(3)
702            ds_child(n)  = Agrif_Child_Gr  % Agrif_dx(3)
703            ds_parent(n) = Agrif_Parent_Gr % Agrif_dx(3)
704!
705            if (root_var % posvar(n)==1) then
706                ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(3)
707            else
708                ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(3) - 1
709                s_child(n)  = s_child(n)  + 0.5*ds_child(n)
710                s_parent(n) = s_parent(n) + 0.5*ds_parent(n)
711            endif
712!
713        case('N') ! No space dimension
714!
715!       The next coefficients are calculated in order to do a simple copy of
716!       values of the grid variable when the interpolation routine is
717!       called for this dimension.
718!
719            if (interp) then
720                call Agrif_get_var_bounds(parent, lb_child(n), ub_child(n), n)
721                nb_child(n) = parent % ub(n) - parent % lb(n)
722            else
723                call Agrif_get_var_bounds(child,  lb_child(n), ub_child(n), n)
724                nb_child(n) = child % ub(n) - child % lb(n)
725            endif
726!
727!           No interpolation but only a copy of the values of the grid variable
728            lb_parent(n) = lb_child(n)
729            s_child(n)   = 0.
730            s_parent(n)  = 0.
731            ds_child(n)  = 1.
732            ds_parent(n) = 1.
733!
734        end select
735!
736    enddo
737!---------------------------------------------------------------------------------------------------
738end subroutine PreProcessToInterpOrUpdate
739!===================================================================================================
740!
741#if defined AGRIF_MPI
742!===================================================================================================
743!   subroutine Agrif_GetLocalBoundaries
744!---------------------------------------------------------------------------------------------------
745subroutine Agrif_GetLocalBoundaries ( tab1, tab2, coord, lb, ub, deb, fin )
746!---------------------------------------------------------------------------------------------------
747    integer, intent(in)     ::  tab1
748    integer, intent(in)     ::  tab2
749    integer, intent(in)     ::  coord
750    integer, intent(in)     ::  lb, ub
751    integer, intent(out)    ::  deb, fin
752!
753    integer :: imin, imax
754    integer :: i1, i2
755!
756    call Agrif_InvLoc(lb, AGRIF_ProcRank, coord, imin)
757    call Agrif_InvLoc(ub, AGRIF_ProcRank, coord, imax)
758!
759    if ( imin > tab2 ) then
760        i1 = imax - imin
761    else
762        i1 = max(tab1 - imin,0)
763    endif
764!
765    if (imax < tab1) then
766        i2 = -(imax - imin)
767    else
768        i2 = min(tab2 - imax,0)
769    endif
770!
771    deb = lb + i1
772    fin = ub + i2
773!---------------------------------------------------------------------------------------------------
774end subroutine Agrif_GetLocalBoundaries
775!===================================================================================================
776!
777!===================================================================================================
778!  subroutine Agrif_GlobalToLocalBounds
779!
780!> For a global index located on the current processor, tabarray gives the corresponding local index
781!---------------------------------------------------------------------------------------------------
782subroutine Agrif_GlobalToLocalBounds ( locbounds, lb_var, ub_var, lb_glob, ub_glob,    &
783                                       coords, nbdim, rank, member )
784!---------------------------------------------------------------------------------------------------
785    integer, dimension(nbdim,2,2), intent(out)   :: locbounds   !< Local values of \b lb_glob and \b ub_glob
786    integer, dimension(nbdim),     intent(in)    :: lb_var      !< Local lower boundary on the current processor
787    integer, dimension(nbdim),     intent(in)    :: ub_var      !< Local upper boundary on the current processor
788    integer, dimension(nbdim),     intent(in)    :: lb_glob     !< Global lower boundary
789    integer, dimension(nbdim),     intent(in)    :: ub_glob     !< Global upper boundary
790    integer, dimension(nbdim),     intent(in)    :: coords
791    integer,                       intent(in)    :: nbdim       !< Dimension of the array
792    integer,                       intent(in)    :: rank        !< Rank of the processor
793    logical,                       intent(out)   :: member
794!
795    integer     :: i, i1, k
796    integer     :: nbloc(nbdim)
797!
798    locbounds(:,1,:) =  HUGE(1)
799    locbounds(:,2,:) = -HUGE(1)
800!
801    nbloc = 0
802!
803    do i = 1,nbdim
804!
805        call Agrif_InvLoc(lb_var(i), rank, coords(i), i1)
806!
807        do k = lb_glob(i)+lb_var(i)-i1,ub_glob(i)+lb_var(i)-i1
808!
809            if ( (k >= lb_var(i)) .AND. (k <= ub_var(i)) ) then
810                nbloc(i) = 1
811                locbounds(i,1,1) = min(locbounds(i,1,1),k-lb_var(i)+i1)
812                locbounds(i,2,1) = max(locbounds(i,2,1),k-lb_var(i)+i1)
813
814                locbounds(i,1,2) = min(locbounds(i,1,2),k)
815                locbounds(i,2,2) = max(locbounds(i,2,2),k)
816            endif
817        enddo
818    enddo
819
820    member = ( sum(nbloc) == nbdim )
821!---------------------------------------------------------------------------------------------------
822end subroutine Agrif_GlobalToLocalBounds
823!===================================================================================================
824#endif
825!
826end module Agrif_Arrays
Note: See TracBrowser for help on using the repository browser.