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 vendors/AGRIF/CMEMS_2020/AGRIF_FILES – NEMO

source: vendors/AGRIF/CMEMS_2020/AGRIF_FILES/modarrays.F90 @ 10725

Last change on this file since 10725 was 10725, checked in by rblod, 5 years ago

Update agrif library and conv see ticket #2129

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