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.
modupdate.F90 in branches/2017/dev_r8126_UKMO_AGRIF_vert_interp/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES – NEMO

source: branches/2017/dev_r8126_UKMO_AGRIF_vert_interp/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modupdate.F90 @ 8138

Last change on this file since 8138 was 8138, checked in by timgraham, 7 years ago

Modifications to AGRIF_FILES as received from Laurent (need to check that these are definitely needed)

  • Property svn:keywords set to Id
File size: 101.8 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!> Module Agrif_Update
24!>
25!> This module  contains procedures to update a parent grid from its child grids.
26!
27module Agrif_Update
28!
29    use Agrif_UpdateBasic
30    use Agrif_Arrays
31    use Agrif_CurgridFunctions
32    use Agrif_Mask
33#if defined AGRIF_MPI
34    use Agrif_Mpp
35#endif
36!
37    implicit none
38!
39    logical, private :: precomputedone(7) = .FALSE.
40!
41contains
42!
43!===================================================================================================
44!  subroutine Agrif_UpdateVariable
45!
46!> subroutine to set arguments for Agrif_UpdatenD
47!---------------------------------------------------------------------------------------------------
48subroutine Agrif_UpdateVariable ( parent, child, updateinf, updatesup, procname )
49!---------------------------------------------------------------------------------------------------
50    type(Agrif_Variable),  pointer      :: parent       !< Variable on the parent grid
51    type(Agrif_Variable),  pointer      :: child        !< Variable on the child grid
52    integer, dimension(6), intent(in)   :: updateinf    !< First positions where interpolations are calculated
53    integer, dimension(6), intent(in)   :: updatesup    !< Last  positions where interpolations are calculated
54    procedure()                         :: procname     !< Data recovery procedure
55!---------------------------------------------------------------------------------------------------
56    integer, dimension(6) :: nb_child           ! Number of cells on the child grid
57    integer, dimension(6) :: lb_child
58    integer, dimension(6) :: ub_child
59    integer, dimension(6) :: lb_parent
60    real   , dimension(6) ::  s_child           ! Child  grid position (s_root = 0)
61    real   , dimension(6) ::  s_parent          ! Parent grid position (s_root = 0)
62    real   , dimension(6) :: ds_child           ! Child  grid dx (ds_root = 1)
63    real   , dimension(6) :: ds_parent          ! Parent grid dx (ds_root = 1)
64    logical, dimension(6) :: do_update          ! Indicates if we perform update for each dimension
65    integer, dimension(6) :: posvar             ! Position of the variable on the cell (1 or 2)
66    integer, dimension(6) :: oldparentlbound, oldparentubound
67    integer               :: n, nbdim
68    logical               :: wholeupdate
69    type(Agrif_Variable), pointer :: root   ! Variable on the root grid
70!
71    root => child % root_var
72    nbdim = root % nbdim
73!
74    call PreProcessToInterpOrUpdate( parent,   child,       &
75                                     nb_child, ub_child,    &
76                                     lb_child, lb_parent,   &
77                                      s_child,  s_parent,   &
78                                     ds_child, ds_parent, nbdim, interp=.false. )
79!
80    do_update(:) = .true.
81    posvar(1:nbdim) = root % posvar(1:nbdim)
82!
83    do n = 1,nbdim
84!
85        if ( root % interptab(n) == 'N' ) then
86            posvar(n) = 1
87            do_update(n) = .false.
88            oldparentlbound(n) = parent % lb(n)
89            oldparentubound(n) = parent % ub(n)
90            parent % lb(n) = child % lb(n)
91            parent % ub(n) = child % ub(n)
92        end if
93!
94    enddo
95
96    wholeupdate = .FALSE.
97!
98    do n = 1,nbdim
99        if ( do_update(n) ) then
100            if ( (updateinf(n) > updatesup(n)) .OR. &
101                ((updateinf(n) == -99) .AND. (updatesup(n) == -99)) &
102               ) then
103                wholeupdate = .TRUE.
104            endif
105        endif
106    enddo
107!
108    IF (wholeupdate) THEN
109        call Agrif_UpdateWhole(parent, child,               &
110                updateinf(1:nbdim), updatesup(1:nbdim),     &
111                lb_child(1:nbdim), lb_parent(1:nbdim),      &
112                nb_child(1:nbdim), posvar(1:nbdim),         &
113                do_update(1:nbdim),                         &
114                s_child(1:nbdim),   s_parent(1:nbdim),      &
115                ds_child(1:nbdim), ds_parent(1:nbdim), nbdim, procname)
116    ELSE
117        call Agrif_UpdateBcnD(parent, child,                &
118                updateinf(1:nbdim), updatesup(1:nbdim),     &
119                lb_child(1:nbdim), lb_parent(1:nbdim),      &
120                nb_child(1:nbdim), posvar(1:nbdim),         &
121                do_update(1:nbdim),                         &
122                s_child(1:nbdim),   s_parent(1:nbdim),      &
123                ds_child(1:nbdim), ds_parent(1:nbdim), nbdim, procname)
124    ENDIF
125!
126    do n = 1,nbdim
127!
128        if ( root % interptab(n) == 'N' ) then  ! No space DIMENSION
129            parent % lb(n) = oldparentlbound(n)
130            parent % ub(n) = oldparentubound(n)
131        end if
132!
133    enddo
134!---------------------------------------------------------------------------------------------------
135end subroutine Agrif_UpdateVariable
136!===================================================================================================
137!
138!===================================================================================================
139!  subroutine Agrif_UpdateWhole
140!---------------------------------------------------------------------------------------------------
141subroutine Agrif_UpdateWhole ( parent, child, uinf, usup,       &
142                               lb_child, lb_parent,             &
143                               nb_child, posvar,                &
144                               do_update,                       &
145                               s_child,   s_parent,             &
146                               ds_child, ds_parent, nbdim, procname )
147!---------------------------------------------------------------------------------------------------
148#if defined AGRIF_MPI
149    include 'mpif.h'
150#endif
151!
152    type(Agrif_Variable),      pointer    :: parent         !< Variable on the parent grid
153    type(Agrif_Variable),      pointer    :: child          !< Variable on the child grid
154    integer, dimension(nbdim), intent(in) :: uinf           !< First positions where interpolations are calculated
155    integer, dimension(nbdim), intent(in) :: usup           !< Last  positions where interpolations are calculated
156    integer,                   intent(in) :: nbdim          !< Number of dimensions of the grid variable
157    integer, dimension(nbdim), intent(in) :: lb_child       !< Index of the first point inside the domain for the parent grid variable
158    integer, dimension(nbdim), intent(in) :: lb_parent      !< Index of the first point inside the domain for the child grid variable
159    integer, dimension(nbdim), intent(in) :: nb_child       !< Number of cells of the child grid
160    integer, dimension(nbdim), intent(in) :: posvar         !< Position of the variable on the cell (1 or 2)
161    logical, dimension(nbdim), intent(in) :: do_update      !< Indicates if we update for each dimension
162    real,    dimension(nbdim), intent(in) :: s_child        !< Positions of the child grid
163    real,    dimension(nbdim), intent(in) :: s_parent       !< Positions of the parent grid
164    real,    dimension(nbdim), intent(in) :: ds_child       !< Space steps of the child grid
165    real,    dimension(nbdim), intent(in) :: ds_parent      !< Space steps of the parent grid
166    procedure()                           :: procname       !< Data recovery procedure
167!
168    integer, dimension(nbdim)     :: type_update ! Type of update (copy or average)
169    integer, dimension(nbdim,2)   :: lubglob
170    integer, dimension(nbdim,2,2) :: indtab      ! limits of the child grid that will be used in the update scheme
171    integer, dimension(nbdim,2,2) :: indtruetab  ! grid variable where boundary conditions are
172    integer :: coeffraf, i
173    integer :: uinfloc, usuploc
174!
175    type_update = child % root_var % type_update(1:nbdim)
176!
177    do i = 1, nbdim
178!
179        if ( do_update(i) ) then
180!
181            coeffraf = nint(ds_parent(i)/ds_child(i))
182            uinfloc = 0
183            usuploc = nb_child(i)/coeffraf - 1
184
185            IF (posvar(i) == 1) THEN
186                usuploc = usuploc - 1
187            ENDIF
188
189            IF (uinf(i) > usup(i)) THEN
190                uinfloc = uinf(i)
191                usuploc = usuploc - uinf(i)
192            ENDIF
193
194            indtab(i,1,1) = lb_child(i) + (uinfloc + 1) * coeffraf
195            indtab(i,1,2) = lb_child(i) + (usuploc + 1) * coeffraf
196
197            IF ( posvar(i) == 1 ) THEN
198                IF ( type_update(i) == Agrif_Update_Full_Weighting ) THEN
199                    indtab(i,1,1) = indtab(i,1,1) - (coeffraf - 1)
200                    indtab(i,1,2) = indtab(i,1,2) + (coeffraf - 1)
201                ELSE IF ( type_update(i) /= Agrif_Update_Copy ) THEN
202                    indtab(i,1,1) = indtab(i,1,1) - coeffraf / 2
203                    indtab(i,1,2) = indtab(i,1,2) + coeffraf / 2
204                ENDIF
205            ELSE
206                indtab(i,1,1) = indtab(i,1,1) - coeffraf
207                indtab(i,1,2) = indtab(i,1,2) - 1
208    ! at this point, indices are OK for an average
209                IF ( type_update(i) == Agrif_Update_Full_Weighting ) THEN
210                    indtab(i,1,1) = indtab(i,1,1) - coeffraf / 2
211                    indtab(i,1,2) = indtab(i,1,2) + coeffraf / 2
212                ENDIF
213            ENDIF
214!
215        else    ! IF ( .not.do_update(i) ) THEN
216!
217            if ( posvar(i) == 1 ) then
218                indtab(i,1,1) = lb_child(i)
219                indtab(i,1,2) = lb_child(i) + nb_child(i)
220            else
221                indtab(i,1,1) = lb_child(i)
222                indtab(i,1,2) = lb_child(i) + nb_child(i) - 1
223            endif
224!
225        endif
226    enddo
227
228! lubglob contains the global lbound and ubound of the child array
229! lubglob(:,1) : global lbound for each dimension
230! lubglob(:,2) : global lbound for each dimension
231!
232    call Agrif_get_var_global_bounds(child, lubglob, nbdim)
233!
234    indtruetab(1:nbdim,1,1) = max(indtab(1:nbdim,1,1), lubglob(1:nbdim,1))
235    indtruetab(1:nbdim,1,2) = min(indtab(1:nbdim,1,2), lubglob(1:nbdim,2))
236!
237    call Agrif_UpdatenD(type_update, parent, child,                         &
238                        indtruetab(1:nbdim,1,1), indtruetab(1:nbdim,1,2),   &
239                        lb_child(1:nbdim), lb_parent(1:nbdim),              &
240                        s_child(1:nbdim), s_parent(1:nbdim),                &
241                        ds_child(1:nbdim), ds_parent(1:nbdim),              &
242#if defined AGRIF_MPI
243                        posvar, do_update,                                  &
244#endif
245                        nbdim, procname)
246!---------------------------------------------------------------------------------------------------
247end subroutine Agrif_UpdateWhole
248!===================================================================================================
249!
250!===================================================================================================
251!  subroutine Agrif_UpdateBcnd
252!---------------------------------------------------------------------------------------------------
253subroutine Agrif_UpdateBcnd ( parent, child, uinf, usup,        &
254                              lb_child, lb_parent,              &
255                              nb_child, posvar,                 &
256                              do_update,                        &
257                              s_child,  s_parent,               &
258                              ds_child, ds_parent, nbdim, procname )
259!---------------------------------------------------------------------------------------------------
260#if defined AGRIF_MPI
261    include 'mpif.h'
262#endif
263!
264    type(Agrif_Variable), pointer           :: parent       !< Variable on the parent grid
265    type(Agrif_Variable), pointer           :: child        !< Variable on the child grid
266    integer, dimension(nbdim), intent(in)   :: uinf         !< First positions where interpolations are calculated
267    integer, dimension(nbdim), intent(in)   :: usup         !< Last  positions where interpolations are calculated
268    integer                                 :: nbdim        !< Number of dimensions of the grid variable
269    integer, dimension(nbdim), intent(in)   :: lb_child     !< Index of the first point inside the domain for
270                                                            !!   the parent grid variable
271    integer, dimension(nbdim), intent(in)   :: lb_parent    !< Index of the first point inside the domain for
272                                                            !!   the child grid variable
273    integer, dimension(nbdim), intent(in)   :: nb_child     !< Number of cells of the child grid
274    integer, dimension(nbdim), intent(in)   :: posvar       !< Position of the variable on the cell (1 or 2)
275    logical, dimension(nbdim), intent(in)   :: do_update    !< Indicates if we update for each dimension
276    real,    dimension(nbdim), intent(in)   :: s_child      !< Positions of the child grid
277    real,    dimension(nbdim), intent(in)   :: s_parent     !< Positions of the parent grid
278    real,    dimension(nbdim), intent(in)   :: ds_child     !< Space steps of the child grid
279    real,    dimension(nbdim), intent(in)   :: ds_parent    !< Space steps of the parent grid
280    procedure()                             :: procname     !< Data recovery procedure
281!
282    integer,dimension(nbdim)     :: type_update ! Type of update (copy or average)
283    integer,dimension(nbdim,2)   :: lubglob
284    integer                      :: i
285    integer,dimension(nbdim,2,2) :: indtab         ! Arrays indicating the limits of the child
286    integer,dimension(nbdim,2,2) :: indtruetab     ! grid variable where boundary conditions are
287    integer,dimension(nbdim,2,2,nbdim)   :: ptres  ! calculated
288    integer                      :: nb, ndir
289    integer :: coeffraf
290!
291    type_update = child % root_var % type_update(1:nbdim)
292!
293    DO i = 1, nbdim
294        coeffraf = nint(ds_parent(i)/ds_child(i))
295        indtab(i,1,1) = lb_child(i) + (uinf(i) + 1) * coeffraf
296        indtab(i,1,2) = lb_child(i) + (usup(i) + 1) * coeffraf
297
298        indtab(i,2,1) = lb_child(i) + nb_child(i) - (usup(i)+1) *  coeffraf
299        indtab(i,2,2) = lb_child(i) + nb_child(i) - (uinf(i)+1) *  coeffraf
300
301        IF (posvar(i) == 1) THEN
302            IF (type_update(i) == Agrif_Update_Full_Weighting) THEN
303                indtab(i,:,1) = indtab(i,:,1) - (coeffraf - 1)
304                indtab(i,:,2) = indtab(i,:,2) + (coeffraf - 1)
305            ELSE IF (type_update(i) /= Agrif_Update_Copy) THEN
306                indtab(i,:,1) = indtab(i,:,1) - coeffraf / 2
307                indtab(i,:,2) = indtab(i,:,2) + coeffraf / 2
308            ENDIF
309        ELSE
310            indtab(i,1,1) = indtab(i,1,1) - coeffraf
311            indtab(i,1,2) = indtab(i,1,2) - 1
312            indtab(i,2,2) = indtab(i,2,2) + coeffraf - 1
313            IF (type_update(i) == Agrif_Update_Full_Weighting) THEN
314                indtab(i,1,1) = indtab(i,1,1) - coeffraf/2
315                indtab(i,1,2) = indtab(i,1,2) + coeffraf/2
316                indtab(i,2,1) = indtab(i,2,1) - coeffraf/2
317                indtab(i,2,2) = indtab(i,2,2) + coeffraf/2
318            ENDIF
319        ENDIF
320    ENDDO
321!
322    call Agrif_get_var_global_bounds(child,lubglob,nbdim)
323!
324    indtruetab(1:nbdim,1,1) = max(indtab(1:nbdim,1,1),lubglob(1:nbdim,1))
325    indtruetab(1:nbdim,1,2) = max(indtab(1:nbdim,1,2),lubglob(1:nbdim,1))
326    indtruetab(1:nbdim,2,1) = min(indtab(1:nbdim,2,1),lubglob(1:nbdim,2))
327    indtruetab(1:nbdim,2,2) = min(indtab(1:nbdim,2,2),lubglob(1:nbdim,2))
328!
329    do nb = 1,nbdim
330        if ( do_update(nb) ) then
331            do ndir = 1,2
332                ptres(nb,1,ndir,nb) = indtruetab(nb,ndir,1)
333                ptres(nb,2,ndir,nb) = indtruetab(nb,ndir,2)
334                do i = 1,nbdim
335                    if ( i /= nb ) then
336                        if ( do_update(i) ) then
337                            ptres(i,1,ndir,nb) = indtruetab(i,1,1)
338                            ptres(i,2,ndir,nb) = indtruetab(i,2,2)
339                        else
340                            if (posvar(i) == 1) then
341                                ptres(i,1,ndir,nb) = lb_child(i)
342                                ptres(i,2,ndir,nb) = lb_child(i) + nb_child(i)
343                            else
344                                ptres(i,1,ndir,nb) = lb_child(i)
345                                ptres(i,2,ndir,nb) = lb_child(i) + nb_child(i) - 1
346                            endif
347                        endif
348                    endif
349                enddo
350            enddo
351        endif
352    enddo
353!
354    do nb = 1,nbdim
355        if ( do_update(nb) ) then
356            do ndir = 1,2
357                call Agrif_UpdatenD(type_update, parent, child,             &
358                        ptres(1:nbdim,1,ndir,nb),ptres(1:nbdim,2,ndir,nb),  &
359                        lb_child(1:nbdim),lb_parent(1:nbdim),               &
360                        s_child(1:nbdim),s_parent(1:nbdim),                 &
361                        ds_child(1:nbdim),ds_parent(1:nbdim),               &
362#if defined AGRIF_MPI
363                        posvar,do_update,                                   &
364#endif
365                        nbdim,procname,nb,ndir)
366            enddo
367        endif
368    enddo
369!---------------------------------------------------------------------------------------------------
370end subroutine Agrif_UpdateBcnd
371!===================================================================================================
372!
373!===================================================================================================
374!  subroutine Agrif_UpdatenD
375!
376!> updates a 2D grid variable on the parent grid of the current grid
377!---------------------------------------------------------------------------------------------------
378subroutine Agrif_UpdatenD ( type_update, parent, child,     &
379                            pttab, petab,                   &
380                            lb_child, lb_parent,            &
381                            s_child,  s_parent,             &
382                            ds_child, ds_parent,            &
383#if defined AGRIF_MPI
384                            posvar, do_update,              &
385#endif
386                            nbdim, procname, nb, ndir )
387!---------------------------------------------------------------------------------------------------
388#if defined AGRIF_MPI
389    include 'mpif.h'
390#endif
391!
392    integer, dimension(6),     intent(in)   :: type_update  !< Type of update (copy or average)
393    type(Agrif_Variable), pointer           :: parent       !< Variable of the parent grid
394    type(Agrif_Variable), pointer           :: child        !< Variable of the child grid
395    integer,                   intent(in)   :: nbdim
396    integer, dimension(nbdim), intent(in)   :: pttab        !< Index of the first point inside the domain
397    integer, dimension(nbdim), intent(in)   :: petab        !< Index of the first point inside the domain
398    integer, dimension(nbdim), intent(in)   :: lb_child  !< Index of the first point inside the domain for the child
399                                                            !!    grid variable
400    integer, dimension(nbdim), intent(in)   :: lb_parent !< Index of the first point inside the domain for the parent
401                                                            !!    grid variable
402    real,    dimension(nbdim), intent(in)   :: s_child      !< Positions of the child grid
403    real,    dimension(nbdim), intent(in)   :: s_parent     !< Positions of the parent grid
404    real,    dimension(nbdim), intent(in)   :: ds_child     !< Space steps of the child grid
405    real,    dimension(nbdim), intent(in)   :: ds_parent    !< Space steps of the parent grid
406    procedure()                             :: procname     !< Data recovery procedure
407    integer, optional,         intent(in)   :: nb, ndir
408!---------------------------------------------------------------------------------------------------
409    integer, dimension(nbdim)       :: pttruetab, cetruetab
410#if defined AGRIF_MPI
411    integer, dimension(nbdim)       :: posvar  !< Position of the variable on the cell (1 or 2)
412    logical, dimension(nbdim)       :: do_update
413#endif
414    integer, dimension(nbdim)       :: coords
415    integer, dimension(nbdim)       :: indmin, indmax
416    integer, dimension(nbdim)       :: indminglob, indmaxglob
417    real   , dimension(nbdim)       :: s_Child_temp, s_Parent_temp
418    integer, dimension(nbdim)       :: lowerbound,upperbound
419    integer, dimension(nbdim)       :: pttruetabwhole, cetruetabwhole
420    integer, dimension(nbdim,2,2)   :: childarray
421    integer, dimension(nbdim,2,2)   :: parentarray
422    integer,dimension(nbdim)     :: type_update_temp
423    logical :: memberin, member
424    integer :: nbin, ndirin
425!
426#if defined AGRIF_MPI
427!
428    integer,dimension(nbdim)    :: indminglob2,indmaxglob2
429    logical, dimension(0:Agrif_Nbprocs-1) :: sendtoproc1,recvfromproc1
430    logical, dimension(0:Agrif_Nbprocs-1) :: sendtoproc2,recvfromproc2
431    integer                               :: code, local_proc
432    integer                               :: i,j,k
433    integer, dimension(nbdim,4)           :: tab3
434    integer, dimension(nbdim,4,0:Agrif_Nbprocs-1) :: tab4
435    integer, dimension(nbdim,0:Agrif_Nbprocs-1,8) :: tab4t
436    integer, dimension(nbdim,0:Agrif_Nbprocs-1,8) :: tab5t
437    logical :: find_list_update
438    logical, dimension(0:Agrif_Nbprocs-1) :: memberinall, memberinall2
439    logical, dimension(1) :: memberin1
440!
441#endif
442!
443    type(Agrif_Variable), pointer, save :: tempC => NULL()       ! Temporary child grid variable
444    type(Agrif_Variable), pointer, save :: tempP => NULL()       ! Temporary parent grid variable
445    type(Agrif_Variable), pointer, save :: tempCextend => NULL() ! Temporary child
446    type(Agrif_Variable), pointer, save :: tempPextend => NULL() ! Temporary parent
447    type(Agrif_Variable), pointer :: tempP_indic, tempP_average
448    type(Agrif_Variable), pointer :: tempC_indic
449    logical :: compute_average
450    real :: coeff_multi
451    integer :: nb_dimensions
452!
453!   Get local lower and upper bound of the child variable
454    call Agrif_get_var_bounds_array(child, lowerbound, upperbound, nbdim)
455
456! here pttab and petab corresponds to the (global) indices of the points needed in the update
457! pttruetab and cetruetab contains only indices that are present on the local processor
458!
459    coords = child % root_var % coords
460!
461    call Agrif_Childbounds( nbdim, lowerbound, upperbound, pttab, petab, Agrif_Procrank,    &
462                            coords, pttruetab, cetruetab, memberin )
463    call Agrif_Prtbounds( nbdim, indminglob, indmaxglob, s_Parent_temp, s_Child_temp,       &
464                         s_child, ds_child, s_parent, ds_parent,                            &
465                         pttab, petab, lb_child, lb_parent                                  &
466#if defined AGRIF_MPI
467                       , posvar, type_update, do_update, pttruetabwhole, cetruetabwhole     &
468#endif
469            )
470
471#if defined AGRIF_MPI
472!
473    IF (memberin) THEN
474        call Agrif_GlobalToLocalBounds(childarray,lowerbound,upperbound,    &
475                                       pttruetab,cetruetab, coords,         &
476                                       nbdim, Agrif_Procrank, member)
477    ENDIF
478
479    call Agrif_Prtbounds(nbdim, indmin, indmax,                     &
480                         s_Parent_temp, s_Child_temp,               &
481                         s_child, ds_child, s_parent, ds_parent,    &
482                         pttruetab, cetruetab, lb_child, lb_parent, &
483                         posvar, type_update, do_update,            &
484                         pttruetabwhole, cetruetabwhole)
485!
486#else
487    indmin = indminglob
488    indmax = indmaxglob
489    pttruetabwhole = pttruetab
490    cetruetabwhole = cetruetab
491    childarray(:,1,2) = pttruetab
492    childarray(:,2,2) = cetruetab
493#endif
494
495    IF (.not.present(nb)) THEN
496        nbin=0
497        ndirin=0
498    ELSE
499        nbin = nb
500        ndirin = ndir
501    ENDIF
502
503    IF (memberin) THEN
504!
505        IF ( .not.associated(tempC) )  allocate(tempC)
506!
507        call Agrif_array_allocate(tempC,pttruetab,cetruetab,nbdim)
508        call Agrif_var_set_array_tozero(tempC,nbdim)
509
510        SELECT CASE (nbdim)
511        CASE(1)
512            CALL procname(tempC%array1,                 &
513                      childarray(1,1,2),childarray(1,2,2),.TRUE.,nbin,ndirin)
514        CASE(2)
515            CALL procname(tempC%array2,                 &
516                      childarray(1,1,2),childarray(1,2,2),  &
517                      childarray(2,1,2),childarray(2,2,2),.TRUE.,nbin,ndirin)
518        CASE(3)
519            CALL procname(tempC%array3,                 &
520                      childarray(1,1,2),childarray(1,2,2),  &
521                      childarray(2,1,2),childarray(2,2,2),  &
522                      childarray(3,1,2),childarray(3,2,2),.TRUE.,nbin,ndirin)
523        CASE(4)
524            CALL procname(tempC%array4,                 &
525                      childarray(1,1,2),childarray(1,2,2),  &
526                      childarray(2,1,2),childarray(2,2,2),  &
527                      childarray(3,1,2),childarray(3,2,2),  &
528                      childarray(4,1,2),childarray(4,2,2),.TRUE.,nbin,ndirin)
529        CASE(5)
530            CALL procname(tempC%array5,                 &
531                      childarray(1,1,2),childarray(1,2,2),  &
532                      childarray(2,1,2),childarray(2,2,2),  &
533                      childarray(3,1,2),childarray(3,2,2),  &
534                      childarray(4,1,2),childarray(4,2,2),  &
535                      childarray(5,1,2),childarray(5,2,2),.TRUE.,nbin,ndirin)
536        CASE(6)
537            CALL procname(tempC%array6,                 &
538                      childarray(1,1,2),childarray(1,2,2),  &
539                      childarray(2,1,2),childarray(2,2,2),  &
540                      childarray(3,1,2),childarray(3,2,2),  &
541                      childarray(4,1,2),childarray(4,2,2),  &
542                      childarray(5,1,2),childarray(5,2,2),  &
543                      childarray(6,1,2),childarray(6,2,2),.TRUE.,nbin,ndirin)
544        END SELECT
545!
546    ENDIF
547!
548#if defined AGRIF_MPI
549!
550!     tab2 contains the necessary limits of the parent grid for each processor
551
552    if (Associated(child%list_update)) then
553        call Agrif_Find_list_update(child%list_update,pttab,petab,                      &
554                                    lb_child,lb_parent,nbdim,                         &
555                                    find_list_update,tab4t,tab5t,memberinall,memberinall2,  &
556                                    sendtoproc1,recvfromproc1,sendtoproc2,recvfromproc2)
557    else
558        find_list_update = .FALSE.
559    endif
560
561    if (.not.find_list_update) then
562        tab3(:,1) = pttruetab(:)
563        tab3(:,2) = cetruetab(:)
564        tab3(:,3) = pttruetabwhole(:)
565        tab3(:,4) = cetruetabwhole(:)
566!
567        call MPI_ALLGATHER(tab3,4*nbdim,MPI_INTEGER,tab4,4*nbdim,MPI_INTEGER,Agrif_mpi_comm,code)
568
569        if ( .not.associated(tempCextend) ) allocate(tempCextend)
570        do k=0,Agrif_Nbprocs-1
571            do j=1,4
572                do i=1,nbdim
573                    tab4t(i,k,j) = tab4(i,j,k)
574                enddo
575            enddo
576        enddo
577
578        memberin1(1) = memberin
579        call MPI_ALLGATHER(memberin1,1,MPI_LOGICAL,memberinall,1,MPI_LOGICAL,Agrif_mpi_comm,code)
580
581        call Get_External_Data_first(tab4t(:,:,1),tab4t(:,:,2),tab4t(:,:,3),tab4t(:,:,4),   &
582                                     nbdim, memberinall, coords,                            &
583                                     sendtoproc1,recvfromproc1,                             &
584                                     tab4t(:,:,5),tab4t(:,:,6),tab4t(:,:,7),tab4t(:,:,8),   &
585                                     tab4t(:,:,1),tab4t(:,:,2))
586    endif
587
588    call ExchangeSameLevel(sendtoproc1,recvfromproc1,nbdim,         &
589            tab4t(:,:,3),tab4t(:,:,4),tab4t(:,:,5),tab4t(:,:,6),    &
590            tab4t(:,:,7),tab4t(:,:,8),memberin,tempC,tempCextend)
591
592#else
593    tempCextend => tempC
594#endif
595!
596!     Update of the parent grid (tempP) from the child grid (tempC)
597!
598    IF (memberin) THEN
599!
600        IF ( .not.associated(tempP) ) allocate(tempP)
601!
602        call Agrif_array_allocate(tempP,indmin,indmax,nbdim)
603!
604        if ( nbdim == 1 ) then
605            tempP % array1 = 0.
606            call Agrif_Update_1D_Recursive( type_update(1),  &
607                                            tempP%array1,   &
608                                            tempCextend%array1, &
609                                            indmin(1), indmax(1),   &
610                                            pttruetabwhole(1), cetruetabwhole(1),   &
611                                            s_Child_temp(1), s_Parent_temp(1),      &
612                                            ds_child(1), ds_parent(1) )
613                                           
614            IF (Agrif_UseSpecialValueInUpdate) THEN
615            allocate(tempC_indic)
616            allocate(tempP_indic)
617            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array1),ubound(tempCextend%array1),nbdim)
618            call Agrif_array_allocate(tempP_indic,lbound(tempP%array1),ubound(tempP%array1),nbdim)
619
620            compute_average = .FALSE.
621            type_update_temp(1:nbdim) = type_update(1:nbdim)
622            IF (ANY(type_update(1:nbdim) == Agrif_Update_Full_Weighting)) THEN
623              compute_average = .TRUE.
624              allocate(tempP_average)
625              call Agrif_array_allocate(tempP_average,lbound(tempP%array1),ubound(tempP%array1),nbdim)
626              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
627                type_update_temp(1:nbdim) = Agrif_Update_Average
628              END WHERE
629              call Agrif_Update_1D_Recursive( type_update_temp(1),   &
630                                            tempP_average%array1,       &
631                                            tempCextend%array1, &
632                                            indmin(1), indmax(1),   &
633                                            pttruetabwhole(1), cetruetabwhole(1),   &
634                                            s_Child_temp(1), s_Parent_temp(1),      &
635                                            ds_child(1), ds_parent(1) )
636              coeff_multi = 1.
637              do nb_dimensions=1,nbdim
638                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
639              enddo
640            ENDIF
641           
642            WHERE (tempCextend%array1 == Agrif_SpecialValueFineGrid)
643              tempC_indic%array1 = 0.
644            ELSEWHERE
645              tempC_indic%array1 = 1.
646            END WHERE
647           
648            Agrif_UseSpecialValueInUpdate = .FALSE.
649            Agrif_Update_Weights = .TRUE.
650 
651             call Agrif_Update_1D_Recursive( type_update_temp(1),   &
652                                            tempP_indic%array1,       &
653                                            tempC_indic%array1, &
654                                            indmin(1), indmax(1),   &
655                                            pttruetabwhole(1), cetruetabwhole(1),   &
656                                            s_Child_temp(1), s_Parent_temp(1),      &
657                                            ds_child(1), ds_parent(1) )
658
659           Agrif_UseSpecialValueInUpdate = .TRUE.
660           Agrif_Update_Weights = .FALSE.
661
662           IF (compute_average) THEN
663               WHERE (tempP_indic%array1 == 0.)
664                  tempP%array1 = Agrif_SpecialValueFineGrid
665               ELSEWHERE ((tempP_indic%array1 == coeff_multi).AND.(tempP%array1 /= Agrif_SpecialValueFineGrid))
666                  tempP%array1 = tempP%array1 /tempP_indic%array1
667               ELSEWHERE
668                  tempP%array1 = tempP_average%array1 /tempP_indic%array1
669               END WHERE
670
671           ELSE
672               WHERE (tempP_indic%array1 == 0.)
673                  tempP%array1 = Agrif_SpecialValueFineGrid
674               ELSEWHERE
675                  tempP%array1 = tempP%array1 /tempP_indic%array1
676               END WHERE
677            ENDIF
678           
679            deallocate(tempP_indic%array1)
680            deallocate(tempC_indic%array1)
681            deallocate(tempC_indic)
682            deallocate(tempP_indic)
683            IF (compute_average) THEN
684              deallocate(tempP_average%array1)
685              deallocate(tempP_average)
686            ENDIF
687            ENDIF
688           
689        endif
690        if ( nbdim == 2 ) then
691            call Agrif_Update_2D_Recursive( type_update(1:2),   &
692                                            tempP%array2,       &
693                                            tempCextend%array2, &
694                                            indmin(1:2), indmax(1:2),   &
695                                            pttruetabwhole(1:2), cetruetabwhole(1:2),   &
696                                            s_Child_temp(1:2), s_Parent_temp(1:2),      &
697                                            ds_child(1:2), ds_parent(1:2) )
698
699            IF (Agrif_UseSpecialValueInUpdate) THEN
700            allocate(tempC_indic)
701            allocate(tempP_indic)
702            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array2),ubound(tempCextend%array2),nbdim)
703            call Agrif_array_allocate(tempP_indic,lbound(tempP%array2),ubound(tempP%array2),nbdim)
704 
705            compute_average = .FALSE.
706            type_update_temp(1:nbdim) = type_update(1:nbdim)
707            IF (ANY(type_update == Agrif_Update_Full_Weighting)) THEN
708              compute_average = .TRUE.
709              allocate(tempP_average)
710              call Agrif_array_allocate(tempP_average,lbound(tempP%array2),ubound(tempP%array2),nbdim)
711              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
712                type_update_temp(1:nbdim) = Agrif_Update_Average
713              END WHERE
714              call Agrif_Update_2D_Recursive( type_update_temp(1:2),   &
715                                            tempP_average%array2,       &
716                                            tempCextend%array2, &
717                                            indmin(1:2), indmax(1:2),   &
718                                            pttruetabwhole(1:2), cetruetabwhole(1:2),   &
719                                            s_Child_temp(1:2), s_Parent_temp(1:2),      &
720                                            ds_child(1:2), ds_parent(1:2) )
721              coeff_multi = 1.
722              do nb_dimensions=1,nbdim
723                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
724              enddo
725            ENDIF
726           
727            WHERE (tempCextend%array2 == Agrif_SpecialValueFineGrid)
728              tempC_indic%array2 = 0.
729            ELSEWHERE
730              tempC_indic%array2 = 1.
731            END WHERE
732           
733            Agrif_UseSpecialValueInUpdate = .FALSE.
734            Agrif_Update_Weights = .TRUE.
735           
736            call Agrif_Update_2D_Recursive( type_update_temp(1:2),   &
737                                            tempP_indic%array2,       &
738                                            tempC_indic%array2, &
739                                            indmin(1:2), indmax(1:2),   &
740                                            pttruetabwhole(1:2), cetruetabwhole(1:2),   &
741                                            s_Child_temp(1:2), s_Parent_temp(1:2),      &
742                                            ds_child(1:2), ds_parent(1:2) )
743
744           Agrif_UseSpecialValueInUpdate = .TRUE.
745           Agrif_Update_Weights = .FALSE.
746
747           IF (compute_average) THEN
748               WHERE (tempP_indic%array2 == 0.)
749                  tempP%array2 = Agrif_SpecialValueFineGrid
750               ELSEWHERE ((tempP_indic%array2 == coeff_multi).AND.(tempP%array2 /= Agrif_SpecialValueFineGrid))
751                  tempP%array2 = tempP%array2 /tempP_indic%array2
752               ELSEWHERE
753                  tempP%array2 = tempP_average%array2 /tempP_indic%array2
754               END WHERE
755
756           ELSE
757               WHERE (tempP_indic%array2 == 0.)
758                  tempP%array2 = Agrif_SpecialValueFineGrid
759               ELSEWHERE
760                  tempP%array2 = tempP%array2 /tempP_indic%array2
761               END WHERE
762            ENDIF
763           
764            deallocate(tempP_indic%array2)
765            deallocate(tempC_indic%array2)
766            deallocate(tempC_indic)
767            deallocate(tempP_indic)
768            IF (compute_average) THEN
769              deallocate(tempP_average%array2)
770              deallocate(tempP_average)
771            ENDIF
772            ENDIF
773           
774        endif
775        if ( nbdim == 3 ) then
776            call Agrif_Update_3D_Recursive( type_update(1:3),   &
777                                            tempP%array3,       &
778                                            tempCextend%array3, &
779                                            indmin(1:3), indmax(1:3),   &
780                                            pttruetabwhole(1:3), cetruetabwhole(1:3),   &
781                                            s_Child_temp(1:3), s_Parent_temp(1:3),      &
782                                            ds_child(1:3), ds_parent(1:3) )
783                                           
784            IF (Agrif_UseSpecialValueInUpdate) THEN
785            allocate(tempC_indic)
786            allocate(tempP_indic)
787            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array3),ubound(tempCextend%array3),nbdim)
788            call Agrif_array_allocate(tempP_indic,lbound(tempP%array3),ubound(tempP%array3),nbdim)
789
790            compute_average = .FALSE.
791            type_update_temp(1:nbdim) = type_update(1:nbdim)
792            IF (ANY(type_update == Agrif_Update_Full_Weighting)) THEN
793              compute_average = .TRUE.
794              allocate(tempP_average)
795              call Agrif_array_allocate(tempP_average,lbound(tempP%array3),ubound(tempP%array3),nbdim)
796              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
797                type_update_temp(1:nbdim) = Agrif_Update_Average
798              END WHERE
799              call Agrif_Update_3D_Recursive( type_update_temp(1:3),   &
800                                            tempP_average%array3,       &
801                                            tempCextend%array3, &
802                                            indmin(1:3), indmax(1:3),   &
803                                            pttruetabwhole(1:3), cetruetabwhole(1:3),   &
804                                            s_Child_temp(1:3), s_Parent_temp(1:3),      &
805                                            ds_child(1:3), ds_parent(1:3) )
806              coeff_multi = 1.
807              do nb_dimensions=1,nbdim
808                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
809              enddo
810            ENDIF
811           
812            WHERE (tempCextend%array3 == Agrif_SpecialValueFineGrid)
813              tempC_indic%array3 = 0.
814            ELSEWHERE
815              tempC_indic%array3 = 1.
816            END WHERE
817           
818            Agrif_UseSpecialValueInUpdate = .FALSE.
819            Agrif_Update_Weights = .TRUE.
820 
821             call Agrif_Update_3D_Recursive( type_update_temp(1:3),   &
822                                            tempP_indic%array3,       &
823                                            tempC_indic%array3, &
824                                            indmin(1:3), indmax(1:3),   &
825                                            pttruetabwhole(1:3), cetruetabwhole(1:3),   &
826                                            s_Child_temp(1:3), s_Parent_temp(1:3),      &
827                                            ds_child(1:3), ds_parent(1:3) )
828
829           Agrif_UseSpecialValueInUpdate = .TRUE.
830           Agrif_Update_Weights = .FALSE.
831
832           IF (compute_average) THEN
833               WHERE (tempP_indic%array3 == 0.)
834                  tempP%array3 = Agrif_SpecialValueFineGrid
835               ELSEWHERE ((tempP_indic%array3 == coeff_multi).AND.(tempP%array3 /= Agrif_SpecialValueFineGrid))
836                  tempP%array3 = tempP%array3 /tempP_indic%array3
837               ELSEWHERE
838                  tempP%array3 = tempP_average%array3 /tempP_indic%array3
839               END WHERE
840
841           ELSE
842               WHERE (tempP_indic%array3 == 0.)
843                  tempP%array3 = Agrif_SpecialValueFineGrid
844               ELSEWHERE
845                  tempP%array3 = tempP%array3 /tempP_indic%array3
846               END WHERE
847            ENDIF
848           
849            deallocate(tempP_indic%array3)
850            deallocate(tempC_indic%array3)
851            deallocate(tempC_indic)
852            deallocate(tempP_indic)
853            IF (compute_average) THEN
854              deallocate(tempP_average%array3)
855              deallocate(tempP_average)
856            ENDIF
857            ENDIF
858           
859        endif
860        if ( nbdim == 4 ) then
861            call Agrif_Update_4D_Recursive( type_update(1:4),   &
862                                            tempP%array4,       &
863                                            tempCextend%array4, &
864                                            indmin(1:4), indmax(1:4),   &
865                                            pttruetabwhole(1:4), cetruetabwhole(1:4),   &
866                                            s_Child_temp(1:4), s_Parent_temp(1:4),      &
867                                            ds_child(1:4), ds_parent(1:4) )
868                                           
869            IF (Agrif_UseSpecialValueInUpdate) THEN
870           
871            allocate(tempC_indic)
872            allocate(tempP_indic)
873            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array4),ubound(tempCextend%array4),nbdim)
874            call Agrif_array_allocate(tempP_indic,lbound(tempP%array4),ubound(tempP%array4),nbdim)
875           
876            compute_average = .FALSE.
877            type_update_temp(1:nbdim) = type_update(1:nbdim)
878            IF (ANY(type_update == Agrif_Update_Full_Weighting)) THEN
879              compute_average = .TRUE.
880              allocate(tempP_average)
881              call Agrif_array_allocate(tempP_average,lbound(tempP%array4),ubound(tempP%array4),nbdim)
882              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
883                type_update_temp(1:nbdim) = Agrif_Update_Average
884              END WHERE
885              call Agrif_Update_4D_Recursive( type_update_temp(1:4),   &
886                                            tempP_average%array4,       &
887                                            tempCextend%array4, &
888                                            indmin(1:4), indmax(1:4),   &
889                                            pttruetabwhole(1:4), cetruetabwhole(1:4),   &
890                                            s_Child_temp(1:4), s_Parent_temp(1:4),      &
891                                            ds_child(1:4), ds_parent(1:4) )
892              coeff_multi = 1.
893              do nb_dimensions=1,nbdim
894                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
895              enddo
896            ENDIF
897           
898            WHERE (tempCextend%array4 == Agrif_SpecialValueFineGrid)
899              tempC_indic%array4 = 0.
900            ELSEWHERE
901              tempC_indic%array4 = 1.
902            END WHERE
903           
904            Agrif_UseSpecialValueInUpdate = .FALSE.
905            Agrif_Update_Weights = .TRUE.
906 
907             call Agrif_Update_4D_Recursive( type_update_temp(1:4),   &
908                                            tempP_indic%array4,       &
909                                            tempC_indic%array4, &
910                                            indmin(1:4), indmax(1:4),   &
911                                            pttruetabwhole(1:4), cetruetabwhole(1:4),   &
912                                            s_Child_temp(1:4), s_Parent_temp(1:4),      &
913                                            ds_child(1:4), ds_parent(1:4) )
914
915           Agrif_UseSpecialValueInUpdate = .TRUE.
916           Agrif_Update_Weights = .FALSE.
917           
918           IF (compute_average) THEN
919               WHERE (tempP_indic%array4 == 0.)
920                  tempP%array4 = Agrif_SpecialValueFineGrid
921               ELSEWHERE ((tempP_indic%array4 == coeff_multi).AND.(tempP%array4 /= Agrif_SpecialValueFineGrid))
922                  tempP%array4 = tempP%array4 /tempP_indic%array4
923               ELSEWHERE
924                  tempP%array4 = tempP_average%array4 /tempP_indic%array4
925               END WHERE
926
927           ELSE
928               WHERE (tempP_indic%array4 == 0.)
929                  tempP%array4 = Agrif_SpecialValueFineGrid
930               ELSEWHERE
931                  tempP%array4 = tempP%array4 /tempP_indic%array4
932               END WHERE
933            ENDIF
934            deallocate(tempP_indic%array4)
935            deallocate(tempC_indic%array4)
936            deallocate(tempC_indic)
937            deallocate(tempP_indic)
938            IF (compute_average) THEN
939              deallocate(tempP_average%array4)
940              deallocate(tempP_average)
941            ENDIF
942            ENDIF
943           
944        endif
945        if ( nbdim == 5 ) then
946            call Agrif_Update_5D_Recursive( type_update(1:5),   &
947                                            tempP%array5,       &
948                                            tempCextend%array5, &
949                                            indmin(1:5), indmax(1:5),   &
950                                            pttruetabwhole(1:5), cetruetabwhole(1:5),   &
951                                            s_Child_temp(1:5), s_Parent_temp(1:5),      &
952                                            ds_child(1:5), ds_parent(1:5) )
953                                           
954            IF (Agrif_UseSpecialValueInUpdate) THEN
955            allocate(tempC_indic)
956            allocate(tempP_indic)
957            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array5),ubound(tempCextend%array5),nbdim)
958            call Agrif_array_allocate(tempP_indic,lbound(tempP%array5),ubound(tempP%array5),nbdim)
959
960            compute_average = .FALSE.
961            type_update_temp(1:nbdim) = type_update(1:nbdim)
962            IF (ANY(type_update == Agrif_Update_Full_Weighting)) THEN
963              compute_average = .TRUE.
964              allocate(tempP_average)
965              call Agrif_array_allocate(tempP_average,lbound(tempP%array5),ubound(tempP%array5),nbdim)
966              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
967                type_update_temp(1:nbdim) = Agrif_Update_Average
968              END WHERE
969              call Agrif_Update_5D_Recursive( type_update_temp(1:5),   &
970                                            tempP_average%array5,       &
971                                            tempCextend%array5, &
972                                            indmin(1:5), indmax(1:5),   &
973                                            pttruetabwhole(1:5), cetruetabwhole(1:5),   &
974                                            s_Child_temp(1:5), s_Parent_temp(1:5),      &
975                                            ds_child(1:5), ds_parent(1:5) )
976              coeff_multi = 1.
977              do nb_dimensions=1,nbdim
978                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
979              enddo
980            ENDIF
981           
982            WHERE (tempCextend%array5 == Agrif_SpecialValueFineGrid)
983              tempC_indic%array5 = 0.
984            ELSEWHERE
985              tempC_indic%array5 = 1.
986            END WHERE
987           
988            Agrif_UseSpecialValueInUpdate = .FALSE.
989            Agrif_Update_Weights = .TRUE.
990 
991             call Agrif_Update_5D_Recursive( type_update_temp(1:5),   &
992                                            tempP_indic%array5,       &
993                                            tempC_indic%array5, &
994                                            indmin(1:5), indmax(1:5),   &
995                                            pttruetabwhole(1:5), cetruetabwhole(1:5),   &
996                                            s_Child_temp(1:5), s_Parent_temp(1:5),      &
997                                            ds_child(1:5), ds_parent(1:5) )
998
999           Agrif_UseSpecialValueInUpdate = .TRUE.
1000           Agrif_Update_Weights = .FALSE.
1001
1002           IF (compute_average) THEN
1003               WHERE (tempP_indic%array5 == 0.)
1004                  tempP%array5 = Agrif_SpecialValueFineGrid
1005               ELSEWHERE ((tempP_indic%array5 == coeff_multi).AND.(tempP%array5 /= Agrif_SpecialValueFineGrid))
1006                  tempP%array5 = tempP%array5 /tempP_indic%array5
1007               ELSEWHERE
1008                  tempP%array5 = tempP_average%array5 /tempP_indic%array5
1009               END WHERE
1010
1011           ELSE
1012               WHERE (tempP_indic%array5 == 0.)
1013                  tempP%array5 = Agrif_SpecialValueFineGrid
1014               ELSEWHERE
1015                  tempP%array5 = tempP%array5 /tempP_indic%array5
1016               END WHERE
1017            ENDIF
1018           
1019            deallocate(tempP_indic%array5)
1020            deallocate(tempC_indic%array5)
1021            deallocate(tempC_indic)
1022            deallocate(tempP_indic)
1023            IF (compute_average) THEN
1024              deallocate(tempP_average%array5)
1025              deallocate(tempP_average)
1026            ENDIF
1027            ENDIF
1028           
1029        endif
1030        if ( nbdim == 6 ) then
1031            call Agrif_Update_6D_Recursive( type_update(1:6),   &
1032                                            tempP%array6,       &
1033                                            tempCextend%array6, &
1034                                            indmin(1:6), indmax(1:6),   &
1035                                            pttruetabwhole(1:6), cetruetabwhole(1:6),   &
1036                                            s_Child_temp(1:6), s_Parent_temp(1:6),      &
1037                                            ds_child(1:6), ds_parent(1:6) )
1038            IF (Agrif_UseSpecialValueInUpdate) THEN
1039            allocate(tempC_indic)
1040            allocate(tempP_indic)
1041            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array6),ubound(tempCextend%array6),nbdim)
1042            call Agrif_array_allocate(tempP_indic,lbound(tempP%array6),ubound(tempP%array6),nbdim)
1043
1044            compute_average = .FALSE.
1045            type_update_temp(1:nbdim) = type_update(1:nbdim)
1046            IF (ANY(type_update == Agrif_Update_Full_Weighting)) THEN
1047              compute_average = .TRUE.
1048              allocate(tempP_average)
1049              call Agrif_array_allocate(tempP_average,lbound(tempP%array6),ubound(tempP%array6),nbdim)
1050              type_update_temp(1:nbdim) = type_update
1051              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
1052                type_update_temp(1:nbdim) = Agrif_Update_Average
1053              END WHERE
1054              call Agrif_Update_6D_Recursive( type_update_temp(1:6),   &
1055                                            tempP_average%array6,       &
1056                                            tempCextend%array6, &
1057                                            indmin(1:6), indmax(1:6),   &
1058                                            pttruetabwhole(1:6), cetruetabwhole(1:6),   &
1059                                            s_Child_temp(1:6), s_Parent_temp(1:6),      &
1060                                            ds_child(1:6), ds_parent(1:6) )
1061              coeff_multi = 1.
1062              do nb_dimensions=1,nbdim
1063                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
1064              enddo
1065            ENDIF
1066
1067           IF (compute_average) THEN
1068               WHERE (tempP_indic%array6 == 0.)
1069                  tempP%array6 = Agrif_SpecialValueFineGrid
1070               ELSEWHERE ((tempP_indic%array6 == coeff_multi).AND.(tempP%array6 /= Agrif_SpecialValueFineGrid))
1071                  tempP%array6 = tempP%array6 /tempP_indic%array6
1072               ELSEWHERE
1073                  tempP%array6 = tempP_average%array6 /tempP_indic%array6
1074               END WHERE
1075
1076           ELSE
1077               WHERE (tempP_indic%array6 == 0.)
1078                  tempP%array6 = Agrif_SpecialValueFineGrid
1079               ELSEWHERE
1080                  tempP%array6 = tempP%array6 /tempP_indic%array6
1081               END WHERE
1082            ENDIF
1083           
1084            Agrif_UseSpecialValueInUpdate = .FALSE.
1085            Agrif_Update_Weights = .TRUE.
1086 
1087             call Agrif_Update_6D_Recursive( type_update_temp(1:6),   &
1088                                            tempP_indic%array6,       &
1089                                            tempC_indic%array6, &
1090                                            indmin(1:6), indmax(1:6),   &
1091                                            pttruetabwhole(1:6), cetruetabwhole(1:6),   &
1092                                            s_Child_temp(1:6), s_Parent_temp(1:6),      &
1093                                            ds_child(1:6), ds_parent(1:6) )
1094
1095           Agrif_UseSpecialValueInUpdate = .TRUE.
1096           Agrif_Update_Weights = .FALSE.
1097           
1098            WHERE (tempP_indic%array6 == 0.)
1099              tempP%array6 = Agrif_SpecialValueFineGrid
1100            ELSEWHERE
1101              tempP%array6 = tempP%array6 /tempP_indic%array6
1102            END WHERE
1103           
1104            deallocate(tempP_indic%array6)
1105            deallocate(tempC_indic%array6)
1106            deallocate(tempC_indic)
1107            deallocate(tempP_indic)
1108            IF (compute_average) THEN
1109              deallocate(tempP_average%array6)
1110              deallocate(tempP_average)
1111            ENDIF
1112            ENDIF
1113        endif
1114!
1115        call Agrif_array_deallocate(tempCextend,nbdim)
1116!
1117    ENDIF
1118
1119#if defined AGRIF_MPI
1120    local_proc = Agrif_Procrank
1121    call Agrif_get_var_bounds_array(parent,lowerbound,upperbound,nbdim)
1122    call Agrif_ChildGrid_to_ParentGrid()
1123    call Agrif_Childbounds(nbdim, lowerbound, upperbound,                   &
1124                           indminglob,  indmaxglob,  local_proc, coords,    &
1125                           indminglob2, indmaxglob2, member)
1126!
1127    IF (member) THEN
1128        call Agrif_GlobalToLocalBounds(parentarray, lowerbound, upperbound, &
1129                                       indminglob2, indmaxglob2, coords,    &
1130                                       nbdim, local_proc, member)
1131    ENDIF
1132
1133    call Agrif_ParentGrid_to_ChildGrid()
1134
1135    if (.not.find_list_update) then
1136        tab3(:,1) = indmin(:)
1137        tab3(:,2) = indmax(:)
1138        tab3(:,3) = indminglob2(:)
1139        tab3(:,4) = indmaxglob2(:)
1140!
1141        call MPI_ALLGATHER(tab3,4*nbdim,MPI_INTEGER,tab4,4*nbdim,MPI_INTEGER,Agrif_mpi_comm,code)
1142
1143        IF ( .not.associated(tempPextend) ) allocate(tempPextend)
1144        DO k=0,Agrif_Nbprocs-1
1145            do j=1,4
1146                do i=1,nbdim
1147                    tab5t(i,k,j) = tab4(i,j,k)
1148                enddo
1149            enddo
1150        enddo
1151
1152        memberin1(1) = member
1153        call MPI_ALLGATHER(memberin1,1,MPI_LOGICAL,memberinall2,1,MPI_LOGICAL,Agrif_mpi_comm,code)
1154        call Get_External_Data_first(tab5t(:,:,1),tab5t(:,:,2),tab5t(:,:,3),tab5t(:,:,4),   &
1155                                     nbdim, memberinall2, coords,                           &
1156                                     sendtoproc2, recvfromproc2,                            &
1157                                     tab5t(:,:,5),tab5t(:,:,6),tab5t(:,:,7),tab5t(:,:,8),   &
1158                                     tab5t(:,:,1),tab5t(:,:,2))
1159
1160        call Agrif_Addto_list_update(child%list_update,pttab,petab,lb_child,lb_parent,      &
1161                                     nbdim,tab4t,tab5t,memberinall,memberinall2,            &
1162                                     sendtoproc1,recvfromproc1,sendtoproc2,recvfromproc2)
1163
1164    endif
1165
1166    call ExchangeSameLevel(sendtoproc2,recvfromproc2,nbdim,                     &
1167                            tab5t(:,:,3),tab5t(:,:,4),tab5t(:,:,5),tab5t(:,:,6),&
1168                            tab5t(:,:,7),tab5t(:,:,8),member,tempP,tempPextend)
1169#else
1170    tempPextend => tempP
1171    parentarray(:,1,1) = indmin
1172    parentarray(:,2,1) = indmax
1173    parentarray(:,1,2) = indmin
1174    parentarray(:,2,2) = indmax
1175    member = .TRUE.
1176#endif
1177!
1178!   Special values on the child grid
1179    if ( Agrif_UseSpecialValueFineGrid ) then
1180!
1181!cc         noraftab(1:nbdim) =
1182!cc     &    child % root_var % interptab(1:nbdim) == 'N'
1183!
1184#if defined AGRIF_MPI
1185!
1186!          allocate(childvalues% var)
1187!
1188!          Call Agrif_array_allocate(childvalues%var,
1189!     &                      pttruetab,cetruetab,nbdim)
1190!          Call Agrif_var_full_copy_array(childvalues% var,
1191!     &                                tempC,
1192!     &                                nbdim)
1193!          Call Agrif_CheckMasknD(tempC,childvalues,
1194!     &                           pttruetab(1:nbdim),cetruetab(1:nbdim),
1195!     &                           pttruetab(1:nbdim),cetruetab(1:nbdim),
1196!     &                           noraftab(1:nbdim),nbdim)
1197!          Call Agrif_array_deallocate(childvalues% var,nbdim)
1198!         Deallocate(childvalues % var)
1199!
1200#else
1201!
1202!          Call Agrif_get_var_bounds_array(child,
1203!     &                              lowerbound,upperbound,nbdim)
1204!          Call Agrif_CheckMasknD(tempC,child,
1205!     &                           pttruetab(1:nbdim),cetruetab(1:nbdim),
1206!     &                           lowerbound,
1207!     &                           upperbound,
1208!     &                           noraftab(1:nbdim),nbdim)
1209!
1210#endif
1211!
1212    endif
1213!
1214!   Special values on the parent grid
1215    if (Agrif_UseSpecialValue) then
1216!
1217#if defined AGRIF_MPI
1218!
1219!          Call GiveAgrif_SpecialValueToTab_mpi(parent,tempP,
1220!     &                 parentarray,
1221!     &                 Agrif_SpecialValue,nbdim)
1222!
1223!
1224#else
1225!
1226!          Call GiveAgrif_SpecialValueToTab(parent,tempP,
1227!     &                  indmin,indmax,
1228!     &                  Agrif_SpecialValue,nbdim)
1229!
1230#endif
1231!
1232    endif
1233!
1234    IF (member) THEN
1235
1236        call Agrif_ChildGrid_to_ParentGrid()
1237!
1238        SELECT CASE(nbdim)
1239        CASE(1)
1240            call procname( tempPextend % array1(            &
1241                    parentarray(1,1,1):parentarray(1,2,1)), &
1242                    parentarray(1,1,2),parentarray(1,2,2),.FALSE.,nbin,ndirin)
1243        CASE(2)
1244            call procname( tempPextend % array2(            &
1245                    parentarray(1,1,1):parentarray(1,2,1),  &
1246                    parentarray(2,1,1):parentarray(2,2,1)), &
1247                    parentarray(1,1,2),parentarray(1,2,2),  &
1248                    parentarray(2,1,2),parentarray(2,2,2),.FALSE.,nbin,ndirin)
1249        CASE(3)
1250            call procname( tempPextend % array3(            &
1251                    parentarray(1,1,1):parentarray(1,2,1),  &
1252                    parentarray(2,1,1):parentarray(2,2,1),  &
1253                    parentarray(3,1,1):parentarray(3,2,1)), &
1254                    parentarray(1,1,2),parentarray(1,2,2),  &
1255                    parentarray(2,1,2),parentarray(2,2,2),  &
1256                    parentarray(3,1,2),parentarray(3,2,2),.FALSE.,nbin,ndirin)
1257        CASE(4)
1258            call procname( tempPextend % array4(            &
1259                    parentarray(1,1,1):parentarray(1,2,1),  &
1260                    parentarray(2,1,1):parentarray(2,2,1),  &
1261                    parentarray(3,1,1):parentarray(3,2,1),  &
1262                    parentarray(4,1,1):parentarray(4,2,1)), &
1263                    parentarray(1,1,2),parentarray(1,2,2),  &
1264                    parentarray(2,1,2),parentarray(2,2,2),  &
1265                    parentarray(3,1,2),parentarray(3,2,2),  &
1266                    parentarray(4,1,2),parentarray(4,2,2),.FALSE.,nbin,ndirin)
1267        CASE(5)
1268            call procname( tempPextend % array5(            &
1269                    parentarray(1,1,1):parentarray(1,2,1),  &
1270                    parentarray(2,1,1):parentarray(2,2,1),  &
1271                    parentarray(3,1,1):parentarray(3,2,1),  &
1272                    parentarray(4,1,1):parentarray(4,2,1),  &
1273                    parentarray(5,1,1):parentarray(5,2,1)), &
1274                    parentarray(1,1,2),parentarray(1,2,2),  &
1275                    parentarray(2,1,2),parentarray(2,2,2),  &
1276                    parentarray(3,1,2),parentarray(3,2,2),  &
1277                    parentarray(4,1,2),parentarray(4,2,2),  &
1278                    parentarray(5,1,2),parentarray(5,2,2),.FALSE.,nbin,ndirin)
1279        CASE(6)
1280            call procname( tempPextend % array6(            &
1281                    parentarray(1,1,1):parentarray(1,2,1),  &
1282                    parentarray(2,1,1):parentarray(2,2,1),  &
1283                    parentarray(3,1,1):parentarray(3,2,1),  &
1284                    parentarray(4,1,1):parentarray(4,2,1),  &
1285                    parentarray(5,1,1):parentarray(5,2,1),  &
1286                    parentarray(6,1,1):parentarray(6,2,1)), &
1287                    parentarray(1,1,2),parentarray(1,2,2),  &
1288                    parentarray(2,1,2),parentarray(2,2,2),  &
1289                    parentarray(3,1,2),parentarray(3,2,2),  &
1290                    parentarray(4,1,2),parentarray(4,2,2),  &
1291                    parentarray(5,1,2),parentarray(5,2,2),  &
1292                    parentarray(6,1,2),parentarray(6,2,2),.FALSE.,nbin,ndirin)
1293        END SELECT
1294!
1295        call Agrif_ParentGrid_to_ChildGrid()
1296!
1297        call Agrif_array_deallocate(tempPextend,nbdim)
1298!
1299    ENDIF
1300!
1301#if defined AGRIF_MPI
1302    IF (memberin) THEN
1303        call Agrif_array_deallocate(tempP,nbdim)
1304        call Agrif_array_deallocate(tempC,nbdim)
1305    ENDIF
1306#endif
1307!---------------------------------------------------------------------------------------------------
1308end subroutine Agrif_UpdatenD
1309!===================================================================================================
1310!
1311!===================================================================================================
1312!  subroutine Agrif_Prtbounds
1313!
1314!> calculates the bounds of the parent grid to be updated by the child grid
1315!---------------------------------------------------------------------------------------------------
1316subroutine Agrif_Prtbounds ( nbdim, indmin, indmax, s_Parent_temp, s_Child_temp,        &
1317                             s_child, ds_child, s_parent, ds_parent,                    &
1318                             pttruetab, cetruetab, lb_child, lb_parent            &
1319#if defined AGRIF_MPI
1320                            ,posvar, type_update, do_update,                &
1321                             pttruetabwhole, cetruetabwhole                             &
1322#endif
1323                )
1324!---------------------------------------------------------------------------------------------------
1325    integer,                   intent(in)   :: nbdim
1326    integer, dimension(nbdim), intent(out)  :: indmin, indmax
1327    real,    dimension(nbdim), intent(out)  :: s_Parent_temp, s_Child_temp
1328    real,    dimension(nbdim), intent(in)   :: s_child,  ds_child
1329    real,    dimension(nbdim), intent(in)   :: s_parent, ds_parent
1330    integer, dimension(nbdim), intent(in)   :: pttruetab, cetruetab
1331    integer, dimension(nbdim), intent(in)   :: lb_child, lb_parent
1332#if defined AGRIF_MPI
1333    integer, dimension(nbdim), intent(in)   :: posvar   !< Position of the variable on the cell (1 or 2)
1334    integer, dimension(nbdim), intent(in)   :: type_update
1335    logical, dimension(nbdim), intent(in)   :: do_update
1336    integer,dimension(nbdim), intent(out)   :: pttruetabwhole, cetruetabwhole
1337#endif
1338!
1339    real,dimension(nbdim) :: dim_newmin,dim_newmax
1340    integer :: i
1341#if defined AGRIF_MPI
1342    real    :: positionmin, positionmax
1343    integer :: imin, imax
1344    integer :: coeffraf
1345#endif
1346!
1347    do i = 1,nbdim
1348!
1349        dim_newmin(i) = s_child(i) + (pttruetab(i) - lb_child(i)) * ds_child(i)
1350        dim_newmax(i) = s_child(i) + (cetruetab(i) - lb_child(i)) * ds_child(i)
1351!
1352        indmin(i) = lb_parent(i) + agrif_ceiling((dim_newmin(i)-s_parent(i))/ds_parent(i))
1353        indmax(i) = lb_parent(i) + agrif_int(    (dim_newmax(i)-s_parent(i))/ds_parent(i))
1354!
1355#if defined AGRIF_MPI
1356        positionmin = s_parent(i) + (indmin(i)-lb_parent(i))*ds_parent(i)
1357        IF ( do_update(i) ) THEN
1358            IF (posvar(i) == 1) THEN
1359                IF      (type_update(i) == Agrif_Update_Average) THEN
1360                    positionmin = positionmin - ds_parent(i)/2.
1361                ELSE IF (type_update(i) == Agrif_Update_Full_Weighting) THEN
1362                    positionmin = positionmin - (ds_parent(i)-ds_child(i))
1363                ENDIF
1364            ELSE
1365                IF (type_update(i) /= Agrif_Update_Full_Weighting) THEN
1366                    positionmin = positionmin - ds_parent(i)/2.
1367                ELSE
1368                    coeffraf = nint(ds_parent(i)/ds_child(i))
1369                    if (mod(coeffraf,2) == 1) then
1370                        positionmin = positionmin - (ds_parent(i)-ds_child(i))
1371                    else
1372                        positionmin = positionmin - (ds_parent(i)-ds_child(i))-ds_child(i)/2.
1373                    endif
1374                ENDIF
1375            ENDIF
1376        ENDIF
1377!
1378        imin = lb_child(i) + agrif_ceiling((positionmin-s_child(i))/ds_child(i))
1379        positionmin = s_child(i)  + (imin - lb_child(i)) * ds_child(i)
1380        positionmax = s_parent(i) + (indmax(i)-lb_parent(i))*ds_parent(i)
1381        pttruetabwhole(i) = imin
1382
1383        IF ( do_update(i) ) THEN
1384            IF (posvar(i) == 1) THEN
1385                IF      (type_update(i) == Agrif_Update_Average) THEN
1386                    positionmax = positionmax  + ds_parent(i)/2.
1387                ELSE IF (type_update(i) == Agrif_Update_Full_Weighting) THEN
1388                    positionmax = positionmax  + (ds_parent(i)-ds_child(i))
1389                ENDIF
1390            ELSE
1391                IF (type_update(i) /= Agrif_Update_Full_Weighting) THEN
1392                    positionmax = positionmax  + ds_parent(i)/2.
1393                ELSE
1394                    coeffraf = nint(ds_parent(i)/ds_child(i))
1395                    if (mod(coeffraf,2) == 1) then
1396                        positionmax = positionmax + (ds_parent(i)-ds_child(i))
1397                    else
1398                        positionmax = positionmax + (ds_parent(i)-ds_child(i)) + ds_child(i)/2.
1399                    endif
1400                ENDIF
1401            ENDIF
1402        ENDIF
1403
1404        imax = lb_child(i) +agrif_int((positionmax-s_child(i))/ds_child(i))
1405        positionmax = s_child(i) + (imax - lb_child(i)) * ds_child(i)
1406        cetruetabwhole(i) = imax
1407#endif
1408!
1409        s_Parent_temp(i) = s_parent(i) + (indmin(i) - lb_parent(i)) * ds_parent(i)
1410        s_Child_temp(i)  = dim_newmin(i)
1411
1412#if defined AGRIF_MPI
1413        s_Child_temp(i) = positionmin
1414#endif
1415!
1416    enddo
1417!---------------------------------------------------------------------------------------------------
1418end subroutine Agrif_Prtbounds
1419!===================================================================================================
1420!
1421!===================================================================================================
1422!  subroutine Agrif_Update_1D_Recursive
1423!
1424!> Updates a 1D grid variable on the parent grid
1425!---------------------------------------------------------------------------------------------------
1426subroutine Agrif_Update_1D_Recursive ( type_update,         &
1427                                       tempP, tempC,        &
1428                                       indmin, indmax,      &
1429                                       lb_child, ub_child,  &
1430                                       s_child,  s_parent,  &
1431                                       ds_child, ds_parent )
1432!---------------------------------------------------------------------------------------------------
1433    integer,                            intent(in)  :: type_update            !< Type of update (copy or average)
1434    integer,                            intent(in)  :: indmin, indmax
1435    integer,                            intent(in)  :: lb_child, ub_child
1436    real,                               intent(in)  ::  s_child,  s_parent
1437    real,                               intent(in)  :: ds_child, ds_parent
1438    real, dimension(indmin:indmax),     intent(out) :: tempP
1439    real, dimension(lb_child:ub_child), intent(in)  :: tempC
1440!---------------------------------------------------------------------------------------------------
1441    call Agrif_UpdateBase(type_update,              &
1442                          tempP(indmin:indmax),     &
1443                          tempC(lb_child:ub_child), &
1444                          indmin, indmax,           &
1445                          lb_child, ub_child,       &
1446                          s_parent, s_child,        &
1447                          ds_parent, ds_child)
1448!---------------------------------------------------------------------------------------------------
1449end subroutine Agrif_Update_1D_Recursive
1450!===================================================================================================
1451!
1452!===================================================================================================
1453!  subroutine Agrif_Update_2D_Recursive
1454!
1455!> updates a 2D grid variable on the parent grid.
1456!! Calls #Agrif_Update_1D_Recursive and #Agrif_UpdateBase
1457!---------------------------------------------------------------------------------------------------
1458subroutine Agrif_Update_2D_Recursive ( type_update,         &
1459                                       tempP, tempC,        &
1460                                       indmin, indmax,      &
1461                                       lb_child, ub_child,  &
1462                                        s_child,  s_parent, &
1463                                       ds_child, ds_parent )
1464!---------------------------------------------------------------------------------------------------
1465    integer, dimension(2),          intent(in)  :: type_update            !< Type of update (copy or average)
1466    integer, dimension(2),          intent(in)  :: indmin, indmax
1467    integer, dimension(2),          intent(in)  :: lb_child, ub_child
1468    real,    dimension(2),          intent(in)  ::  s_child,  s_parent
1469    real,    dimension(2),          intent(in)  :: ds_child, ds_parent
1470    real,    dimension(          &
1471        indmin(1):indmax(1),     &
1472        indmin(2):indmax(2)),       intent(out) :: tempP
1473    real,    dimension(:,:),        intent(in)  :: tempC
1474!---------------------------------------------------------------------------------------------------
1475    real, dimension(indmin(1):indmax(1), lb_child(2):ub_child(2)) :: tabtemp
1476    real, dimension(indmin(2):indmax(2), indmin(1):indmax(1))     :: tempP_trsp
1477    real, dimension(lb_child(2):ub_child(2), indmin(1):indmax(1)) :: tabtemp_trsp
1478    integer :: i, j
1479    integer :: coeffraf
1480!
1481    tabtemp = 0.
1482    coeffraf = nint ( ds_parent(1) / ds_child(1) )
1483!
1484    IF((type_update(1) == Agrif_Update_Average) .AND. (coeffraf /= 1 )) THEN
1485!---CDIR NEXPAND
1486        if ( .NOT. precomputedone(1) ) then
1487            call Average1dPrecompute( ub_child(2)-lb_child(2)+1,    &
1488                                      indmax(1)-indmin(1)+1,              &
1489                                      ub_child(1)-lb_child(1)+1,    &
1490                                      s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1491!            precomputedone(1) = .TRUE.
1492        endif
1493!---CDIR NEXPAND
1494        call Average1dAfterCompute( tabtemp, tempC, size(tabtemp), size(tempC), &
1495                    s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1496!
1497    ELSE IF ((type_update(1) == Agrif_Update_Copy) .AND. (coeffraf /= 1 ))THEN
1498!---CDIR NEXPAND
1499         if ( .NOT. precomputedone(1) ) then
1500            call Agrif_basicupdate_copy1d_before( ub_child(2)-lb_child(2)+1, &
1501                                   indmax(1)-indmin(1)+1,           &
1502                                   ub_child(1)-lb_child(1)+1, &
1503                                   s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1504!            precomputedone(1) = .TRUE.
1505        endif
1506!---CDIR NEXPAND
1507        call Agrif_basicupdate_copy1d_after(tabtemp,tempC,size(tabtemp),size(tempC),1)
1508!
1509    ELSE
1510        do j = lb_child(2),ub_child(2)
1511!
1512!---CDIR NEXPAND
1513            call Agrif_Update_1D_Recursive( type_update(1),             &
1514                                            tabtemp(:,j),               &
1515                                            tempC(:,j-lb_child(2)+1),   &
1516                                            indmin(1), indmax(1),       &
1517                                            lb_child(1),ub_child(1),    &
1518                                            s_child(1),  s_parent(1),   &
1519                                            ds_child(1),ds_parent(1))
1520        enddo
1521    ENDIF
1522!
1523    tabtemp_trsp = TRANSPOSE(tabtemp)
1524    coeffraf = nint(ds_parent(2)/ds_child(2))
1525!
1526    tempP_trsp = 0.
1527!
1528    IF((type_update(2) == Agrif_Update_Average) .AND. (coeffraf /= 1 )) THEN
1529!---CDIR NEXPAND
1530        if ( .NOT. precomputedone(2) ) then
1531            call Average1dPrecompute( indmax(1)-indmin(1)+1,          &
1532                                      indmax(2)-indmin(2)+1,          &
1533                                      ub_child(2)-lb_child(2)+1,&
1534                                      s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1535!            precomputedone(2) = .TRUE.
1536        endif
1537!---CDIR NEXPAND
1538        call Average1dAfterCompute( tempP_trsp, tabtemp_trsp, size(tempP_trsp), size(tabtemp_trsp),&
1539                s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1540!
1541    ELSE IF ((type_update(2) == Agrif_Update_Copy) .AND. (coeffraf /= 1 )) THEN
1542!---CDIR NEXPAND
1543        if ( .NOT. precomputedone(2) ) then
1544            call Agrif_basicupdate_copy1d_before( indmax(1)-indmin(1)+1,             &
1545                                   indmax(2)-indmin(2)+1,             &
1546                                   ub_child(2)-lb_child(2)+1,   &
1547                                   s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1548!            precomputedone(2) = .TRUE.
1549        endif
1550!---CDIR NEXPAND
1551        call Agrif_basicupdate_copy1d_after( tempP_trsp, tabtemp_trsp, size(tempP_trsp), size(tabtemp_trsp),2)
1552!
1553    ELSE
1554        do i = indmin(1),indmax(1)
1555!
1556!---CDIR NEXPAND
1557            call Agrif_UpdateBase(type_update(2),                                       &
1558                                  tempP_trsp(indmin(2):indmax(2),i),            &
1559                                  tabtemp_trsp(lb_child(2):ub_child(2),i),&
1560                                  indmin(2),indmax(2),                          &
1561                                  lb_child(2),ub_child(2),                &
1562                                  s_parent(2),s_child(2),                       &
1563                                  ds_parent(2),ds_child(2))
1564!
1565        enddo
1566    ENDIF
1567!
1568    tempP = TRANSPOSE(tempP_trsp)
1569!---------------------------------------------------------------------------------------------------
1570end subroutine Agrif_Update_2D_Recursive
1571!===================================================================================================
1572!
1573subroutine Agrif_Update_2D_Recursive_ok ( type_update, &
1574                                        tempP, tempC, &
1575                                        indmin, indmax,   &
1576                                       lb_child, ub_child,                    &
1577                                       s_child, s_parent, ds_child, ds_parent )
1578!---------------------------------------------------------------------------------------------------
1579    INTEGER, DIMENSION(2), intent(in)   :: type_update            !< Type of update (copy or average)
1580    INTEGER, DIMENSION(2), intent(in)   :: indmin, indmax
1581    INTEGER, DIMENSION(2), intent(in)   :: lb_child, ub_child
1582    REAL,    DIMENSION(2), intent(in)   :: s_child,  s_parent
1583    REAL,    DIMENSION(2), intent(in)   :: ds_child, ds_parent
1584    REAL,    DIMENSION(                 &
1585                indmin(1):indmax(1),    &
1586                indmin(2):indmax(2)),           intent(out) :: tempP
1587    REAL, DIMENSION(                            &
1588                lb_child(1):ub_child(1),  &
1589                lb_child(2):ub_child(2)), intent(in)  :: tempC
1590!
1591    REAL, DIMENSION(indmin(1):indmax(1), lb_child(2):ub_child(2)) :: tabtemp
1592    INTEGER :: i
1593!
1594    do i = lb_child(2),ub_child(2)
1595        call Agrif_Update_1D_Recursive(type_update(1),                              &
1596                                       tabtemp(:, i),          &
1597                                       tempC(:,i),  &
1598                                       indmin(1),indmax(1),                 &
1599                                       lb_child(1),ub_child(1),       &
1600                                       s_child(1), s_parent(1),             &
1601                                      ds_child(1),ds_parent(1))
1602    enddo
1603!
1604    tempP = 0.
1605!
1606    do i = indmin(1),indmax(1)
1607        call Agrif_UpdateBase(type_update(2),                                       &
1608                              tempP(i,:),             &
1609                              tabtemp(i,:), &
1610                              indmin(2),indmax(2),                          &
1611                              lb_child(2),ub_child(2),                &
1612                              s_parent(2),s_child(2),                       &
1613                             ds_parent(2),ds_child(2))
1614    enddo
1615!---------------------------------------------------------------------------------------------------
1616end subroutine Agrif_Update_2D_Recursive_ok
1617!===================================================================================================
1618
1619!
1620!===================================================================================================
1621!  subroutine Agrif_Update_3D_Recursive
1622!
1623!> Updates a 3D grid variable on the parent grid.
1624!! Calls #Agrif_Update_2D_Recursive and #Agrif_UpdateBase.
1625!---------------------------------------------------------------------------------------------------
1626subroutine Agrif_Update_3D_Recursive ( type_update,         &
1627                                       tempP, tempC,        &
1628                                       indmin, indmax,      &
1629                                       lb_child, ub_child,  &
1630                                        s_child,  s_parent, &
1631                                       ds_child, ds_parent )
1632!---------------------------------------------------------------------------------------------------
1633    integer, dimension(3),          intent(in)  :: type_update            !< Type of update (copy or average)
1634    integer, dimension(3),          intent(in)  :: indmin, indmax
1635    integer, dimension(3),          intent(in)  :: lb_child, ub_child
1636    real,    dimension(3),          intent(in)  ::  s_child,  s_parent
1637    real,    dimension(3),          intent(in)  :: ds_child, ds_parent
1638    real,    dimension(          &
1639        indmin(1):indmax(1),     &
1640        indmin(2):indmax(2),     &
1641        indmin(3):indmax(3)),       intent(out) :: tempP
1642    real, dimension(             &
1643        lb_child(1):ub_child(1), &
1644        lb_child(2):ub_child(2), &
1645        lb_child(3):ub_child(3)),   intent(in)  :: tempC
1646!---------------------------------------------------------------------------------------------------
1647    real, dimension(            &
1648        indmin(1):indmax(1),    &
1649        indmin(2):indmax(2),    &
1650        lb_child(3):ub_child(3))    :: tabtemp
1651    integer :: i,j,k
1652    integer :: coeffraf,locind_child_left
1653    integer :: kuinf
1654!
1655    coeffraf = nint ( ds_parent(1) / ds_child(1) )
1656!
1657    if ((type_update(1) == Agrif_Update_Average) .AND. (coeffraf /= 1 )) then
1658!---CDIR NEXPAND
1659        call Average1dPrecompute(ub_child(2)-lb_child(2)+1,&
1660                                 indmax(1)-indmin(1)+1,&
1661                                 ub_child(1)-lb_child(1)+1,&
1662                                 s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1663        precomputedone(1) = .TRUE.
1664    else if ((type_update(1) == Agrif_Update_Copy) .AND. (coeffraf /= 1 )) then
1665!---CDIR NEXPAND
1666         call Agrif_basicupdate_copy1d_before(ub_child(2)-lb_child(2)+1, &
1667                               indmax(1)-indmin(1)+1,           &
1668                               ub_child(1)-lb_child(1)+1, &
1669                               s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1670        precomputedone(1) = .TRUE.
1671    endif
1672!
1673    coeffraf = nint ( ds_parent(2) / ds_child(2) )
1674!
1675    if ((type_update(2) == Agrif_Update_Average) .AND. (coeffraf /= 1 )) then
1676!---CDIR NEXPAND
1677        call Average1dPrecompute(indmax(1)-indmin(1)+1,&
1678                                 indmax(2)-indmin(2)+1,&
1679                                 ub_child(2)-lb_child(2)+1,&
1680                                 s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1681        precomputedone(2) = .TRUE.
1682    else if ((type_update(2) == Agrif_Update_Copy) .AND. (coeffraf /= 1 )) then
1683!---CDIR NEXPAND
1684        call Agrif_basicupdate_copy1d_before( indmax(1)-indmin(1)+1,           &
1685                               indmax(2)-indmin(2)+1,           &
1686                               ub_child(2)-lb_child(2)+1, &
1687                               s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1688        precomputedone(2) = .TRUE.
1689    endif
1690!
1691    do k = lb_child(3),ub_child(3)
1692        call Agrif_Update_2D_Recursive( type_update(1:2),tabtemp(:,:,k),tempC(:,:,k), &
1693                                        indmin(1:2),indmax(1:2),                &
1694                                        lb_child(1:2),ub_child(1:2),      &
1695                                        s_child(1:2),s_parent(1:2),             &
1696                                        ds_child(1:2),ds_parent(1:2) )
1697    enddo
1698!
1699    precomputedone(1) = .FALSE.
1700    precomputedone(2) = .FALSE.
1701!
1702    coeffraf = nint ( ds_parent(3) / ds_child(3) )
1703    locind_child_left = 1 + agrif_int((s_parent(3)-s_child(3))/ds_child(3))
1704!
1705    if (coeffraf == 1) then
1706        kuinf = lb_child(3)+locind_child_left-2
1707        do k=indmin(3),indmax(3)
1708            kuinf = kuinf + 1
1709            do j = indmin(2),indmax(2)
1710            do i = indmin(1),indmax(1)
1711                tempP(i,j,k) = tabtemp(i,j,kuinf)
1712            enddo
1713            enddo
1714        enddo
1715    else
1716        tempP = 0.
1717        do j = indmin(2),indmax(2)
1718        do i = indmin(1),indmax(1)
1719            call Agrif_UpdateBase(type_update(3),tempP(i,j,:),tabtemp(i,j,:),   &
1720                                  indmin(3),indmax(3),                  &
1721                                  lb_child(3),ub_child(3),        &
1722                                  s_parent(3),s_child(3),               &
1723                                  ds_parent(3),ds_child(3))
1724!
1725        enddo
1726        enddo
1727    endif
1728!---------------------------------------------------------------------------------------------------
1729end subroutine Agrif_Update_3D_Recursive
1730!===================================================================================================
1731!
1732!===================================================================================================
1733!  subroutine Agrif_Update_4D_Recursive
1734!
1735!> Updates a 4D grid variable on the parent grid.
1736!! Calls #Agrif_Update_3D_Recursive and #Agrif_UpdateBase.
1737!---------------------------------------------------------------------------------------------------
1738subroutine Agrif_Update_4D_Recursive ( type_update,         &
1739                                       tempP, tempC,        &
1740                                       indmin, indmax,      &
1741                                       lb_child, ub_child,  &
1742                                        s_child,  s_parent, &
1743                                       ds_child, ds_parent )
1744!---------------------------------------------------------------------------------------------------
1745    integer, dimension(4),          intent(in)  :: type_update            !< Type of update (copy or average)
1746    integer, dimension(4),          intent(in)  :: indmin, indmax
1747    integer, dimension(4),          intent(in)  :: lb_child, ub_child
1748    real,    dimension(4),          intent(in)  ::  s_child,  s_parent
1749    real,    dimension(4),          intent(in)  :: ds_child, ds_parent
1750    real,    dimension(          &
1751        indmin(1):indmax(1),     &
1752        indmin(2):indmax(2),     &
1753        indmin(3):indmax(3),     &
1754        indmin(4):indmax(4)),       intent(out) :: tempP
1755    real, dimension(             &
1756        lb_child(1):ub_child(1), &
1757        lb_child(2):ub_child(2), &
1758        lb_child(3):ub_child(3), &
1759        lb_child(4):ub_child(4)),   intent(in)  :: tempC
1760!---------------------------------------------------------------------------------------------------
1761    real, dimension(:,:,:,:), allocatable       :: tabtemp
1762    integer :: i,j,k,l
1763!
1764    allocate(tabtemp(indmin(1):indmax(1), &
1765                     indmin(2):indmax(2), &
1766                     indmin(3):indmax(3), &
1767                     lb_child(4):ub_child(4)))
1768!
1769    do l = lb_child(4), ub_child(4)
1770        call Agrif_Update_3D_Recursive(type_update(1:3),                    &
1771                                       tabtemp(indmin(1):indmax(1),         &
1772                                               indmin(2):indmax(2),         &
1773                                               indmin(3):indmax(3), l),     &
1774                                       tempC(lb_child(1):ub_child(1),       &
1775                                             lb_child(2):ub_child(2),       &
1776                                             lb_child(3):ub_child(3), l),   &
1777                                       indmin(1:3), indmax(1:3),            &
1778                                       lb_child(1:3), ub_child(1:3),        &
1779                                        s_child(1:3),  s_parent(1:3),       &
1780                                       ds_child(1:3), ds_parent(1:3))
1781    enddo
1782!
1783    tempP = 0.
1784!
1785    do k = indmin(3), indmax(3)
1786    do j = indmin(2), indmax(2)
1787    do i = indmin(1), indmax(1)
1788        call Agrif_UpdateBase(type_update(4),                               &
1789                              tempP(i,j,k,indmin(4):indmax(4)),             &
1790                              tabtemp(i,j,k,lb_child(4):ub_child(4)),       &
1791                              indmin(4), indmax(4),                         &
1792                              lb_child(4), ub_child(4),                     &
1793                               s_parent(4), s_child(4),                     &
1794                              ds_parent(4),ds_child(4) )
1795    enddo
1796    enddo
1797    enddo
1798!
1799    deallocate(tabtemp)
1800!---------------------------------------------------------------------------------------------------
1801end subroutine Agrif_Update_4D_Recursive
1802!===================================================================================================
1803!
1804!===================================================================================================
1805!  subroutine Agrif_Update_5D_Recursive
1806!
1807!> Updates a 5D grid variable on the parent grid.
1808!! Calls #Agrif_Update_4D_Recursive and #Agrif_UpdateBase.
1809!---------------------------------------------------------------------------------------------------
1810subroutine Agrif_Update_5D_Recursive ( type_update,         &
1811                                       tempP, tempC,        &
1812                                       indmin, indmax,      &
1813                                       lb_child, ub_child,  &
1814                                        s_child,  s_parent, &
1815                                       ds_child, ds_parent )
1816!---------------------------------------------------------------------------------------------------
1817    integer, dimension(5),          intent(in)  :: type_update            !< Type of update (copy or average)
1818    integer, dimension(5),          intent(in)  :: indmin, indmax
1819    integer, dimension(5),          intent(in)  :: lb_child, ub_child
1820    real,    dimension(5),          intent(in)  ::  s_child,  s_parent
1821    real,    dimension(5),          intent(in)  :: ds_child, ds_parent
1822    real,    dimension(          &
1823        indmin(1):indmax(1),     &
1824        indmin(2):indmax(2),     &
1825        indmin(3):indmax(3),     &
1826        indmin(4):indmax(4),     &
1827        indmin(5):indmax(5)),       intent(out) :: tempP
1828    real, dimension(             &
1829        lb_child(1):ub_child(1), &
1830        lb_child(2):ub_child(2), &
1831        lb_child(3):ub_child(3), &
1832        lb_child(4):ub_child(4), &
1833        lb_child(5):ub_child(5)),   intent(in)  :: tempC
1834!---------------------------------------------------------------------------------------------------
1835    real, dimension(:,:,:,:,:), allocatable     :: tabtemp
1836    integer :: i,j,k,l,m
1837!
1838    allocate(tabtemp(indmin(1):indmax(1), &
1839                     indmin(2):indmax(2), &
1840                     indmin(3):indmax(3), &
1841                     indmin(4):indmax(4), &
1842                     lb_child(5):ub_child(5)))
1843!
1844    do m = lb_child(5), ub_child(5)
1845        call Agrif_Update_4D_Recursive(type_update(1:4),                    &
1846                                       tabtemp(indmin(1):indmax(1),         &
1847                                               indmin(2):indmax(2),         &
1848                                               indmin(3):indmax(3),         &
1849                                               indmin(4):indmax(4), m),     &
1850                                       tempC(lb_child(1):ub_child(1),       &
1851                                             lb_child(2):ub_child(2),       &
1852                                             lb_child(3):ub_child(3),       &
1853                                             lb_child(4):ub_child(4), m),   &
1854                                       indmin(1:4),indmax(1:4),             &
1855                                       lb_child(1:4), ub_child(1:4),        &
1856                                        s_child(1:4), s_parent(1:4),        &
1857                                       ds_child(1:4), ds_parent(1:4))
1858    enddo
1859!
1860    tempP = 0.
1861!
1862    do l = indmin(4), indmax(4)
1863    do k = indmin(3), indmax(3)
1864    do j = indmin(2), indmax(2)
1865    do i = indmin(1), indmax(1)
1866        call Agrif_UpdateBase( type_update(5),                              &
1867                               tempP(i,j,k,l,indmin(5):indmax(5)),          &
1868                               tabtemp(i,j,k,l,lb_child(5):ub_child(5)),    &
1869                               indmin(5), indmax(5),                        &
1870                               lb_child(5), ub_child(5),                    &
1871                                s_parent(5), s_child(5),                    &
1872                               ds_parent(5),ds_child(5) )
1873    enddo
1874    enddo
1875    enddo
1876    enddo
1877!
1878    deallocate(tabtemp)
1879!---------------------------------------------------------------------------------------------------
1880end subroutine Agrif_Update_5D_Recursive
1881!===================================================================================================
1882!
1883!===================================================================================================
1884!  subroutine Agrif_Update_6D_Recursive
1885!
1886!> Updates a 6D grid variable on the parent grid.
1887!! Calls #Agrif_Update_5D_Recursive and #Agrif_UpdateBase.
1888!---------------------------------------------------------------------------------------------------
1889subroutine Agrif_Update_6D_Recursive ( type_update,         &
1890                                       tempP, tempC,        &
1891                                       indmin, indmax,      &
1892                                       lb_child, ub_child,  &
1893                                        s_child,  s_parent, &
1894                                       ds_child, ds_parent )
1895!---------------------------------------------------------------------------------------------------
1896    integer, dimension(6),          intent(in)  :: type_update            !< Type of update (copy or average)
1897    integer, dimension(6),          intent(in)  :: indmin, indmax
1898    integer, dimension(6),          intent(in)  :: lb_child, ub_child
1899    real,    dimension(6),          intent(in)  ::  s_child,  s_parent
1900    real,    dimension(6),          intent(in)  :: ds_child, ds_parent
1901    real,    dimension(          &
1902        indmin(1):indmax(1),     &
1903        indmin(2):indmax(2),     &
1904        indmin(3):indmax(3),     &
1905        indmin(4):indmax(4),     &
1906        indmin(5):indmax(5),     &
1907        indmin(6):indmax(6)),       intent(out) :: tempP
1908    real, dimension(             &
1909        lb_child(1):ub_child(1), &
1910        lb_child(2):ub_child(2), &
1911        lb_child(3):ub_child(3), &
1912        lb_child(4):ub_child(4), &
1913        lb_child(5):ub_child(5), &
1914        lb_child(6):ub_child(6)),   intent(in)  :: tempC
1915!---------------------------------------------------------------------------------------------------
1916    real, dimension(:,:,:,:,:,:), allocatable   :: tabtemp
1917    integer :: i,j,k,l,m,n
1918!
1919    allocate(tabtemp(indmin(1):indmax(1), &
1920                     indmin(2):indmax(2), &
1921                     indmin(3):indmax(3), &
1922                     indmin(4):indmax(4), &
1923                     indmin(5):indmax(5), &
1924                     lb_child(6):ub_child(6)))
1925!
1926    do n = lb_child(6),ub_child(6)
1927        call Agrif_Update_5D_Recursive(type_update(1:5),                    &
1928                                       tabtemp(indmin(1):indmax(1),         &
1929                                               indmin(2):indmax(2),         &
1930                                               indmin(3):indmax(3),         &
1931                                               indmin(4):indmax(4),         &
1932                                               indmin(5):indmax(5), n),     &
1933                                       tempC(lb_child(1):ub_child(1),       &
1934                                             lb_child(2):ub_child(2),       &
1935                                             lb_child(3):ub_child(3),       &
1936                                             lb_child(4):ub_child(4),       &
1937                                             lb_child(5):ub_child(5), n),   &
1938                                       indmin(1:5), indmax(1:5),            &
1939                                       lb_child(1:5),ub_child(1:5),         &
1940                                       s_child(1:5), s_parent(1:5),         &
1941                                      ds_child(1:5),ds_parent(1:5))
1942    enddo
1943!
1944    tempP = 0.
1945!
1946    do m = indmin(5), indmax(5)
1947    do l = indmin(4), indmax(4)
1948    do k = indmin(3), indmax(3)
1949    do j = indmin(2), indmax(2)
1950    do i = indmin(1), indmax(1)
1951        call Agrif_UpdateBase( type_update(6),                              &
1952                               tempP(i,j,k,l,m,indmin(6):indmax(6)),        &
1953                               tabtemp(i,j,k,l,m,lb_child(6):ub_child(6)),  &
1954                               indmin(6), indmax(6),                        &
1955                               lb_child(6),  ub_child(6),                   &
1956                                s_parent(6),  s_child(6),                   &
1957                               ds_parent(6), ds_child(6) )
1958    enddo
1959    enddo
1960    enddo
1961    enddo
1962    enddo
1963!
1964    deallocate(tabtemp)
1965!---------------------------------------------------------------------------------------------------
1966end subroutine Agrif_Update_6D_Recursive
1967!===================================================================================================
1968!
1969!===================================================================================================
1970!  subroutine Agrif_UpdateBase
1971!
1972!> Calls the updating method chosen by the user (copy, average or full-weighting).
1973!---------------------------------------------------------------------------------------------------
1974subroutine Agrif_UpdateBase ( type_update,              &
1975                              parent_tab, child_tab,    &
1976                              indmin, indmax,           &
1977                              lb_child, ub_child,       &
1978                              s_parent, s_child,        &
1979                              ds_parent, ds_child )
1980!---------------------------------------------------------------------------------------------------
1981    integer,                            intent(in) :: type_update
1982    integer,                            intent(in) :: indmin, indmax
1983    integer,                            intent(in) :: lb_child, ub_child
1984    real, dimension(indmin:indmax),     intent(out):: parent_tab
1985    real, dimension(lb_child:ub_child), intent(in) :: child_tab
1986    real,                               intent(in) :: s_parent,  s_child
1987    real,                               intent(in) :: ds_parent, ds_child
1988!---------------------------------------------------------------------------------------------------
1989    integer :: np       ! Length of parent array
1990    integer :: nc       ! Length of child  array
1991!
1992    np = indmax - indmin + 1
1993    nc = ub_child - lb_child + 1
1994!
1995    if     ( type_update == Agrif_Update_Copy ) then
1996!
1997        call Agrif_basicupdate_copy1d(              &
1998                    parent_tab, child_tab,          &
1999                    np,         nc,                 &
2000                    s_parent,  s_child,             &
2001                    ds_parent, ds_child )
2002!
2003    elseif ( type_update == Agrif_Update_Average ) then
2004!
2005        call Agrif_basicupdate_average1d(           &
2006                    parent_tab, child_tab,          &
2007                    np,         nc,                 &
2008                     s_parent,  s_child,            &
2009                    ds_parent, ds_child )
2010!
2011    elseif ( type_update == Agrif_Update_Full_Weighting ) then
2012!
2013        call Agrif_basicupdate_full_weighting1D(    &
2014                    parent_tab, child_tab,          &
2015                    np,         nc,                 &
2016                    s_parent,  s_child,             &
2017                    ds_parent, ds_child )
2018!
2019    endif
2020!---------------------------------------------------------------------------------------------------
2021end subroutine Agrif_UpdateBase
2022!===================================================================================================
2023!
2024#if defined AGRIF_MPI
2025!===================================================================================================
2026!  subroutine Agrif_Find_list_update
2027!---------------------------------------------------------------------------------------------------
2028subroutine Agrif_Find_list_update ( list_update, pttab, petab, lb_child, lb_parent, nbdim, &
2029                                    find_list_update, tab4t, tab5t, memberinall, memberinall2,   &
2030                                    sendtoproc1, recvfromproc1, sendtoproc2, recvfromproc2 )
2031!---------------------------------------------------------------------------------------------------
2032    TYPE(Agrif_List_Interp_Loc),                   pointer     :: list_update
2033    INTEGER,                                       intent(in)  :: nbdim
2034    INTEGER, DIMENSION(nbdim),                     intent(in)  :: pttab, petab
2035    INTEGER, DIMENSION(nbdim),                     intent(in)  :: lb_child, lb_parent
2036    LOGICAL,                                       intent(out) :: find_list_update
2037    INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1,8), intent(out) :: tab4t
2038    INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1,8), intent(out) :: tab5t
2039    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(out) :: memberinall,memberinall2
2040    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(out) :: sendtoproc1,recvfromproc1
2041    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(out) :: sendtoproc2,recvfromproc2
2042!
2043    Type(Agrif_List_Interp_Loc), Pointer :: parcours
2044    INTEGER :: i
2045!
2046    find_list_update = .FALSE.
2047!
2048    parcours => list_update
2049
2050    Find_loop :  do while ( associated(parcours) )
2051        do i = 1,nbdim
2052            IF ((pttab(i) /= parcours%interp_loc%pttab(i)) .OR. &
2053                (petab(i) /= parcours%interp_loc%petab(i)) .OR. &
2054                (lb_child(i)  /= parcours%interp_loc%pttab_child(i)) .OR. &
2055                (lb_parent(i) /= parcours%interp_loc%pttab_parent(i))) THEN
2056                parcours => parcours%suiv
2057                cycle Find_loop
2058            ENDIF
2059        enddo
2060!
2061        tab4t = parcours%interp_loc%tab4t(1:nbdim,0:Agrif_Nbprocs-1,1:8)
2062        tab5t = parcours%interp_loc%tab5t(1:nbdim,0:Agrif_Nbprocs-1,1:8)
2063        memberinall  =  parcours%interp_loc%memberinall(0:Agrif_Nbprocs-1)
2064        memberinall2 =  parcours%interp_loc%memberinall2(0:Agrif_Nbprocs-1)
2065        sendtoproc1 =   parcours%interp_loc%sendtoproc1(0:Agrif_Nbprocs-1)
2066        sendtoproc2 =   parcours%interp_loc%sendtoproc2(0:Agrif_Nbprocs-1)
2067        recvfromproc1 = parcours%interp_loc%recvfromproc1(0:Agrif_Nbprocs-1)
2068        recvfromproc2 = parcours%interp_loc%recvfromproc2(0:Agrif_Nbprocs-1)
2069!
2070        find_list_update = .TRUE.
2071        exit Find_loop
2072!
2073    enddo Find_loop
2074!---------------------------------------------------------------------------------------------------
2075end subroutine Agrif_Find_list_update
2076!===================================================================================================
2077!
2078!===================================================================================================
2079!  subroutine Agrif_AddTo_list_update
2080!---------------------------------------------------------------------------------------------------
2081subroutine Agrif_AddTo_list_update ( list_update, pttab, petab, lb_child, lb_parent,  &
2082                                     nbdim, tab4t, tab5t, memberinall, memberinall2,        &
2083                                     sendtoproc1, recvfromproc1, sendtoproc2, recvfromproc2 )
2084!---------------------------------------------------------------------------------------------------
2085    TYPE(Agrif_List_Interp_Loc), pointer :: list_update
2086    INTEGER,                                        intent(in) :: nbdim
2087    INTEGER, DIMENSION(nbdim),                      intent(in) :: pttab, petab
2088    INTEGER, DIMENSION(nbdim),                      intent(in) :: lb_child, lb_parent
2089    INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1,8),  intent(in) :: tab4t
2090    INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1,8),  intent(in) :: tab5t
2091    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1),          intent(in) :: memberinall, memberinall2
2092    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1),          intent(in) :: sendtoproc1, recvfromproc1
2093    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1),          intent(in) :: sendtoproc2, recvfromproc2
2094!
2095    Type(Agrif_List_Interp_Loc), pointer :: parcours
2096!
2097    allocate(parcours)
2098    allocate(parcours%interp_loc)
2099
2100    parcours%interp_loc%pttab(1:nbdim) = pttab(1:nbdim)
2101    parcours%interp_loc%petab(1:nbdim) = petab(1:nbdim)
2102    parcours%interp_loc%pttab_child(1:nbdim)  = lb_child(1:nbdim)
2103    parcours%interp_loc%pttab_parent(1:nbdim) = lb_parent(1:nbdim)
2104
2105    allocate(parcours%interp_loc%tab4t(nbdim,0:Agrif_Nbprocs-1,8))
2106    allocate(parcours%interp_loc%tab5t(nbdim,0:Agrif_Nbprocs-1,8))
2107
2108    allocate(parcours%interp_loc%memberinall (0:Agrif_Nbprocs-1))
2109    allocate(parcours%interp_loc%memberinall2(0:Agrif_Nbprocs-1))
2110
2111    allocate(parcours%interp_loc%recvfromproc1(0:Agrif_Nbprocs-1))
2112    allocate(parcours%interp_loc%recvfromproc2(0:Agrif_Nbprocs-1))
2113    allocate(parcours%interp_loc%sendtoproc1(0:Agrif_Nbprocs-1))
2114    allocate(parcours%interp_loc%sendtoproc2(0:Agrif_Nbprocs-1))
2115
2116    parcours%interp_loc%tab4t = tab4t
2117    parcours%interp_loc%tab5t = tab5t
2118    parcours%interp_loc%memberinall   = memberinall
2119    parcours%interp_loc%memberinall2  = memberinall2
2120    parcours%interp_loc%sendtoproc1   = sendtoproc1
2121    parcours%interp_loc%sendtoproc2   = sendtoproc2
2122    parcours%interp_loc%recvfromproc1 = recvfromproc1
2123    parcours%interp_loc%recvfromproc2 = recvfromproc2
2124
2125    parcours%suiv => list_update
2126    list_update => parcours
2127!---------------------------------------------------------------------------------------------------
2128end subroutine Agrif_Addto_list_update
2129!===================================================================================================
2130#endif
2131!
2132end module Agrif_Update
Note: See TracBrowser for help on using the repository browser.