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 @ 13144

Last change on this file since 13144 was 13027, checked in by rblod, 4 years ago

New AGRIF library, see ticket #2129

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