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/dev/AGRIF_FILES – NEMO

source: vendors/AGRIF/dev/AGRIF_FILES/modarrays.F90 @ 15504

Last change on this file since 15504 was 14975, checked in by jchanut, 3 years ago

#2638, merge new AGRIF library into trunk

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