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

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

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

Last change on this file since 10087 was 10087, checked in by rblod, 6 years ago

update AGRIF library

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