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/UKMO/r5936_hadgem3_cplfld/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES – NEMO

source: branches/UKMO/r5936_hadgem3_cplfld/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modupdate.F90 @ 7138

Last change on this file since 7138 was 7138, checked in by jcastill, 7 years ago

Remove svn keywords

File size: 101.6 KB
Line 
1!
2! $Id$
3!
4!     AGRIF (Adaptive Grid Refinement In Fortran)
5!
6!     Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
7!                        Christophe Vouland (Christophe.Vouland@imag.fr)
8!
9!     This program is free software; you can redistribute it and/or modify
10!     it under the terms of the GNU General Public License as published by
11!     the Free Software Foundation; either version 2 of the License, or
12!     (at your option) any later version.
13!
14!     This program is distributed in the hope that it will be useful,
15!     but WITHOUT ANY WARRANTY; without even the implied warranty of
16!     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
17!     GNU General Public License for more details.
18!
19!     You should have received a copy of the GNU General Public License
20!     along with this program; if not, write to the Free Software
21!     Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA.
22!
23!> 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    endif
586
587    call ExchangeSameLevel(sendtoproc1,recvfromproc1,nbdim,         &
588            tab4t(:,:,3),tab4t(:,:,4),tab4t(:,:,5),tab4t(:,:,6),    &
589            tab4t(:,:,7),tab4t(:,:,8),memberin,tempC,tempCextend)
590
591#else
592    tempCextend => tempC
593#endif
594!
595!     Update of the parent grid (tempP) from the child grid (tempC)
596!
597    IF (memberin) THEN
598!
599        IF ( .not.associated(tempP) ) allocate(tempP)
600!
601        call Agrif_array_allocate(tempP,indmin,indmax,nbdim)
602!
603        if ( nbdim == 1 ) then
604            tempP % array1 = 0.
605            call Agrif_Update_1D_Recursive( type_update(1),  &
606                                            tempP%array1,   &
607                                            tempCextend%array1, &
608                                            indmin(1), indmax(1),   &
609                                            pttruetabwhole(1), cetruetabwhole(1),   &
610                                            s_Child_temp(1), s_Parent_temp(1),      &
611                                            ds_child(1), ds_parent(1) )
612                                           
613            IF (Agrif_UseSpecialValueInUpdate) THEN
614            allocate(tempC_indic)
615            allocate(tempP_indic)
616            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array1),ubound(tempCextend%array1),nbdim)
617            call Agrif_array_allocate(tempP_indic,lbound(tempP%array1),ubound(tempP%array1),nbdim)
618
619            compute_average = .FALSE.
620            type_update_temp(1:nbdim) = type_update(1:nbdim)
621            IF (ANY(type_update(1:nbdim) == Agrif_Update_Full_Weighting)) THEN
622              compute_average = .TRUE.
623              allocate(tempP_average)
624              call Agrif_array_allocate(tempP_average,lbound(tempP%array1),ubound(tempP%array1),nbdim)
625              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
626                type_update_temp(1:nbdim) = Agrif_Update_Average
627              END WHERE
628              call Agrif_Update_1D_Recursive( type_update_temp(1),   &
629                                            tempP_average%array1,       &
630                                            tempCextend%array1, &
631                                            indmin(1), indmax(1),   &
632                                            pttruetabwhole(1), cetruetabwhole(1),   &
633                                            s_Child_temp(1), s_Parent_temp(1),      &
634                                            ds_child(1), ds_parent(1) )
635              coeff_multi = 1.
636              do nb_dimensions=1,nbdim
637                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
638              enddo
639            ENDIF
640           
641            WHERE (tempCextend%array1 == Agrif_SpecialValueFineGrid)
642              tempC_indic%array1 = 0.
643            ELSEWHERE
644              tempC_indic%array1 = 1.
645            END WHERE
646           
647            Agrif_UseSpecialValueInUpdate = .FALSE.
648            Agrif_Update_Weights = .TRUE.
649 
650             call Agrif_Update_1D_Recursive( type_update_temp(1),   &
651                                            tempP_indic%array1,       &
652                                            tempC_indic%array1, &
653                                            indmin(1), indmax(1),   &
654                                            pttruetabwhole(1), cetruetabwhole(1),   &
655                                            s_Child_temp(1), s_Parent_temp(1),      &
656                                            ds_child(1), ds_parent(1) )
657
658           Agrif_UseSpecialValueInUpdate = .TRUE.
659           Agrif_Update_Weights = .FALSE.
660
661           IF (compute_average) THEN
662               WHERE (tempP_indic%array1 == 0.)
663                  tempP%array1 = Agrif_SpecialValueFineGrid
664               ELSEWHERE ((tempP_indic%array1 == coeff_multi).AND.(tempP%array1 /= Agrif_SpecialValueFineGrid))
665                  tempP%array1 = tempP%array1 /tempP_indic%array1
666               ELSEWHERE
667                  tempP%array1 = tempP_average%array1 /tempP_indic%array1
668               END WHERE
669
670           ELSE
671               WHERE (tempP_indic%array1 == 0.)
672                  tempP%array1 = Agrif_SpecialValueFineGrid
673               ELSEWHERE
674                  tempP%array1 = tempP%array1 /tempP_indic%array1
675               END WHERE
676            ENDIF
677           
678            deallocate(tempP_indic%array1)
679            deallocate(tempC_indic%array1)
680            deallocate(tempC_indic)
681            deallocate(tempP_indic)
682            IF (compute_average) THEN
683              deallocate(tempP_average%array1)
684              deallocate(tempP_average)
685            ENDIF
686            ENDIF
687           
688        endif
689        if ( nbdim == 2 ) then
690            call Agrif_Update_2D_Recursive( type_update(1:2),   &
691                                            tempP%array2,       &
692                                            tempCextend%array2, &
693                                            indmin(1:2), indmax(1:2),   &
694                                            pttruetabwhole(1:2), cetruetabwhole(1:2),   &
695                                            s_Child_temp(1:2), s_Parent_temp(1:2),      &
696                                            ds_child(1:2), ds_parent(1:2) )
697
698            IF (Agrif_UseSpecialValueInUpdate) THEN
699            allocate(tempC_indic)
700            allocate(tempP_indic)
701            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array2),ubound(tempCextend%array2),nbdim)
702            call Agrif_array_allocate(tempP_indic,lbound(tempP%array2),ubound(tempP%array2),nbdim)
703 
704            compute_average = .FALSE.
705            type_update_temp(1:nbdim) = type_update(1:nbdim)
706            IF (ANY(type_update == Agrif_Update_Full_Weighting)) THEN
707              compute_average = .TRUE.
708              allocate(tempP_average)
709              call Agrif_array_allocate(tempP_average,lbound(tempP%array2),ubound(tempP%array2),nbdim)
710              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
711                type_update_temp(1:nbdim) = Agrif_Update_Average
712              END WHERE
713              call Agrif_Update_2D_Recursive( type_update_temp(1:2),   &
714                                            tempP_average%array2,       &
715                                            tempCextend%array2, &
716                                            indmin(1:2), indmax(1:2),   &
717                                            pttruetabwhole(1:2), cetruetabwhole(1:2),   &
718                                            s_Child_temp(1:2), s_Parent_temp(1:2),      &
719                                            ds_child(1:2), ds_parent(1:2) )
720              coeff_multi = 1.
721              do nb_dimensions=1,nbdim
722                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
723              enddo
724            ENDIF
725           
726            WHERE (tempCextend%array2 == Agrif_SpecialValueFineGrid)
727              tempC_indic%array2 = 0.
728            ELSEWHERE
729              tempC_indic%array2 = 1.
730            END WHERE
731           
732            Agrif_UseSpecialValueInUpdate = .FALSE.
733            Agrif_Update_Weights = .TRUE.
734           
735            call Agrif_Update_2D_Recursive( type_update_temp(1:2),   &
736                                            tempP_indic%array2,       &
737                                            tempC_indic%array2, &
738                                            indmin(1:2), indmax(1:2),   &
739                                            pttruetabwhole(1:2), cetruetabwhole(1:2),   &
740                                            s_Child_temp(1:2), s_Parent_temp(1:2),      &
741                                            ds_child(1:2), ds_parent(1:2) )
742
743           Agrif_UseSpecialValueInUpdate = .TRUE.
744           Agrif_Update_Weights = .FALSE.
745
746           IF (compute_average) THEN
747               WHERE (tempP_indic%array2 == 0.)
748                  tempP%array2 = Agrif_SpecialValueFineGrid
749               ELSEWHERE ((tempP_indic%array2 == coeff_multi).AND.(tempP%array2 /= Agrif_SpecialValueFineGrid))
750                  tempP%array2 = tempP%array2 /tempP_indic%array2
751               ELSEWHERE
752                  tempP%array2 = tempP_average%array2 /tempP_indic%array2
753               END WHERE
754
755           ELSE
756               WHERE (tempP_indic%array2 == 0.)
757                  tempP%array2 = Agrif_SpecialValueFineGrid
758               ELSEWHERE
759                  tempP%array2 = tempP%array2 /tempP_indic%array2
760               END WHERE
761            ENDIF
762           
763            deallocate(tempP_indic%array2)
764            deallocate(tempC_indic%array2)
765            deallocate(tempC_indic)
766            deallocate(tempP_indic)
767            IF (compute_average) THEN
768              deallocate(tempP_average%array2)
769              deallocate(tempP_average)
770            ENDIF
771            ENDIF
772           
773        endif
774        if ( nbdim == 3 ) then
775            call Agrif_Update_3D_Recursive( type_update(1:3),   &
776                                            tempP%array3,       &
777                                            tempCextend%array3, &
778                                            indmin(1:3), indmax(1:3),   &
779                                            pttruetabwhole(1:3), cetruetabwhole(1:3),   &
780                                            s_Child_temp(1:3), s_Parent_temp(1:3),      &
781                                            ds_child(1:3), ds_parent(1:3) )
782                                           
783            IF (Agrif_UseSpecialValueInUpdate) THEN
784            allocate(tempC_indic)
785            allocate(tempP_indic)
786            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array3),ubound(tempCextend%array3),nbdim)
787            call Agrif_array_allocate(tempP_indic,lbound(tempP%array3),ubound(tempP%array3),nbdim)
788
789            compute_average = .FALSE.
790            type_update_temp(1:nbdim) = type_update(1:nbdim)
791            IF (ANY(type_update == Agrif_Update_Full_Weighting)) THEN
792              compute_average = .TRUE.
793              allocate(tempP_average)
794              call Agrif_array_allocate(tempP_average,lbound(tempP%array3),ubound(tempP%array3),nbdim)
795              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
796                type_update_temp(1:nbdim) = Agrif_Update_Average
797              END WHERE
798              call Agrif_Update_3D_Recursive( type_update_temp(1:3),   &
799                                            tempP_average%array3,       &
800                                            tempCextend%array3, &
801                                            indmin(1:3), indmax(1:3),   &
802                                            pttruetabwhole(1:3), cetruetabwhole(1:3),   &
803                                            s_Child_temp(1:3), s_Parent_temp(1:3),      &
804                                            ds_child(1:3), ds_parent(1:3) )
805              coeff_multi = 1.
806              do nb_dimensions=1,nbdim
807                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
808              enddo
809            ENDIF
810           
811            WHERE (tempCextend%array3 == Agrif_SpecialValueFineGrid)
812              tempC_indic%array3 = 0.
813            ELSEWHERE
814              tempC_indic%array3 = 1.
815            END WHERE
816           
817            Agrif_UseSpecialValueInUpdate = .FALSE.
818            Agrif_Update_Weights = .TRUE.
819 
820             call Agrif_Update_3D_Recursive( type_update_temp(1:3),   &
821                                            tempP_indic%array3,       &
822                                            tempC_indic%array3, &
823                                            indmin(1:3), indmax(1:3),   &
824                                            pttruetabwhole(1:3), cetruetabwhole(1:3),   &
825                                            s_Child_temp(1:3), s_Parent_temp(1:3),      &
826                                            ds_child(1:3), ds_parent(1:3) )
827
828           Agrif_UseSpecialValueInUpdate = .TRUE.
829           Agrif_Update_Weights = .FALSE.
830
831           IF (compute_average) THEN
832               WHERE (tempP_indic%array3 == 0.)
833                  tempP%array3 = Agrif_SpecialValueFineGrid
834               ELSEWHERE ((tempP_indic%array3 == coeff_multi).AND.(tempP%array3 /= Agrif_SpecialValueFineGrid))
835                  tempP%array3 = tempP%array3 /tempP_indic%array3
836               ELSEWHERE
837                  tempP%array3 = tempP_average%array3 /tempP_indic%array3
838               END WHERE
839
840           ELSE
841               WHERE (tempP_indic%array3 == 0.)
842                  tempP%array3 = Agrif_SpecialValueFineGrid
843               ELSEWHERE
844                  tempP%array3 = tempP%array3 /tempP_indic%array3
845               END WHERE
846            ENDIF
847           
848            deallocate(tempP_indic%array3)
849            deallocate(tempC_indic%array3)
850            deallocate(tempC_indic)
851            deallocate(tempP_indic)
852            IF (compute_average) THEN
853              deallocate(tempP_average%array3)
854              deallocate(tempP_average)
855            ENDIF
856            ENDIF
857           
858        endif
859        if ( nbdim == 4 ) then
860            call Agrif_Update_4D_Recursive( type_update(1:4),   &
861                                            tempP%array4,       &
862                                            tempCextend%array4, &
863                                            indmin(1:4), indmax(1:4),   &
864                                            pttruetabwhole(1:4), cetruetabwhole(1:4),   &
865                                            s_Child_temp(1:4), s_Parent_temp(1:4),      &
866                                            ds_child(1:4), ds_parent(1:4) )
867                                           
868            IF (Agrif_UseSpecialValueInUpdate) THEN
869           
870            allocate(tempC_indic)
871            allocate(tempP_indic)
872            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array4),ubound(tempCextend%array4),nbdim)
873            call Agrif_array_allocate(tempP_indic,lbound(tempP%array4),ubound(tempP%array4),nbdim)
874           
875            compute_average = .FALSE.
876            type_update_temp(1:nbdim) = type_update(1:nbdim)
877            IF (ANY(type_update == Agrif_Update_Full_Weighting)) THEN
878              compute_average = .TRUE.
879              allocate(tempP_average)
880              call Agrif_array_allocate(tempP_average,lbound(tempP%array4),ubound(tempP%array4),nbdim)
881              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
882                type_update_temp(1:nbdim) = Agrif_Update_Average
883              END WHERE
884              call Agrif_Update_4D_Recursive( type_update_temp(1:4),   &
885                                            tempP_average%array4,       &
886                                            tempCextend%array4, &
887                                            indmin(1:4), indmax(1:4),   &
888                                            pttruetabwhole(1:4), cetruetabwhole(1:4),   &
889                                            s_Child_temp(1:4), s_Parent_temp(1:4),      &
890                                            ds_child(1:4), ds_parent(1:4) )
891              coeff_multi = 1.
892              do nb_dimensions=1,nbdim
893                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
894              enddo
895            ENDIF
896           
897            WHERE (tempCextend%array4 == Agrif_SpecialValueFineGrid)
898              tempC_indic%array4 = 0.
899            ELSEWHERE
900              tempC_indic%array4 = 1.
901            END WHERE
902           
903            Agrif_UseSpecialValueInUpdate = .FALSE.
904            Agrif_Update_Weights = .TRUE.
905 
906             call Agrif_Update_4D_Recursive( type_update_temp(1:4),   &
907                                            tempP_indic%array4,       &
908                                            tempC_indic%array4, &
909                                            indmin(1:4), indmax(1:4),   &
910                                            pttruetabwhole(1:4), cetruetabwhole(1:4),   &
911                                            s_Child_temp(1:4), s_Parent_temp(1:4),      &
912                                            ds_child(1:4), ds_parent(1:4) )
913
914           Agrif_UseSpecialValueInUpdate = .TRUE.
915           Agrif_Update_Weights = .FALSE.
916           
917           IF (compute_average) THEN
918               WHERE (tempP_indic%array4 == 0.)
919                  tempP%array4 = Agrif_SpecialValueFineGrid
920               ELSEWHERE ((tempP_indic%array4 == coeff_multi).AND.(tempP%array4 /= Agrif_SpecialValueFineGrid))
921                  tempP%array4 = tempP%array4 /tempP_indic%array4
922               ELSEWHERE
923                  tempP%array4 = tempP_average%array4 /tempP_indic%array4
924               END WHERE
925
926           ELSE
927               WHERE (tempP_indic%array4 == 0.)
928                  tempP%array4 = Agrif_SpecialValueFineGrid
929               ELSEWHERE
930                  tempP%array4 = tempP%array4 /tempP_indic%array4
931               END WHERE
932            ENDIF
933            deallocate(tempP_indic%array4)
934            deallocate(tempC_indic%array4)
935            deallocate(tempC_indic)
936            deallocate(tempP_indic)
937            IF (compute_average) THEN
938              deallocate(tempP_average%array4)
939              deallocate(tempP_average)
940            ENDIF
941            ENDIF
942           
943        endif
944        if ( nbdim == 5 ) then
945            call Agrif_Update_5D_Recursive( type_update(1:5),   &
946                                            tempP%array5,       &
947                                            tempCextend%array5, &
948                                            indmin(1:5), indmax(1:5),   &
949                                            pttruetabwhole(1:5), cetruetabwhole(1:5),   &
950                                            s_Child_temp(1:5), s_Parent_temp(1:5),      &
951                                            ds_child(1:5), ds_parent(1:5) )
952                                           
953            IF (Agrif_UseSpecialValueInUpdate) THEN
954            allocate(tempC_indic)
955            allocate(tempP_indic)
956            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array5),ubound(tempCextend%array5),nbdim)
957            call Agrif_array_allocate(tempP_indic,lbound(tempP%array5),ubound(tempP%array5),nbdim)
958
959            compute_average = .FALSE.
960            type_update_temp(1:nbdim) = type_update(1:nbdim)
961            IF (ANY(type_update == Agrif_Update_Full_Weighting)) THEN
962              compute_average = .TRUE.
963              allocate(tempP_average)
964              call Agrif_array_allocate(tempP_average,lbound(tempP%array5),ubound(tempP%array5),nbdim)
965              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
966                type_update_temp(1:nbdim) = Agrif_Update_Average
967              END WHERE
968              call Agrif_Update_5D_Recursive( type_update_temp(1:5),   &
969                                            tempP_average%array5,       &
970                                            tempCextend%array5, &
971                                            indmin(1:5), indmax(1:5),   &
972                                            pttruetabwhole(1:5), cetruetabwhole(1:5),   &
973                                            s_Child_temp(1:5), s_Parent_temp(1:5),      &
974                                            ds_child(1:5), ds_parent(1:5) )
975              coeff_multi = 1.
976              do nb_dimensions=1,nbdim
977                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
978              enddo
979            ENDIF
980           
981            WHERE (tempCextend%array5 == Agrif_SpecialValueFineGrid)
982              tempC_indic%array5 = 0.
983            ELSEWHERE
984              tempC_indic%array5 = 1.
985            END WHERE
986           
987            Agrif_UseSpecialValueInUpdate = .FALSE.
988            Agrif_Update_Weights = .TRUE.
989 
990             call Agrif_Update_5D_Recursive( type_update_temp(1:5),   &
991                                            tempP_indic%array5,       &
992                                            tempC_indic%array5, &
993                                            indmin(1:5), indmax(1:5),   &
994                                            pttruetabwhole(1:5), cetruetabwhole(1:5),   &
995                                            s_Child_temp(1:5), s_Parent_temp(1:5),      &
996                                            ds_child(1:5), ds_parent(1:5) )
997
998           Agrif_UseSpecialValueInUpdate = .TRUE.
999           Agrif_Update_Weights = .FALSE.
1000
1001           IF (compute_average) THEN
1002               WHERE (tempP_indic%array5 == 0.)
1003                  tempP%array5 = Agrif_SpecialValueFineGrid
1004               ELSEWHERE ((tempP_indic%array5 == coeff_multi).AND.(tempP%array5 /= Agrif_SpecialValueFineGrid))
1005                  tempP%array5 = tempP%array5 /tempP_indic%array5
1006               ELSEWHERE
1007                  tempP%array5 = tempP_average%array5 /tempP_indic%array5
1008               END WHERE
1009
1010           ELSE
1011               WHERE (tempP_indic%array5 == 0.)
1012                  tempP%array5 = Agrif_SpecialValueFineGrid
1013               ELSEWHERE
1014                  tempP%array5 = tempP%array5 /tempP_indic%array5
1015               END WHERE
1016            ENDIF
1017           
1018            deallocate(tempP_indic%array5)
1019            deallocate(tempC_indic%array5)
1020            deallocate(tempC_indic)
1021            deallocate(tempP_indic)
1022            IF (compute_average) THEN
1023              deallocate(tempP_average%array5)
1024              deallocate(tempP_average)
1025            ENDIF
1026            ENDIF
1027           
1028        endif
1029        if ( nbdim == 6 ) then
1030            call Agrif_Update_6D_Recursive( type_update(1:6),   &
1031                                            tempP%array6,       &
1032                                            tempCextend%array6, &
1033                                            indmin(1:6), indmax(1:6),   &
1034                                            pttruetabwhole(1:6), cetruetabwhole(1:6),   &
1035                                            s_Child_temp(1:6), s_Parent_temp(1:6),      &
1036                                            ds_child(1:6), ds_parent(1:6) )
1037            IF (Agrif_UseSpecialValueInUpdate) THEN
1038            allocate(tempC_indic)
1039            allocate(tempP_indic)
1040            call Agrif_array_allocate(tempC_indic,lbound(tempCextend%array6),ubound(tempCextend%array6),nbdim)
1041            call Agrif_array_allocate(tempP_indic,lbound(tempP%array6),ubound(tempP%array6),nbdim)
1042
1043            compute_average = .FALSE.
1044            type_update_temp(1:nbdim) = type_update(1:nbdim)
1045            IF (ANY(type_update == Agrif_Update_Full_Weighting)) THEN
1046              compute_average = .TRUE.
1047              allocate(tempP_average)
1048              call Agrif_array_allocate(tempP_average,lbound(tempP%array6),ubound(tempP%array6),nbdim)
1049              type_update_temp(1:nbdim) = type_update
1050              WHERE (type_update(1:nbdim) == Agrif_Update_Full_Weighting)
1051                type_update_temp(1:nbdim) = Agrif_Update_Average
1052              END WHERE
1053              call Agrif_Update_6D_Recursive( type_update_temp(1:6),   &
1054                                            tempP_average%array6,       &
1055                                            tempCextend%array6, &
1056                                            indmin(1:6), indmax(1:6),   &
1057                                            pttruetabwhole(1:6), cetruetabwhole(1:6),   &
1058                                            s_Child_temp(1:6), s_Parent_temp(1:6),      &
1059                                            ds_child(1:6), ds_parent(1:6) )
1060              coeff_multi = 1.
1061              do nb_dimensions=1,nbdim
1062                coeff_multi = coeff_multi * nint(ds_parent(nb_dimensions)/ds_child(nb_dimensions))
1063              enddo
1064            ENDIF
1065
1066           IF (compute_average) THEN
1067               WHERE (tempP_indic%array6 == 0.)
1068                  tempP%array6 = Agrif_SpecialValueFineGrid
1069               ELSEWHERE ((tempP_indic%array6 == coeff_multi).AND.(tempP%array6 /= Agrif_SpecialValueFineGrid))
1070                  tempP%array6 = tempP%array6 /tempP_indic%array6
1071               ELSEWHERE
1072                  tempP%array6 = tempP_average%array6 /tempP_indic%array6
1073               END WHERE
1074
1075           ELSE
1076               WHERE (tempP_indic%array6 == 0.)
1077                  tempP%array6 = Agrif_SpecialValueFineGrid
1078               ELSEWHERE
1079                  tempP%array6 = tempP%array6 /tempP_indic%array6
1080               END WHERE
1081            ENDIF
1082           
1083            Agrif_UseSpecialValueInUpdate = .FALSE.
1084            Agrif_Update_Weights = .TRUE.
1085 
1086             call Agrif_Update_6D_Recursive( type_update_temp(1:6),   &
1087                                            tempP_indic%array6,       &
1088                                            tempC_indic%array6, &
1089                                            indmin(1:6), indmax(1:6),   &
1090                                            pttruetabwhole(1:6), cetruetabwhole(1:6),   &
1091                                            s_Child_temp(1:6), s_Parent_temp(1:6),      &
1092                                            ds_child(1:6), ds_parent(1:6) )
1093
1094           Agrif_UseSpecialValueInUpdate = .TRUE.
1095           Agrif_Update_Weights = .FALSE.
1096           
1097            WHERE (tempP_indic%array6 == 0.)
1098              tempP%array6 = Agrif_SpecialValueFineGrid
1099            ELSEWHERE
1100              tempP%array6 = tempP%array6 /tempP_indic%array6
1101            END WHERE
1102           
1103            deallocate(tempP_indic%array6)
1104            deallocate(tempC_indic%array6)
1105            deallocate(tempC_indic)
1106            deallocate(tempP_indic)
1107            IF (compute_average) THEN
1108              deallocate(tempP_average%array6)
1109              deallocate(tempP_average)
1110            ENDIF
1111            ENDIF
1112        endif
1113!
1114        call Agrif_array_deallocate(tempCextend,nbdim)
1115!
1116    ENDIF
1117
1118#if defined AGRIF_MPI
1119    local_proc = Agrif_Procrank
1120    call Agrif_get_var_bounds_array(parent,lowerbound,upperbound,nbdim)
1121    call Agrif_ChildGrid_to_ParentGrid()
1122    call Agrif_Childbounds(nbdim, lowerbound, upperbound,                   &
1123                           indminglob,  indmaxglob,  local_proc, coords,    &
1124                           indminglob2, indmaxglob2, member)
1125!
1126    IF (member) THEN
1127        call Agrif_GlobalToLocalBounds(parentarray, lowerbound, upperbound, &
1128                                       indminglob2, indmaxglob2, coords,    &
1129                                       nbdim, local_proc, member)
1130    ENDIF
1131
1132    call Agrif_ParentGrid_to_ChildGrid()
1133
1134    if (.not.find_list_update) then
1135        tab3(:,1) = indmin(:)
1136        tab3(:,2) = indmax(:)
1137        tab3(:,3) = indminglob2(:)
1138        tab3(:,4) = indmaxglob2(:)
1139!
1140        call MPI_ALLGATHER(tab3,4*nbdim,MPI_INTEGER,tab4,4*nbdim,MPI_INTEGER,Agrif_mpi_comm,code)
1141
1142        IF ( .not.associated(tempPextend) ) allocate(tempPextend)
1143        DO k=0,Agrif_Nbprocs-1
1144            do j=1,4
1145                do i=1,nbdim
1146                    tab5t(i,k,j) = tab4(i,j,k)
1147                enddo
1148            enddo
1149        enddo
1150
1151        memberin1(1) = member
1152        call MPI_ALLGATHER(memberin1,1,MPI_LOGICAL,memberinall2,1,MPI_LOGICAL,Agrif_mpi_comm,code)
1153        call Get_External_Data_first(tab5t(:,:,1),tab5t(:,:,2),tab5t(:,:,3),tab5t(:,:,4),   &
1154                                     nbdim, memberinall2, coords,                           &
1155                                     sendtoproc2, recvfromproc2,                            &
1156                                     tab5t(:,:,5),tab5t(:,:,6),tab5t(:,:,7),tab5t(:,:,8))
1157
1158        call Agrif_Addto_list_update(child%list_update,pttab,petab,lb_child,lb_parent,      &
1159                                     nbdim,tab4t,tab5t,memberinall,memberinall2,            &
1160                                     sendtoproc1,recvfromproc1,sendtoproc2,recvfromproc2)
1161
1162    endif
1163
1164    call ExchangeSameLevel(sendtoproc2,recvfromproc2,nbdim,                     &
1165                            tab5t(:,:,3),tab5t(:,:,4),tab5t(:,:,5),tab5t(:,:,6),&
1166                            tab5t(:,:,7),tab5t(:,:,8),member,tempP,tempPextend)
1167#else
1168    tempPextend => tempP
1169    parentarray(:,1,1) = indmin
1170    parentarray(:,2,1) = indmax
1171    parentarray(:,1,2) = indmin
1172    parentarray(:,2,2) = indmax
1173    member = .TRUE.
1174#endif
1175!
1176!   Special values on the child grid
1177    if ( Agrif_UseSpecialValueFineGrid ) then
1178!
1179!cc         noraftab(1:nbdim) =
1180!cc     &    child % root_var % interptab(1:nbdim) == 'N'
1181!
1182#if defined AGRIF_MPI
1183!
1184!          allocate(childvalues% var)
1185!
1186!          Call Agrif_array_allocate(childvalues%var,
1187!     &                      pttruetab,cetruetab,nbdim)
1188!          Call Agrif_var_full_copy_array(childvalues% var,
1189!     &                                tempC,
1190!     &                                nbdim)
1191!          Call Agrif_CheckMasknD(tempC,childvalues,
1192!     &                           pttruetab(1:nbdim),cetruetab(1:nbdim),
1193!     &                           pttruetab(1:nbdim),cetruetab(1:nbdim),
1194!     &                           noraftab(1:nbdim),nbdim)
1195!          Call Agrif_array_deallocate(childvalues% var,nbdim)
1196!         Deallocate(childvalues % var)
1197!
1198#else
1199!
1200!          Call Agrif_get_var_bounds_array(child,
1201!     &                              lowerbound,upperbound,nbdim)
1202!          Call Agrif_CheckMasknD(tempC,child,
1203!     &                           pttruetab(1:nbdim),cetruetab(1:nbdim),
1204!     &                           lowerbound,
1205!     &                           upperbound,
1206!     &                           noraftab(1:nbdim),nbdim)
1207!
1208#endif
1209!
1210    endif
1211!
1212!   Special values on the parent grid
1213    if (Agrif_UseSpecialValue) then
1214!
1215#if defined AGRIF_MPI
1216!
1217!          Call GiveAgrif_SpecialValueToTab_mpi(parent,tempP,
1218!     &                 parentarray,
1219!     &                 Agrif_SpecialValue,nbdim)
1220!
1221!
1222#else
1223!
1224!          Call GiveAgrif_SpecialValueToTab(parent,tempP,
1225!     &                  indmin,indmax,
1226!     &                  Agrif_SpecialValue,nbdim)
1227!
1228#endif
1229!
1230    endif
1231!
1232    IF (member) THEN
1233
1234        call Agrif_ChildGrid_to_ParentGrid()
1235!
1236        SELECT CASE(nbdim)
1237        CASE(1)
1238            call procname( tempPextend % array1(            &
1239                    parentarray(1,1,1):parentarray(1,2,1)), &
1240                    parentarray(1,1,2),parentarray(1,2,2),.FALSE.,nbin,ndirin)
1241        CASE(2)
1242            call procname( tempPextend % array2(            &
1243                    parentarray(1,1,1):parentarray(1,2,1),  &
1244                    parentarray(2,1,1):parentarray(2,2,1)), &
1245                    parentarray(1,1,2),parentarray(1,2,2),  &
1246                    parentarray(2,1,2),parentarray(2,2,2),.FALSE.,nbin,ndirin)
1247        CASE(3)
1248            call procname( tempPextend % array3(            &
1249                    parentarray(1,1,1):parentarray(1,2,1),  &
1250                    parentarray(2,1,1):parentarray(2,2,1),  &
1251                    parentarray(3,1,1):parentarray(3,2,1)), &
1252                    parentarray(1,1,2),parentarray(1,2,2),  &
1253                    parentarray(2,1,2),parentarray(2,2,2),  &
1254                    parentarray(3,1,2),parentarray(3,2,2),.FALSE.,nbin,ndirin)
1255        CASE(4)
1256            call procname( tempPextend % array4(            &
1257                    parentarray(1,1,1):parentarray(1,2,1),  &
1258                    parentarray(2,1,1):parentarray(2,2,1),  &
1259                    parentarray(3,1,1):parentarray(3,2,1),  &
1260                    parentarray(4,1,1):parentarray(4,2,1)), &
1261                    parentarray(1,1,2),parentarray(1,2,2),  &
1262                    parentarray(2,1,2),parentarray(2,2,2),  &
1263                    parentarray(3,1,2),parentarray(3,2,2),  &
1264                    parentarray(4,1,2),parentarray(4,2,2),.FALSE.,nbin,ndirin)
1265        CASE(5)
1266            call procname( tempPextend % array5(            &
1267                    parentarray(1,1,1):parentarray(1,2,1),  &
1268                    parentarray(2,1,1):parentarray(2,2,1),  &
1269                    parentarray(3,1,1):parentarray(3,2,1),  &
1270                    parentarray(4,1,1):parentarray(4,2,1),  &
1271                    parentarray(5,1,1):parentarray(5,2,1)), &
1272                    parentarray(1,1,2),parentarray(1,2,2),  &
1273                    parentarray(2,1,2),parentarray(2,2,2),  &
1274                    parentarray(3,1,2),parentarray(3,2,2),  &
1275                    parentarray(4,1,2),parentarray(4,2,2),  &
1276                    parentarray(5,1,2),parentarray(5,2,2),.FALSE.,nbin,ndirin)
1277        CASE(6)
1278            call procname( tempPextend % array6(            &
1279                    parentarray(1,1,1):parentarray(1,2,1),  &
1280                    parentarray(2,1,1):parentarray(2,2,1),  &
1281                    parentarray(3,1,1):parentarray(3,2,1),  &
1282                    parentarray(4,1,1):parentarray(4,2,1),  &
1283                    parentarray(5,1,1):parentarray(5,2,1),  &
1284                    parentarray(6,1,1):parentarray(6,2,1)), &
1285                    parentarray(1,1,2),parentarray(1,2,2),  &
1286                    parentarray(2,1,2),parentarray(2,2,2),  &
1287                    parentarray(3,1,2),parentarray(3,2,2),  &
1288                    parentarray(4,1,2),parentarray(4,2,2),  &
1289                    parentarray(5,1,2),parentarray(5,2,2),  &
1290                    parentarray(6,1,2),parentarray(6,2,2),.FALSE.,nbin,ndirin)
1291        END SELECT
1292!
1293        call Agrif_ParentGrid_to_ChildGrid()
1294!
1295        call Agrif_array_deallocate(tempPextend,nbdim)
1296!
1297    ENDIF
1298!
1299#if defined AGRIF_MPI
1300    IF (memberin) THEN
1301        call Agrif_array_deallocate(tempP,nbdim)
1302        call Agrif_array_deallocate(tempC,nbdim)
1303    ENDIF
1304#endif
1305!---------------------------------------------------------------------------------------------------
1306end subroutine Agrif_UpdatenD
1307!===================================================================================================
1308!
1309!===================================================================================================
1310!  subroutine Agrif_Prtbounds
1311!
1312!> calculates the bounds of the parent grid to be updated by the child grid
1313!---------------------------------------------------------------------------------------------------
1314subroutine Agrif_Prtbounds ( nbdim, indmin, indmax, s_Parent_temp, s_Child_temp,        &
1315                             s_child, ds_child, s_parent, ds_parent,                    &
1316                             pttruetab, cetruetab, lb_child, lb_parent            &
1317#if defined AGRIF_MPI
1318                            ,posvar, type_update, do_update,                &
1319                             pttruetabwhole, cetruetabwhole                             &
1320#endif
1321                )
1322!---------------------------------------------------------------------------------------------------
1323    integer,                   intent(in)   :: nbdim
1324    integer, dimension(nbdim), intent(out)  :: indmin, indmax
1325    real,    dimension(nbdim), intent(out)  :: s_Parent_temp, s_Child_temp
1326    real,    dimension(nbdim), intent(in)   :: s_child,  ds_child
1327    real,    dimension(nbdim), intent(in)   :: s_parent, ds_parent
1328    integer, dimension(nbdim), intent(in)   :: pttruetab, cetruetab
1329    integer, dimension(nbdim), intent(in)   :: lb_child, lb_parent
1330#if defined AGRIF_MPI
1331    integer, dimension(nbdim), intent(in)   :: posvar   !< Position of the variable on the cell (1 or 2)
1332    integer, dimension(nbdim), intent(in)   :: type_update
1333    logical, dimension(nbdim), intent(in)   :: do_update
1334    integer,dimension(nbdim), intent(out)   :: pttruetabwhole, cetruetabwhole
1335#endif
1336!
1337    real,dimension(nbdim) :: dim_newmin,dim_newmax
1338    integer :: i
1339#if defined AGRIF_MPI
1340    real    :: positionmin, positionmax
1341    integer :: imin, imax
1342    integer :: coeffraf
1343#endif
1344!
1345    do i = 1,nbdim
1346!
1347        dim_newmin(i) = s_child(i) + (pttruetab(i) - lb_child(i)) * ds_child(i)
1348        dim_newmax(i) = s_child(i) + (cetruetab(i) - lb_child(i)) * ds_child(i)
1349!
1350        indmin(i) = lb_parent(i) + agrif_ceiling((dim_newmin(i)-s_parent(i))/ds_parent(i))
1351        indmax(i) = lb_parent(i) + agrif_int(    (dim_newmax(i)-s_parent(i))/ds_parent(i))
1352!
1353#if defined AGRIF_MPI
1354        positionmin = s_parent(i) + (indmin(i)-lb_parent(i))*ds_parent(i)
1355        IF ( do_update(i) ) THEN
1356            IF (posvar(i) == 1) THEN
1357                IF      (type_update(i) == Agrif_Update_Average) THEN
1358                    positionmin = positionmin - ds_parent(i)/2.
1359                ELSE IF (type_update(i) == Agrif_Update_Full_Weighting) THEN
1360                    positionmin = positionmin - (ds_parent(i)-ds_child(i))
1361                ENDIF
1362            ELSE
1363                IF (type_update(i) /= Agrif_Update_Full_Weighting) THEN
1364                    positionmin = positionmin - ds_parent(i)/2.
1365                ELSE
1366                    coeffraf = nint(ds_parent(i)/ds_child(i))
1367                    if (mod(coeffraf,2) == 1) then
1368                        positionmin = positionmin - (ds_parent(i)-ds_child(i))
1369                    else
1370                        positionmin = positionmin - (ds_parent(i)-ds_child(i))-ds_child(i)/2.
1371                    endif
1372                ENDIF
1373            ENDIF
1374        ENDIF
1375!
1376        imin = lb_child(i) + agrif_ceiling((positionmin-s_child(i))/ds_child(i))
1377        positionmin = s_child(i)  + (imin - lb_child(i)) * ds_child(i)
1378        positionmax = s_parent(i) + (indmax(i)-lb_parent(i))*ds_parent(i)
1379        pttruetabwhole(i) = imin
1380
1381        IF ( do_update(i) ) THEN
1382            IF (posvar(i) == 1) THEN
1383                IF      (type_update(i) == Agrif_Update_Average) THEN
1384                    positionmax = positionmax  + ds_parent(i)/2.
1385                ELSE IF (type_update(i) == Agrif_Update_Full_Weighting) THEN
1386                    positionmax = positionmax  + (ds_parent(i)-ds_child(i))
1387                ENDIF
1388            ELSE
1389                IF (type_update(i) /= Agrif_Update_Full_Weighting) THEN
1390                    positionmax = positionmax  + ds_parent(i)/2.
1391                ELSE
1392                    coeffraf = nint(ds_parent(i)/ds_child(i))
1393                    if (mod(coeffraf,2) == 1) then
1394                        positionmax = positionmax + (ds_parent(i)-ds_child(i))
1395                    else
1396                        positionmax = positionmax + (ds_parent(i)-ds_child(i)) + ds_child(i)/2.
1397                    endif
1398                ENDIF
1399            ENDIF
1400        ENDIF
1401
1402        imax = lb_child(i) +agrif_int((positionmax-s_child(i))/ds_child(i))
1403        positionmax = s_child(i) + (imax - lb_child(i)) * ds_child(i)
1404        cetruetabwhole(i) = imax
1405#endif
1406!
1407        s_Parent_temp(i) = s_parent(i) + (indmin(i) - lb_parent(i)) * ds_parent(i)
1408        s_Child_temp(i)  = dim_newmin(i)
1409
1410#if defined AGRIF_MPI
1411        s_Child_temp(i) = positionmin
1412#endif
1413!
1414    enddo
1415!---------------------------------------------------------------------------------------------------
1416end subroutine Agrif_Prtbounds
1417!===================================================================================================
1418!
1419!===================================================================================================
1420!  subroutine Agrif_Update_1D_Recursive
1421!
1422!> Updates a 1D grid variable on the parent grid
1423!---------------------------------------------------------------------------------------------------
1424subroutine Agrif_Update_1D_Recursive ( type_update,         &
1425                                       tempP, tempC,        &
1426                                       indmin, indmax,      &
1427                                       lb_child, ub_child,  &
1428                                       s_child,  s_parent,  &
1429                                       ds_child, ds_parent )
1430!---------------------------------------------------------------------------------------------------
1431    integer,                            intent(in)  :: type_update            !< Type of update (copy or average)
1432    integer,                            intent(in)  :: indmin, indmax
1433    integer,                            intent(in)  :: lb_child, ub_child
1434    real,                               intent(in)  ::  s_child,  s_parent
1435    real,                               intent(in)  :: ds_child, ds_parent
1436    real, dimension(indmin:indmax),     intent(out) :: tempP
1437    real, dimension(lb_child:ub_child), intent(in)  :: tempC
1438!---------------------------------------------------------------------------------------------------
1439    call Agrif_UpdateBase(type_update,              &
1440                          tempP(indmin:indmax),     &
1441                          tempC(lb_child:ub_child), &
1442                          indmin, indmax,           &
1443                          lb_child, ub_child,       &
1444                          s_parent, s_child,        &
1445                          ds_parent, ds_child)
1446!---------------------------------------------------------------------------------------------------
1447end subroutine Agrif_Update_1D_Recursive
1448!===================================================================================================
1449!
1450!===================================================================================================
1451!  subroutine Agrif_Update_2D_Recursive
1452!
1453!> updates a 2D grid variable on the parent grid.
1454!! Calls #Agrif_Update_1D_Recursive and #Agrif_UpdateBase
1455!---------------------------------------------------------------------------------------------------
1456subroutine Agrif_Update_2D_Recursive ( type_update,         &
1457                                       tempP, tempC,        &
1458                                       indmin, indmax,      &
1459                                       lb_child, ub_child,  &
1460                                        s_child,  s_parent, &
1461                                       ds_child, ds_parent )
1462!---------------------------------------------------------------------------------------------------
1463    integer, dimension(2),          intent(in)  :: type_update            !< Type of update (copy or average)
1464    integer, dimension(2),          intent(in)  :: indmin, indmax
1465    integer, dimension(2),          intent(in)  :: lb_child, ub_child
1466    real,    dimension(2),          intent(in)  ::  s_child,  s_parent
1467    real,    dimension(2),          intent(in)  :: ds_child, ds_parent
1468    real,    dimension(          &
1469        indmin(1):indmax(1),     &
1470        indmin(2):indmax(2)),       intent(out) :: tempP
1471    real,    dimension(:,:),        intent(in)  :: tempC
1472!---------------------------------------------------------------------------------------------------
1473    real, dimension(indmin(1):indmax(1), lb_child(2):ub_child(2)) :: tabtemp
1474    real, dimension(indmin(2):indmax(2), indmin(1):indmax(1))     :: tempP_trsp
1475    real, dimension(lb_child(2):ub_child(2), indmin(1):indmax(1)) :: tabtemp_trsp
1476    integer :: i, j
1477    integer :: coeffraf
1478!
1479    tabtemp = 0.
1480    coeffraf = nint ( ds_parent(1) / ds_child(1) )
1481!
1482    IF((type_update(1) == Agrif_Update_Average) .AND. (coeffraf /= 1 )) THEN
1483!---CDIR NEXPAND
1484        if ( .NOT. precomputedone(1) ) then
1485            call Average1dPrecompute( ub_child(2)-lb_child(2)+1,    &
1486                                      indmax(1)-indmin(1)+1,              &
1487                                      ub_child(1)-lb_child(1)+1,    &
1488                                      s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1489!            precomputedone(1) = .TRUE.
1490        endif
1491!---CDIR NEXPAND
1492        call Average1dAfterCompute( tabtemp, tempC, size(tabtemp), size(tempC), &
1493                    s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1494!
1495    ELSE IF ((type_update(1) == Agrif_Update_Copy) .AND. (coeffraf /= 1 ))THEN
1496!---CDIR NEXPAND
1497         if ( .NOT. precomputedone(1) ) then
1498            call Agrif_basicupdate_copy1d_before( ub_child(2)-lb_child(2)+1, &
1499                                   indmax(1)-indmin(1)+1,           &
1500                                   ub_child(1)-lb_child(1)+1, &
1501                                   s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1502!            precomputedone(1) = .TRUE.
1503        endif
1504!---CDIR NEXPAND
1505        call Agrif_basicupdate_copy1d_after(tabtemp,tempC,size(tabtemp),size(tempC),1)
1506!
1507    ELSE
1508        do j = lb_child(2),ub_child(2)
1509!
1510!---CDIR NEXPAND
1511            call Agrif_Update_1D_Recursive( type_update(1),             &
1512                                            tabtemp(:,j),               &
1513                                            tempC(:,j-lb_child(2)+1),   &
1514                                            indmin(1), indmax(1),       &
1515                                            lb_child(1),ub_child(1),    &
1516                                            s_child(1),  s_parent(1),   &
1517                                            ds_child(1),ds_parent(1))
1518        enddo
1519    ENDIF
1520!
1521    tabtemp_trsp = TRANSPOSE(tabtemp)
1522    coeffraf = nint(ds_parent(2)/ds_child(2))
1523!
1524    tempP_trsp = 0.
1525!
1526    IF((type_update(2) == Agrif_Update_Average) .AND. (coeffraf /= 1 )) THEN
1527!---CDIR NEXPAND
1528        if ( .NOT. precomputedone(2) ) then
1529            call Average1dPrecompute( indmax(1)-indmin(1)+1,          &
1530                                      indmax(2)-indmin(2)+1,          &
1531                                      ub_child(2)-lb_child(2)+1,&
1532                                      s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1533!            precomputedone(2) = .TRUE.
1534        endif
1535!---CDIR NEXPAND
1536        call Average1dAfterCompute( tempP_trsp, tabtemp_trsp, size(tempP_trsp), size(tabtemp_trsp),&
1537                s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1538!
1539    ELSE IF ((type_update(2) == Agrif_Update_Copy) .AND. (coeffraf /= 1 )) THEN
1540!---CDIR NEXPAND
1541        if ( .NOT. precomputedone(2) ) then
1542            call Agrif_basicupdate_copy1d_before( indmax(1)-indmin(1)+1,             &
1543                                   indmax(2)-indmin(2)+1,             &
1544                                   ub_child(2)-lb_child(2)+1,   &
1545                                   s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1546!            precomputedone(2) = .TRUE.
1547        endif
1548!---CDIR NEXPAND
1549        call Agrif_basicupdate_copy1d_after( tempP_trsp, tabtemp_trsp, size(tempP_trsp), size(tabtemp_trsp),2)
1550!
1551    ELSE
1552        do i = indmin(1),indmax(1)
1553!
1554!---CDIR NEXPAND
1555            call Agrif_UpdateBase(type_update(2),                                       &
1556                                  tempP_trsp(indmin(2):indmax(2),i),            &
1557                                  tabtemp_trsp(lb_child(2):ub_child(2),i),&
1558                                  indmin(2),indmax(2),                          &
1559                                  lb_child(2),ub_child(2),                &
1560                                  s_parent(2),s_child(2),                       &
1561                                  ds_parent(2),ds_child(2))
1562!
1563        enddo
1564    ENDIF
1565!
1566    tempP = TRANSPOSE(tempP_trsp)
1567!---------------------------------------------------------------------------------------------------
1568end subroutine Agrif_Update_2D_Recursive
1569!===================================================================================================
1570!
1571subroutine Agrif_Update_2D_Recursive_ok ( type_update, &
1572                                        tempP, tempC, &
1573                                        indmin, indmax,   &
1574                                       lb_child, ub_child,                    &
1575                                       s_child, s_parent, ds_child, ds_parent )
1576!---------------------------------------------------------------------------------------------------
1577    INTEGER, DIMENSION(2), intent(in)   :: type_update            !< Type of update (copy or average)
1578    INTEGER, DIMENSION(2), intent(in)   :: indmin, indmax
1579    INTEGER, DIMENSION(2), intent(in)   :: lb_child, ub_child
1580    REAL,    DIMENSION(2), intent(in)   :: s_child,  s_parent
1581    REAL,    DIMENSION(2), intent(in)   :: ds_child, ds_parent
1582    REAL,    DIMENSION(                 &
1583                indmin(1):indmax(1),    &
1584                indmin(2):indmax(2)),           intent(out) :: tempP
1585    REAL, DIMENSION(                            &
1586                lb_child(1):ub_child(1),  &
1587                lb_child(2):ub_child(2)), intent(in)  :: tempC
1588!
1589    REAL, DIMENSION(indmin(1):indmax(1), lb_child(2):ub_child(2)) :: tabtemp
1590    INTEGER :: i
1591!
1592    do i = lb_child(2),ub_child(2)
1593        call Agrif_Update_1D_Recursive(type_update(1),                              &
1594                                       tabtemp(:, i),          &
1595                                       tempC(:,i),  &
1596                                       indmin(1),indmax(1),                 &
1597                                       lb_child(1),ub_child(1),       &
1598                                       s_child(1), s_parent(1),             &
1599                                      ds_child(1),ds_parent(1))
1600    enddo
1601!
1602    tempP = 0.
1603!
1604    do i = indmin(1),indmax(1)
1605        call Agrif_UpdateBase(type_update(2),                                       &
1606                              tempP(i,:),             &
1607                              tabtemp(i,:), &
1608                              indmin(2),indmax(2),                          &
1609                              lb_child(2),ub_child(2),                &
1610                              s_parent(2),s_child(2),                       &
1611                             ds_parent(2),ds_child(2))
1612    enddo
1613!---------------------------------------------------------------------------------------------------
1614end subroutine Agrif_Update_2D_Recursive_ok
1615!===================================================================================================
1616
1617!
1618!===================================================================================================
1619!  subroutine Agrif_Update_3D_Recursive
1620!
1621!> Updates a 3D grid variable on the parent grid.
1622!! Calls #Agrif_Update_2D_Recursive and #Agrif_UpdateBase.
1623!---------------------------------------------------------------------------------------------------
1624subroutine Agrif_Update_3D_Recursive ( type_update,         &
1625                                       tempP, tempC,        &
1626                                       indmin, indmax,      &
1627                                       lb_child, ub_child,  &
1628                                        s_child,  s_parent, &
1629                                       ds_child, ds_parent )
1630!---------------------------------------------------------------------------------------------------
1631    integer, dimension(3),          intent(in)  :: type_update            !< Type of update (copy or average)
1632    integer, dimension(3),          intent(in)  :: indmin, indmax
1633    integer, dimension(3),          intent(in)  :: lb_child, ub_child
1634    real,    dimension(3),          intent(in)  ::  s_child,  s_parent
1635    real,    dimension(3),          intent(in)  :: ds_child, ds_parent
1636    real,    dimension(          &
1637        indmin(1):indmax(1),     &
1638        indmin(2):indmax(2),     &
1639        indmin(3):indmax(3)),       intent(out) :: tempP
1640    real, dimension(             &
1641        lb_child(1):ub_child(1), &
1642        lb_child(2):ub_child(2), &
1643        lb_child(3):ub_child(3)),   intent(in)  :: tempC
1644!---------------------------------------------------------------------------------------------------
1645    real, dimension(            &
1646        indmin(1):indmax(1),    &
1647        indmin(2):indmax(2),    &
1648        lb_child(3):ub_child(3))    :: tabtemp
1649    integer :: i,j,k
1650    integer :: coeffraf,locind_child_left
1651    integer :: kuinf
1652!
1653    coeffraf = nint ( ds_parent(1) / ds_child(1) )
1654!
1655    if ((type_update(1) == Agrif_Update_Average) .AND. (coeffraf /= 1 )) then
1656!---CDIR NEXPAND
1657        call Average1dPrecompute(ub_child(2)-lb_child(2)+1,&
1658                                 indmax(1)-indmin(1)+1,&
1659                                 ub_child(1)-lb_child(1)+1,&
1660                                 s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1661        precomputedone(1) = .TRUE.
1662    else if ((type_update(1) == Agrif_Update_Copy) .AND. (coeffraf /= 1 )) then
1663!---CDIR NEXPAND
1664         call Agrif_basicupdate_copy1d_before(ub_child(2)-lb_child(2)+1, &
1665                               indmax(1)-indmin(1)+1,           &
1666                               ub_child(1)-lb_child(1)+1, &
1667                               s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1668        precomputedone(1) = .TRUE.
1669    endif
1670!
1671    coeffraf = nint ( ds_parent(2) / ds_child(2) )
1672!
1673    if ((type_update(2) == Agrif_Update_Average) .AND. (coeffraf /= 1 )) then
1674!---CDIR NEXPAND
1675        call Average1dPrecompute(indmax(1)-indmin(1)+1,&
1676                                 indmax(2)-indmin(2)+1,&
1677                                 ub_child(2)-lb_child(2)+1,&
1678                                 s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1679        precomputedone(2) = .TRUE.
1680    else if ((type_update(2) == Agrif_Update_Copy) .AND. (coeffraf /= 1 )) then
1681!---CDIR NEXPAND
1682        call Agrif_basicupdate_copy1d_before( indmax(1)-indmin(1)+1,           &
1683                               indmax(2)-indmin(2)+1,           &
1684                               ub_child(2)-lb_child(2)+1, &
1685                               s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1686        precomputedone(2) = .TRUE.
1687    endif
1688!
1689    do k = lb_child(3),ub_child(3)
1690        call Agrif_Update_2D_Recursive( type_update(1:2),tabtemp(:,:,k),tempC(:,:,k), &
1691                                        indmin(1:2),indmax(1:2),                &
1692                                        lb_child(1:2),ub_child(1:2),      &
1693                                        s_child(1:2),s_parent(1:2),             &
1694                                        ds_child(1:2),ds_parent(1:2) )
1695    enddo
1696!
1697    precomputedone(1) = .FALSE.
1698    precomputedone(2) = .FALSE.
1699!
1700    coeffraf = nint ( ds_parent(3) / ds_child(3) )
1701    locind_child_left = 1 + agrif_int((s_parent(3)-s_child(3))/ds_child(3))
1702!
1703    if (coeffraf == 1) then
1704        kuinf = lb_child(3)+locind_child_left-2
1705        do k=indmin(3),indmax(3)
1706            kuinf = kuinf + 1
1707            do j = indmin(2),indmax(2)
1708            do i = indmin(1),indmax(1)
1709                tempP(i,j,k) = tabtemp(i,j,kuinf)
1710            enddo
1711            enddo
1712        enddo
1713    else
1714        tempP = 0.
1715        do j = indmin(2),indmax(2)
1716        do i = indmin(1),indmax(1)
1717            call Agrif_UpdateBase(type_update(3),tempP(i,j,:),tabtemp(i,j,:),   &
1718                                  indmin(3),indmax(3),                  &
1719                                  lb_child(3),ub_child(3),        &
1720                                  s_parent(3),s_child(3),               &
1721                                  ds_parent(3),ds_child(3))
1722!
1723        enddo
1724        enddo
1725    endif
1726!---------------------------------------------------------------------------------------------------
1727end subroutine Agrif_Update_3D_Recursive
1728!===================================================================================================
1729!
1730!===================================================================================================
1731!  subroutine Agrif_Update_4D_Recursive
1732!
1733!> Updates a 4D grid variable on the parent grid.
1734!! Calls #Agrif_Update_3D_Recursive and #Agrif_UpdateBase.
1735!---------------------------------------------------------------------------------------------------
1736subroutine Agrif_Update_4D_Recursive ( type_update,         &
1737                                       tempP, tempC,        &
1738                                       indmin, indmax,      &
1739                                       lb_child, ub_child,  &
1740                                        s_child,  s_parent, &
1741                                       ds_child, ds_parent )
1742!---------------------------------------------------------------------------------------------------
1743    integer, dimension(4),          intent(in)  :: type_update            !< Type of update (copy or average)
1744    integer, dimension(4),          intent(in)  :: indmin, indmax
1745    integer, dimension(4),          intent(in)  :: lb_child, ub_child
1746    real,    dimension(4),          intent(in)  ::  s_child,  s_parent
1747    real,    dimension(4),          intent(in)  :: ds_child, ds_parent
1748    real,    dimension(          &
1749        indmin(1):indmax(1),     &
1750        indmin(2):indmax(2),     &
1751        indmin(3):indmax(3),     &
1752        indmin(4):indmax(4)),       intent(out) :: tempP
1753    real, dimension(             &
1754        lb_child(1):ub_child(1), &
1755        lb_child(2):ub_child(2), &
1756        lb_child(3):ub_child(3), &
1757        lb_child(4):ub_child(4)),   intent(in)  :: tempC
1758!---------------------------------------------------------------------------------------------------
1759    real, dimension(:,:,:,:), allocatable       :: tabtemp
1760    integer :: i,j,k,l
1761!
1762    allocate(tabtemp(indmin(1):indmax(1), &
1763                     indmin(2):indmax(2), &
1764                     indmin(3):indmax(3), &
1765                     lb_child(4):ub_child(4)))
1766!
1767    do l = lb_child(4), ub_child(4)
1768        call Agrif_Update_3D_Recursive(type_update(1:3),                    &
1769                                       tabtemp(indmin(1):indmax(1),         &
1770                                               indmin(2):indmax(2),         &
1771                                               indmin(3):indmax(3), l),     &
1772                                       tempC(lb_child(1):ub_child(1),       &
1773                                             lb_child(2):ub_child(2),       &
1774                                             lb_child(3):ub_child(3), l),   &
1775                                       indmin(1:3), indmax(1:3),            &
1776                                       lb_child(1:3), ub_child(1:3),        &
1777                                        s_child(1:3),  s_parent(1:3),       &
1778                                       ds_child(1:3), ds_parent(1:3))
1779    enddo
1780!
1781    tempP = 0.
1782!
1783    do k = indmin(3), indmax(3)
1784    do j = indmin(2), indmax(2)
1785    do i = indmin(1), indmax(1)
1786        call Agrif_UpdateBase(type_update(4),                               &
1787                              tempP(i,j,k,indmin(4):indmax(4)),             &
1788                              tabtemp(i,j,k,lb_child(4):ub_child(4)),       &
1789                              indmin(4), indmax(4),                         &
1790                              lb_child(4), ub_child(4),                     &
1791                               s_parent(4), s_child(4),                     &
1792                              ds_parent(4),ds_child(4) )
1793    enddo
1794    enddo
1795    enddo
1796!
1797    deallocate(tabtemp)
1798!---------------------------------------------------------------------------------------------------
1799end subroutine Agrif_Update_4D_Recursive
1800!===================================================================================================
1801!
1802!===================================================================================================
1803!  subroutine Agrif_Update_5D_Recursive
1804!
1805!> Updates a 5D grid variable on the parent grid.
1806!! Calls #Agrif_Update_4D_Recursive and #Agrif_UpdateBase.
1807!---------------------------------------------------------------------------------------------------
1808subroutine Agrif_Update_5D_Recursive ( type_update,         &
1809                                       tempP, tempC,        &
1810                                       indmin, indmax,      &
1811                                       lb_child, ub_child,  &
1812                                        s_child,  s_parent, &
1813                                       ds_child, ds_parent )
1814!---------------------------------------------------------------------------------------------------
1815    integer, dimension(5),          intent(in)  :: type_update            !< Type of update (copy or average)
1816    integer, dimension(5),          intent(in)  :: indmin, indmax
1817    integer, dimension(5),          intent(in)  :: lb_child, ub_child
1818    real,    dimension(5),          intent(in)  ::  s_child,  s_parent
1819    real,    dimension(5),          intent(in)  :: ds_child, ds_parent
1820    real,    dimension(          &
1821        indmin(1):indmax(1),     &
1822        indmin(2):indmax(2),     &
1823        indmin(3):indmax(3),     &
1824        indmin(4):indmax(4),     &
1825        indmin(5):indmax(5)),       intent(out) :: tempP
1826    real, dimension(             &
1827        lb_child(1):ub_child(1), &
1828        lb_child(2):ub_child(2), &
1829        lb_child(3):ub_child(3), &
1830        lb_child(4):ub_child(4), &
1831        lb_child(5):ub_child(5)),   intent(in)  :: tempC
1832!---------------------------------------------------------------------------------------------------
1833    real, dimension(:,:,:,:,:), allocatable     :: tabtemp
1834    integer :: i,j,k,l,m
1835!
1836    allocate(tabtemp(indmin(1):indmax(1), &
1837                     indmin(2):indmax(2), &
1838                     indmin(3):indmax(3), &
1839                     indmin(4):indmax(4), &
1840                     lb_child(5):ub_child(5)))
1841!
1842    do m = lb_child(5), ub_child(5)
1843        call Agrif_Update_4D_Recursive(type_update(1:4),                    &
1844                                       tabtemp(indmin(1):indmax(1),         &
1845                                               indmin(2):indmax(2),         &
1846                                               indmin(3):indmax(3),         &
1847                                               indmin(4):indmax(4), m),     &
1848                                       tempC(lb_child(1):ub_child(1),       &
1849                                             lb_child(2):ub_child(2),       &
1850                                             lb_child(3):ub_child(3),       &
1851                                             lb_child(4):ub_child(4), m),   &
1852                                       indmin(1:4),indmax(1:4),             &
1853                                       lb_child(1:4), ub_child(1:4),        &
1854                                        s_child(1:4), s_parent(1:4),        &
1855                                       ds_child(1:4), ds_parent(1:4))
1856    enddo
1857!
1858    tempP = 0.
1859!
1860    do l = indmin(4), indmax(4)
1861    do k = indmin(3), indmax(3)
1862    do j = indmin(2), indmax(2)
1863    do i = indmin(1), indmax(1)
1864        call Agrif_UpdateBase( type_update(5),                              &
1865                               tempP(i,j,k,l,indmin(5):indmax(5)),          &
1866                               tabtemp(i,j,k,l,lb_child(5):ub_child(5)),    &
1867                               indmin(5), indmax(5),                        &
1868                               lb_child(5), ub_child(5),                    &
1869                                s_parent(5), s_child(5),                    &
1870                               ds_parent(5),ds_child(5) )
1871    enddo
1872    enddo
1873    enddo
1874    enddo
1875!
1876    deallocate(tabtemp)
1877!---------------------------------------------------------------------------------------------------
1878end subroutine Agrif_Update_5D_Recursive
1879!===================================================================================================
1880!
1881!===================================================================================================
1882!  subroutine Agrif_Update_6D_Recursive
1883!
1884!> Updates a 6D grid variable on the parent grid.
1885!! Calls #Agrif_Update_5D_Recursive and #Agrif_UpdateBase.
1886!---------------------------------------------------------------------------------------------------
1887subroutine Agrif_Update_6D_Recursive ( type_update,         &
1888                                       tempP, tempC,        &
1889                                       indmin, indmax,      &
1890                                       lb_child, ub_child,  &
1891                                        s_child,  s_parent, &
1892                                       ds_child, ds_parent )
1893!---------------------------------------------------------------------------------------------------
1894    integer, dimension(6),          intent(in)  :: type_update            !< Type of update (copy or average)
1895    integer, dimension(6),          intent(in)  :: indmin, indmax
1896    integer, dimension(6),          intent(in)  :: lb_child, ub_child
1897    real,    dimension(6),          intent(in)  ::  s_child,  s_parent
1898    real,    dimension(6),          intent(in)  :: ds_child, ds_parent
1899    real,    dimension(          &
1900        indmin(1):indmax(1),     &
1901        indmin(2):indmax(2),     &
1902        indmin(3):indmax(3),     &
1903        indmin(4):indmax(4),     &
1904        indmin(5):indmax(5),     &
1905        indmin(6):indmax(6)),       intent(out) :: tempP
1906    real, dimension(             &
1907        lb_child(1):ub_child(1), &
1908        lb_child(2):ub_child(2), &
1909        lb_child(3):ub_child(3), &
1910        lb_child(4):ub_child(4), &
1911        lb_child(5):ub_child(5), &
1912        lb_child(6):ub_child(6)),   intent(in)  :: tempC
1913!---------------------------------------------------------------------------------------------------
1914    real, dimension(:,:,:,:,:,:), allocatable   :: tabtemp
1915    integer :: i,j,k,l,m,n
1916!
1917    allocate(tabtemp(indmin(1):indmax(1), &
1918                     indmin(2):indmax(2), &
1919                     indmin(3):indmax(3), &
1920                     indmin(4):indmax(4), &
1921                     indmin(5):indmax(5), &
1922                     lb_child(6):ub_child(6)))
1923!
1924    do n = lb_child(6),ub_child(6)
1925        call Agrif_Update_5D_Recursive(type_update(1:5),                    &
1926                                       tabtemp(indmin(1):indmax(1),         &
1927                                               indmin(2):indmax(2),         &
1928                                               indmin(3):indmax(3),         &
1929                                               indmin(4):indmax(4),         &
1930                                               indmin(5):indmax(5), n),     &
1931                                       tempC(lb_child(1):ub_child(1),       &
1932                                             lb_child(2):ub_child(2),       &
1933                                             lb_child(3):ub_child(3),       &
1934                                             lb_child(4):ub_child(4),       &
1935                                             lb_child(5):ub_child(5), n),   &
1936                                       indmin(1:5), indmax(1:5),            &
1937                                       lb_child(1:5),ub_child(1:5),         &
1938                                       s_child(1:5), s_parent(1:5),         &
1939                                      ds_child(1:5),ds_parent(1:5))
1940    enddo
1941!
1942    tempP = 0.
1943!
1944    do m = indmin(5), indmax(5)
1945    do l = indmin(4), indmax(4)
1946    do k = indmin(3), indmax(3)
1947    do j = indmin(2), indmax(2)
1948    do i = indmin(1), indmax(1)
1949        call Agrif_UpdateBase( type_update(6),                              &
1950                               tempP(i,j,k,l,m,indmin(6):indmax(6)),        &
1951                               tabtemp(i,j,k,l,m,lb_child(6):ub_child(6)),  &
1952                               indmin(6), indmax(6),                        &
1953                               lb_child(6),  ub_child(6),                   &
1954                                s_parent(6),  s_child(6),                   &
1955                               ds_parent(6), ds_child(6) )
1956    enddo
1957    enddo
1958    enddo
1959    enddo
1960    enddo
1961!
1962    deallocate(tabtemp)
1963!---------------------------------------------------------------------------------------------------
1964end subroutine Agrif_Update_6D_Recursive
1965!===================================================================================================
1966!
1967!===================================================================================================
1968!  subroutine Agrif_UpdateBase
1969!
1970!> Calls the updating method chosen by the user (copy, average or full-weighting).
1971!---------------------------------------------------------------------------------------------------
1972subroutine Agrif_UpdateBase ( type_update,              &
1973                              parent_tab, child_tab,    &
1974                              indmin, indmax,           &
1975                              lb_child, ub_child,       &
1976                              s_parent, s_child,        &
1977                              ds_parent, ds_child )
1978!---------------------------------------------------------------------------------------------------
1979    integer,                            intent(in) :: type_update
1980    integer,                            intent(in) :: indmin, indmax
1981    integer,                            intent(in) :: lb_child, ub_child
1982    real, dimension(indmin:indmax),     intent(out):: parent_tab
1983    real, dimension(lb_child:ub_child), intent(in) :: child_tab
1984    real,                               intent(in) :: s_parent,  s_child
1985    real,                               intent(in) :: ds_parent, ds_child
1986!---------------------------------------------------------------------------------------------------
1987    integer :: np       ! Length of parent array
1988    integer :: nc       ! Length of child  array
1989!
1990    np = indmax - indmin + 1
1991    nc = ub_child - lb_child + 1
1992!
1993    if     ( type_update == Agrif_Update_Copy ) then
1994!
1995        call Agrif_basicupdate_copy1d(              &
1996                    parent_tab, child_tab,          &
1997                    np,         nc,                 &
1998                    s_parent,  s_child,             &
1999                    ds_parent, ds_child )
2000!
2001    elseif ( type_update == Agrif_Update_Average ) then
2002!
2003        call Agrif_basicupdate_average1d(           &
2004                    parent_tab, child_tab,          &
2005                    np,         nc,                 &
2006                     s_parent,  s_child,            &
2007                    ds_parent, ds_child )
2008!
2009    elseif ( type_update == Agrif_Update_Full_Weighting ) then
2010!
2011        call Agrif_basicupdate_full_weighting1D(    &
2012                    parent_tab, child_tab,          &
2013                    np,         nc,                 &
2014                    s_parent,  s_child,             &
2015                    ds_parent, ds_child )
2016!
2017    endif
2018!---------------------------------------------------------------------------------------------------
2019end subroutine Agrif_UpdateBase
2020!===================================================================================================
2021!
2022#if defined AGRIF_MPI
2023!===================================================================================================
2024!  subroutine Agrif_Find_list_update
2025!---------------------------------------------------------------------------------------------------
2026subroutine Agrif_Find_list_update ( list_update, pttab, petab, lb_child, lb_parent, nbdim, &
2027                                    find_list_update, tab4t, tab5t, memberinall, memberinall2,   &
2028                                    sendtoproc1, recvfromproc1, sendtoproc2, recvfromproc2 )
2029!---------------------------------------------------------------------------------------------------
2030    TYPE(Agrif_List_Interp_Loc),                   pointer     :: list_update
2031    INTEGER,                                       intent(in)  :: nbdim
2032    INTEGER, DIMENSION(nbdim),                     intent(in)  :: pttab, petab
2033    INTEGER, DIMENSION(nbdim),                     intent(in)  :: lb_child, lb_parent
2034    LOGICAL,                                       intent(out) :: find_list_update
2035    INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1,8), intent(out) :: tab4t
2036    INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1,8), intent(out) :: tab5t
2037    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(out) :: memberinall,memberinall2
2038    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(out) :: sendtoproc1,recvfromproc1
2039    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(out) :: sendtoproc2,recvfromproc2
2040!
2041    Type(Agrif_List_Interp_Loc), Pointer :: parcours
2042    INTEGER :: i
2043!
2044    find_list_update = .FALSE.
2045!
2046    parcours => list_update
2047
2048    Find_loop :  do while ( associated(parcours) )
2049        do i = 1,nbdim
2050            IF ((pttab(i) /= parcours%interp_loc%pttab(i)) .OR. &
2051                (petab(i) /= parcours%interp_loc%petab(i)) .OR. &
2052                (lb_child(i)  /= parcours%interp_loc%pttab_child(i)) .OR. &
2053                (lb_parent(i) /= parcours%interp_loc%pttab_parent(i))) THEN
2054                parcours => parcours%suiv
2055                cycle Find_loop
2056            ENDIF
2057        enddo
2058!
2059        tab4t = parcours%interp_loc%tab4t(1:nbdim,0:Agrif_Nbprocs-1,1:8)
2060        tab5t = parcours%interp_loc%tab5t(1:nbdim,0:Agrif_Nbprocs-1,1:8)
2061        memberinall  =  parcours%interp_loc%memberinall(0:Agrif_Nbprocs-1)
2062        memberinall2 =  parcours%interp_loc%memberinall2(0:Agrif_Nbprocs-1)
2063        sendtoproc1 =   parcours%interp_loc%sendtoproc1(0:Agrif_Nbprocs-1)
2064        sendtoproc2 =   parcours%interp_loc%sendtoproc2(0:Agrif_Nbprocs-1)
2065        recvfromproc1 = parcours%interp_loc%recvfromproc1(0:Agrif_Nbprocs-1)
2066        recvfromproc2 = parcours%interp_loc%recvfromproc2(0:Agrif_Nbprocs-1)
2067!
2068        find_list_update = .TRUE.
2069        exit Find_loop
2070!
2071    enddo Find_loop
2072!---------------------------------------------------------------------------------------------------
2073end subroutine Agrif_Find_list_update
2074!===================================================================================================
2075!
2076!===================================================================================================
2077!  subroutine Agrif_AddTo_list_update
2078!---------------------------------------------------------------------------------------------------
2079subroutine Agrif_AddTo_list_update ( list_update, pttab, petab, lb_child, lb_parent,  &
2080                                     nbdim, tab4t, tab5t, memberinall, memberinall2,        &
2081                                     sendtoproc1, recvfromproc1, sendtoproc2, recvfromproc2 )
2082!---------------------------------------------------------------------------------------------------
2083    TYPE(Agrif_List_Interp_Loc), pointer :: list_update
2084    INTEGER,                                        intent(in) :: nbdim
2085    INTEGER, DIMENSION(nbdim),                      intent(in) :: pttab, petab
2086    INTEGER, DIMENSION(nbdim),                      intent(in) :: lb_child, lb_parent
2087    INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1,8),  intent(in) :: tab4t
2088    INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1,8),  intent(in) :: tab5t
2089    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1),          intent(in) :: memberinall, memberinall2
2090    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1),          intent(in) :: sendtoproc1, recvfromproc1
2091    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1),          intent(in) :: sendtoproc2, recvfromproc2
2092!
2093    Type(Agrif_List_Interp_Loc), pointer :: parcours
2094!
2095    allocate(parcours)
2096    allocate(parcours%interp_loc)
2097
2098    parcours%interp_loc%pttab(1:nbdim) = pttab(1:nbdim)
2099    parcours%interp_loc%petab(1:nbdim) = petab(1:nbdim)
2100    parcours%interp_loc%pttab_child(1:nbdim)  = lb_child(1:nbdim)
2101    parcours%interp_loc%pttab_parent(1:nbdim) = lb_parent(1:nbdim)
2102
2103    allocate(parcours%interp_loc%tab4t(nbdim,0:Agrif_Nbprocs-1,8))
2104    allocate(parcours%interp_loc%tab5t(nbdim,0:Agrif_Nbprocs-1,8))
2105
2106    allocate(parcours%interp_loc%memberinall (0:Agrif_Nbprocs-1))
2107    allocate(parcours%interp_loc%memberinall2(0:Agrif_Nbprocs-1))
2108
2109    allocate(parcours%interp_loc%recvfromproc1(0:Agrif_Nbprocs-1))
2110    allocate(parcours%interp_loc%recvfromproc2(0:Agrif_Nbprocs-1))
2111    allocate(parcours%interp_loc%sendtoproc1(0:Agrif_Nbprocs-1))
2112    allocate(parcours%interp_loc%sendtoproc2(0:Agrif_Nbprocs-1))
2113
2114    parcours%interp_loc%tab4t = tab4t
2115    parcours%interp_loc%tab5t = tab5t
2116    parcours%interp_loc%memberinall   = memberinall
2117    parcours%interp_loc%memberinall2  = memberinall2
2118    parcours%interp_loc%sendtoproc1   = sendtoproc1
2119    parcours%interp_loc%sendtoproc2   = sendtoproc2
2120    parcours%interp_loc%recvfromproc1 = recvfromproc1
2121    parcours%interp_loc%recvfromproc2 = recvfromproc2
2122
2123    parcours%suiv => list_update
2124    list_update => parcours
2125!---------------------------------------------------------------------------------------------------
2126end subroutine Agrif_Addto_list_update
2127!===================================================================================================
2128#endif
2129!
2130end module Agrif_Update
Note: See TracBrowser for help on using the repository browser.