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

source: vendors/AGRIF/dev_r12970_AGRIF_CMEMS/AGRIF_FILES/modarrays.F90 @ 13370

Last change on this file since 13370 was 13370, checked in by jchanut, 4 years ago

#2222, changes to allow less vertical levels over child grid than parent

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