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.
modinterp.F90 in vendors/AGRIF/CMEMS_2020/AGRIF_FILES – NEMO

source: vendors/AGRIF/CMEMS_2020/AGRIF_FILES/modinterp.F90 @ 10872

Last change on this file since 10872 was 10725, checked in by rblod, 5 years ago

Update agrif library and conv see ticket #2129

  • Property svn:keywords set to Id
File size: 77.5 KB
Line 
1!
2! $Id$
3!
4!     AGRIF (Adaptive Grid Refinement In Fortran)
5!
6!     Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
7!                        Christophe Vouland (Christophe.Vouland@imag.fr)
8!
9!     This program is free software; you can redistribute it and/or modify
10!     it under the terms of the GNU General Public License as published by
11!     the Free Software Foundation; either version 2 of the License, or
12!     (at your option) any later version.
13!
14!     This program is distributed in the hope that it will be useful,
15!     but WITHOUT ANY WARRANTY; without even the implied warranty of
16!     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
17!     GNU General Public License for more details.
18!
19!     You should have received a copy of the GNU General Public License
20!     along with this program; if not, write to the Free Software
21!     Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA.
22!
23!
24!> Module to initialize a fine grid from its parent grid, by using a space interpolation
25!
26module Agrif_Interpolation
27!
28    use Agrif_InterpBasic
29    use Agrif_Arrays
30    use Agrif_Mask
31    use Agrif_CurgridFunctions
32#if defined AGRIF_MPI
33    use Agrif_Mpp
34#endif
35!
36    implicit none
37!
38    logical, private:: precomputedone(7) = .FALSE.
39!
40    private :: Agrif_Parentbounds
41    private :: Agrif_Interp_1D_recursive, Agrif_Interp_2D_recursive, Agrif_Interp_3D_recursive
42    private :: Agrif_Interp_4D_recursive, Agrif_Interp_5D_recursive, Agrif_Interp_6D_recursive
43    private :: Agrif_InterpBase
44    private :: Agrif_Find_list_interp, Agrif_AddTo_list_interp
45!
46contains
47!
48!===================================================================================================
49!  subroutine Agrif_InterpVariable
50!
51!> Sets some arguments of subroutine Agrif_InterpnD, n being the dimension of the grid variable
52!---------------------------------------------------------------------------------------------------
53subroutine Agrif_InterpVariable ( parent, child, torestore, procname )
54!---------------------------------------------------------------------------------------------------
55    type(Agrif_Variable), pointer       :: parent       !< Variable on the parent grid
56    type(Agrif_Variable), pointer       :: child        !< Variable on the child grid
57    logical,              intent(in)    :: torestore    !< .false. indicates that the results of the
58                                                        !! interpolation are applied on the whole current grid
59    procedure()                         :: procname     !< Data recovery procedure
60!---------------------------------------------------------------------------------------------------
61    logical               :: memberin
62    integer               :: nbdim        ! Number of dimensions of the current grid
63    integer, dimension(6) :: type_interp  ! Type of interpolation (linear,spline,...)
64    integer, dimension(6) :: nb_child
65    integer, dimension(6) :: lb_child
66    integer, dimension(6) :: ub_child
67    integer, dimension(6) :: lb_parent
68    real   , dimension(6) :: s_child,   s_parent
69    real   , dimension(6) :: ds_child, ds_parent
70    integer, dimension(child % root_var % nbdim,2,2) :: childarray
71!
72    nbdim       = child % root_var % nbdim
73    type_interp = child % root_var % type_interp
74!
75    call PreProcessToInterpOrUpdate( parent,   child,       &
76                                     nb_child, ub_child,    &
77                                     lb_child, lb_parent,   &
78                                      s_child,  s_parent,   &
79                                     ds_child, ds_parent, nbdim, interp=.true.)
80!
81!   Call to a procedure of interpolation against the number of dimensions of the grid variable
82!
83    call Agrif_InterpnD(type_interp, parent, child, &
84                        lb_child, ub_child,         &
85                        lb_child, lb_parent,        &
86                        s_child,   s_parent,        &
87                        ds_child, ds_parent,        &
88                        child, torestore, nbdim,    &
89                        childarray, memberin,       &
90                        .false., procname, 0, 0)
91!---------------------------------------------------------------------------------------------------
92end subroutine Agrif_InterpVariable
93!===================================================================================================
94!
95!===================================================================================================
96!  subroutine Agrif_InterpnD
97!
98!> Interpolates a nD grid variable from its parent grid, by using a space interpolation
99!---------------------------------------------------------------------------------------------------
100subroutine Agrif_InterpnD ( type_interp, parent, child, pttab, petab, pttab_Child, pttab_Parent, &
101                            s_Child, s_Parent, ds_Child, ds_Parent, restore, torestore,          &
102                            nbdim, childarray, memberin, in_bc, procname, nb, ndir )
103!---------------------------------------------------------------------------------------------------
104#if defined AGRIF_MPI
105    include 'mpif.h'
106#endif
107!
108    INTEGER, DIMENSION(6),     INTENT(in)   :: type_interp  !< Type of interpolation !    (linear,...)
109    TYPE(Agrif_Variable),      pointer      :: parent       !< Variable of the parent grid
110    TYPE(Agrif_Variable),      pointer      :: child        !< Variable of the child grid
111    INTEGER, DIMENSION(nbdim), INTENT(in)   :: pttab        !< Index of the first point inside the domain
112    INTEGER, DIMENSION(nbdim), INTENT(in)   :: petab        !< Index of the first point inside the domain
113    INTEGER, DIMENSION(nbdim), INTENT(in)   :: pttab_Child  !< Index of the first point inside the domain
114                                                            !<    for the child grid variable
115    INTEGER, DIMENSION(nbdim), INTENT(in)   :: pttab_Parent !< Index of the first point inside the domain
116                                                            !<    for the parent grid variable
117    REAL,    DIMENSION(nbdim), INTENT(in)   :: s_Child,s_Parent   !< Positions of the parent and child grids
118    REAL,    DIMENSION(nbdim), INTENT(in)   :: ds_Child,ds_Parent !< Space steps of the parent and child grids
119    TYPE(Agrif_Variable),      pointer      :: restore            !< Indicates points where interpolation
120    LOGICAL,                   INTENT(in)   :: torestore          !< Indicates if the array restore is used
121    INTEGER,                   INTENT(in)   :: nbdim
122    LOGICAL,                   INTENT(out)  :: memberin
123    LOGICAL,                   INTENT(in)   :: in_bc              !< .true. if called from Agrif_CorrectVariable \n
124                                                                  !! .false. if called from Agrif_InterpVariable
125    procedure()                             :: procname           !< Data recovery procedure
126    INTEGER,                   INTENT(in)   :: nb, ndir
127!
128    INTEGER                       :: i,j,k,l,m,n
129    INTEGER, DIMENSION(nbdim)     :: pttruetab,cetruetab
130    INTEGER, DIMENSION(nbdim)     :: indmin,     indmax
131    INTEGER, DIMENSION(nbdim)     :: indminglob, indmaxglob
132#if defined AGRIF_MPI
133    INTEGER, DIMENSION(nbdim)     :: indminglob2,indmaxglob2
134    INTEGER, DIMENSION(nbdim)     :: indminglob3,indmaxglob3
135#endif
136    LOGICAL, DIMENSION(nbdim)     :: noraftab
137    REAL   , DIMENSION(nbdim)     :: s_Child_temp,s_Parent_temp
138    INTEGER, DIMENSION(nbdim)     :: lowerbound, upperbound, coords
139    INTEGER, DIMENSION(nbdim,2,2), INTENT(OUT) :: childarray
140    INTEGER, DIMENSION(nbdim,2,2)              :: parentarray
141    LOGICAL :: member
142    LOGICAL :: find_list_interp
143!
144#if defined AGRIF_MPI
145!
146    INTEGER, PARAMETER          :: etiquette = 100
147    INTEGER                     :: code, local_proc
148    INTEGER, DIMENSION(nbdim,4)                     :: tab3
149    INTEGER, DIMENSION(nbdim,4,0:Agrif_Nbprocs-1)   :: tab4
150    INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1,8)   :: tab4t
151    INTEGER,DIMENSION(nbdim,2) :: tab5
152    INTEGER,DIMENSION(nbdim,2,0:Agrif_Nbprocs-1) :: tab6
153    INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1,2) :: tab5t
154    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1)           :: memberinall
155    LOGICAL, DIMENSION(0:Agrif_Nbprocs-1)           :: sendtoproc1, recvfromproc1
156    LOGICAL, DIMENSION(1)                           :: memberin1
157    LOGICAL                                         :: memberout
158!
159#endif
160!
161    type(Agrif_Variable), pointer, save :: tempC => NULL()        ! Temporary child grid variable
162    type(Agrif_Variable), pointer, save :: tempP => NULL()        ! Temporary parent grid variable
163    type(Agrif_Variable), pointer, save :: tempPextend => NULL()  ! Temporary parent grid variable
164    type(Agrif_Variable), pointer, save :: parentvalues => NULL()
165!
166    coords = child % root_var % coords
167!
168!   Boundaries of the current grid where interpolation is done
169    find_list_interp =                                              &
170        Agrif_Find_list_interp(                                     &
171            child % list_interp,                                    &
172            pttab, petab, pttab_Child, pttab_Parent, nbdim,         &
173            indmin, indmax, indminglob, indmaxglob,                 &
174            pttruetab, cetruetab, memberin                          &
175#if defined AGRIF_MPI
176           ,indminglob2, indmaxglob2, parentarray,                  &
177            member, tab4t,memberinall,  sendtoproc1, recvfromproc1  &
178#endif
179        )
180!
181    if (.not.find_list_interp) then
182!
183        call Agrif_get_var_bounds_array(child, lowerbound, upperbound, nbdim)
184        call Agrif_Childbounds(nbdim, lowerbound, upperbound,               &
185                               pttab, petab, Agrif_Procrank, coords,        &
186                               pttruetab, cetruetab, memberin)
187        call Agrif_Parentbounds(type_interp,nbdim,indminglob,indmaxglob,    &
188                                s_Parent_temp,s_Child_temp,                 &
189                                s_Child,ds_Child,                           &
190                                s_Parent,ds_Parent,                         &
191                                pttab,petab,                                &
192                                pttab_Child,pttab_Parent,                   &
193                                child%root_var % posvar, coords)
194#if defined AGRIF_MPI
195        if (memberin) then
196            call Agrif_Parentbounds(type_interp,nbdim,indmin,indmax,        &
197                                    s_Parent_temp,s_Child_temp,             &
198                                    s_Child,ds_Child,                       &
199                                    s_Parent,ds_Parent,                     &
200                                    pttruetab,cetruetab,                    &
201                                    pttab_Child,pttab_Parent,               &
202                                    child%root_var % posvar, coords)
203        endif
204
205        local_proc = Agrif_Procrank
206        call Agrif_get_var_bounds_array(parent,lowerbound,upperbound,nbdim)
207        call Agrif_ChildGrid_to_ParentGrid()
208!
209        call Agrif_Childbounds(nbdim,lowerbound,upperbound,                 &
210                               indminglob,indmaxglob, local_proc, coords,   &
211                               indminglob2,indmaxglob2,member,              &
212                               indminglob3,indmaxglob3)
213!
214        if (member) then
215            call Agrif_GlobalToLocalBounds(parentarray,                     &
216                                           lowerbound,  upperbound,         &
217                                           indminglob2, indmaxglob2, coords,&
218                                           nbdim, local_proc, member)
219        endif
220
221        call Agrif_ParentGrid_to_ChildGrid()
222#else
223        parentarray(:,1,1) = indminglob
224        parentarray(:,2,1) = indmaxglob
225        parentarray(:,1,2) = indminglob
226        parentarray(:,2,2) = indmaxglob
227        indmin = indminglob
228        indmax = indmaxglob
229        member = .TRUE.
230#endif
231!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
232! Correct for non refined directions
233!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
234        do i=1,nbdim
235          if (coords(i) == 0) then
236             indmin(i) = indminglob(i)
237             indmax(i) = indmaxglob(i)
238             pttruetab(i) = indminglob(i)
239             cetruetab(i) = indmaxglob(i)
240          endif
241        enddo
242
243    else
244
245#if defined AGRIF_MPI
246        s_Parent_temp = s_Parent + (indmin - pttab_Parent) * ds_Parent
247        s_Child_temp  = s_Child + (pttruetab - pttab_Child) * ds_Child
248#else
249        parentarray(:,1,1) = indminglob
250        parentarray(:,2,1) = indmaxglob
251        parentarray(:,1,2) = indminglob
252        parentarray(:,2,2) = indmaxglob
253        indmin = indminglob
254        indmax = indmaxglob
255        member = .TRUE.
256        s_Parent_temp = s_Parent + (indminglob - pttab_Parent) * ds_Parent
257        s_Child_temp  = s_Child + (pttab - pttab_Child) * ds_Child
258#endif
259    endif
260!
261    if (member) then
262        if (.not.associated(tempP)) allocate(tempP)
263!
264        call Agrif_array_allocate(tempP,parentarray(:,1,1),parentarray(:,2,1),nbdim)
265        call Agrif_var_set_array_tozero(tempP,nbdim)
266
267        call Agrif_ChildGrid_to_ParentGrid()
268!
269        select case (nbdim)
270        case(1)
271            call procname(tempP%array1,                         &
272                      parentarray(1,1,2),parentarray(1,2,2),.TRUE.,nb,ndir)
273        case(2)
274            call procname(tempP%array2,                         &
275                      parentarray(1,1,2),parentarray(1,2,2),    &
276                      parentarray(2,1,2),parentarray(2,2,2),.TRUE.,nb,ndir)
277        case(3)
278            call procname(tempP%array3,                         &
279                      parentarray(1,1,2),parentarray(1,2,2),    &
280                      parentarray(2,1,2),parentarray(2,2,2),    &
281                      parentarray(3,1,2),parentarray(3,2,2),.TRUE.,nb,ndir)
282        case(4)
283            call procname(tempP%array4,                         &
284                      parentarray(1,1,2),parentarray(1,2,2),    &
285                      parentarray(2,1,2),parentarray(2,2,2),    &
286                      parentarray(3,1,2),parentarray(3,2,2),    &
287                      parentarray(4,1,2),parentarray(4,2,2),.TRUE.,nb,ndir)
288        case(5)
289            call procname(tempP%array5,                         &
290                      parentarray(1,1,2),parentarray(1,2,2),    &
291                      parentarray(2,1,2),parentarray(2,2,2),    &
292                      parentarray(3,1,2),parentarray(3,2,2),    &
293                      parentarray(4,1,2),parentarray(4,2,2),    &
294                      parentarray(5,1,2),parentarray(5,2,2),.TRUE.,nb,ndir)
295        case(6)
296            call procname(tempP%array6,                         &
297                      parentarray(1,1,2),parentarray(1,2,2),    &
298                      parentarray(2,1,2),parentarray(2,2,2),    &
299                      parentarray(3,1,2),parentarray(3,2,2),    &
300                      parentarray(4,1,2),parentarray(4,2,2),    &
301                      parentarray(5,1,2),parentarray(5,2,2),    &
302                      parentarray(6,1,2),parentarray(6,2,2),.TRUE.,nb,ndir)
303        end select
304!
305        call Agrif_ParentGrid_to_ChildGrid()
306!
307    endif
308
309#if defined AGRIF_MPI
310    if (.not.find_list_interp) then
311!
312        tab3(:,1) = indminglob2(:)
313        tab3(:,2) = indmaxglob2(:)
314        tab3(:,3) = indmin(:)
315        tab3(:,4) = indmax(:)
316        tab5(:,1) = indminglob3(:)
317        tab5(:,2) = indmaxglob3(:)
318!
319        call MPI_ALLGATHER(tab3,4*nbdim,MPI_INTEGER,tab4,4*nbdim,MPI_INTEGER,Agrif_mpi_comm,code)
320        call MPI_ALLGATHER(tab5,2*nbdim,MPI_INTEGER,tab6,2*nbdim,MPI_INTEGER,Agrif_mpi_comm,code)
321        if (.not.associated(tempPextend))   allocate(tempPextend)
322
323        do k=0,Agrif_Nbprocs-1
324            do j=1,4
325                do i=1,nbdim
326                    tab4t(i,k,j) = tab4(i,j,k)
327                enddo
328            enddo
329        enddo
330
331        do k=0,Agrif_Nbprocs-1
332          do j=1,2
333            do i=1,nbdim
334               tab5t(i,k,j) = tab6(i,j,k)
335            enddo
336          enddo
337        enddo
338     
339        memberin1(1) = memberin
340        call MPI_ALLGATHER(memberin1,1,MPI_LOGICAL,memberinall,1,MPI_LOGICAL,Agrif_mpi_comm,code)
341
342        call Get_External_Data_first(tab4t(:,:,1),tab4t(:,:,2),         &
343                                     tab4t(:,:,3),tab4t(:,:,4),         &
344                                     nbdim,memberinall, coords,         &
345                                     sendtoproc1,recvfromproc1,         &
346                                     tab4t(:,:,5),tab4t(:,:,6),         &
347                                     tab4t(:,:,7),tab4t(:,:,8),         &
348                                     tab5t(:,:,1),tab5t(:,:,2))
349    endif
350
351    call ExchangeSameLevel(sendtoproc1,recvfromproc1,nbdim,         &
352            tab4t(:,:,3),tab4t(:,:,4),tab4t(:,:,5),tab4t(:,:,6),    &
353            tab4t(:,:,7),tab4t(:,:,8),memberin,tempP,tempPextend)
354#else
355    tempPextend => tempP
356#endif
357
358    if (.not.find_list_interp) then
359        call Agrif_Addto_list_interp(                           &
360                child%list_interp,pttab,petab,                  &
361                pttab_Child,pttab_Parent,indmin,indmax,         &
362                indminglob,indmaxglob,                          &
363                pttruetab,cetruetab,                            &
364                memberin,nbdim                                  &
365#if defined AGRIF_MPI
366               ,indminglob2,indmaxglob2,                        &
367                parentarray,                                    &
368                member,                                         &
369                tab4t,memberinall,sendtoproc1,recvfromproc1     &
370#endif
371        )
372    endif
373!
374    if (memberin) then
375!
376        if (.not.associated(tempC)) allocate(tempC)
377!
378        call Agrif_array_allocate(tempC,pttruetab,cetruetab,nbdim)
379!
380!       Special values on the parent grid
381        if (Agrif_UseSpecialValue) then
382!
383            noraftab(1:nbdim) = child % root_var % interptab(1:nbdim) == 'N'
384!
385            if (.not.associated(parentvalues))  allocate(parentvalues)
386!
387            call Agrif_array_allocate(parentvalues,indmin,indmax,nbdim)
388            call Agrif_var_full_copy_array(parentvalues,tempPextend,nbdim)
389!
390            call Agrif_CheckMasknD(tempPextend,parentvalues,    &
391                    indmin(1:nbdim),indmax(1:nbdim),            &
392                    indmin(1:nbdim),indmax(1:nbdim),            &
393                    noraftab(1:nbdim),nbdim)
394!
395            call Agrif_array_deallocate(parentvalues,nbdim)
396!
397        endif
398!
399!       Interpolation of the current grid
400!
401        if ( memberin ) then
402            select case(nbdim)
403            case(1)
404                call Agrif_Interp_1D_recursive( type_interp(1),                         &
405                                                tempPextend%array1,                     &
406                                                tempC%array1,                           &
407                                                indmin(1), indmax(1),                   &
408                                                pttruetab(1),    cetruetab(1),          &
409                                                s_Child_temp(1), s_Parent_temp(1),      &
410                                                ds_Child(1),     ds_Parent(1) )
411            case(2)
412                call Agrif_Interp_2D_recursive( type_interp(1:2),                       &
413                                                tempPextend % array2,                   &
414                                                tempC       % array2,                   &
415                                                indmin(1:2), indmax(1:2),               &
416                                                pttruetab(1:2),    cetruetab(1:2),      &
417                                                s_Child_temp(1:2), s_Parent_temp(1:2),  &
418                                                ds_Child(1:2),    ds_Parent(1:2) )
419            case(3)
420                call Agrif_Interp_3D_recursive( type_interp(1:3),                       &
421                                                tempPextend % array3,                   &
422                                                tempC       % array3,                   &
423                                                indmin(1:3), indmax(1:3),               &
424                                                pttruetab(1:3),    cetruetab(1:3),      &
425                                                s_Child_temp(1:3), s_Parent_temp(1:3),  &
426                                                ds_Child(1:3),    ds_Parent(1:3) )
427            case(4)
428                call Agrif_Interp_4D_recursive( type_interp(1:4),                       &
429                                                tempPextend % array4,                   &
430                                                tempC       % array4,                   &
431                                                indmin(1:4), indmax(1:4),               &
432                                                pttruetab(1:4),    cetruetab(1:4),      &
433                                                s_Child_temp(1:4), s_Parent_temp(1:4),  &
434                                                ds_Child(1:4),    ds_Parent(1:4) )
435            case(5)
436                call Agrif_Interp_5D_recursive( type_interp(1:5),                       &
437                                                tempPextend % array5,                   &
438                                                tempC       % array5,                   &
439                                                indmin(1:5), indmax(1:5),               &
440                                                pttruetab(1:5),    cetruetab(1:5),      &
441                                                s_Child_temp(1:5), s_Parent_temp(1:5),  &
442                                                ds_Child(1:5),    ds_Parent(1:5) )
443            case(6)
444                call Agrif_Interp_6D_recursive( type_interp(1:6),                       &
445                                                tempPextend % array6,                   &
446                                                tempC       % array6,                   &
447                                                indmin(1:6), indmax(1:6),               &
448                                                pttruetab(1:6),    cetruetab(1:6),      &
449                                                s_Child_temp(1:6), s_Parent_temp(1:6),  &
450                                                ds_Child(1:6),    ds_Parent(1:6) )
451            end select
452!
453            call Agrif_get_var_bounds_array(child,lowerbound,upperbound,nbdim)
454
455#if defined AGRIF_MPI
456            call Agrif_GlobalToLocalBounds(childarray, lowerbound, upperbound,  &
457                                            pttruetab, cetruetab, coords,       &
458                                            nbdim, Agrif_Procrank, memberout)
459#else
460            childarray(:,1,1) = pttruetab
461            childarray(:,2,1) = cetruetab
462            childarray(:,1,2) = pttruetab
463            childarray(:,2,2) = cetruetab
464!cccccccccccccc       memberout = .TRUE.
465#endif
466!
467!           Special values on the child grid
468            if (Agrif_UseSpecialValueFineGrid) then
469                call GiveAgrif_SpecialValueToTab_mpi( child, tempC, childarray, Agrif_SpecialValueFineGrid,nbdim )
470            endif
471!
472        endif   ! ( memberin )
473!
474        if (torestore) then
475!
476#if defined AGRIF_MPI
477!
478            SELECT CASE (nbdim)
479            CASE (1)
480                do i = pttruetab(1),cetruetab(1)
481!hildarrayAModifier     if (restore%restore1D(i) == 0)      &
482!hildarrayAModifier         child%array1(childarray(i,1,2)) = tempC%array1(i)
483                enddo
484            CASE (2)
485                do i = pttruetab(1),cetruetab(1)
486                do j = pttruetab(2),cetruetab(2)
487!hildarrayAModifier     if (restore%restore2D(i,j) == 0)    &
488!hildarrayAModifier         child%array2(childarray(i,1,2), &
489!hildarrayAModifier                          childarray(j,2,2)) = tempC%array2(i,j)
490                enddo
491                enddo
492            CASE (3)
493                do i = pttruetab(1),cetruetab(1)
494                do j = pttruetab(2),cetruetab(2)
495                do k = pttruetab(3),cetruetab(3)
496!hildarrayAModifier     if (restore%restore3D(i,j,k) == 0)  &
497!hildarrayAModifier         child%array3(childarray(i,1,2), &
498!hildarrayAModifier                          childarray(j,2,2), &
499!hildarrayAModifier                          childarray(k,3,2)) = tempC%array3(i,j,k)
500                enddo
501                enddo
502                enddo
503            CASE (4)
504                do i = pttruetab(1),cetruetab(1)
505                do j = pttruetab(2),cetruetab(2)
506                do k = pttruetab(3),cetruetab(3)
507                do l = pttruetab(4),cetruetab(4)
508!hildarrayAModifier     if (restore%restore4D(i,j,k,l) == 0)    &
509!hildarrayAModifier         child%array4(childarray(i,1,2),     &
510!hildarrayAModifier                          childarray(j,2,2),     &
511!hildarrayAModifier                          childarray(k,3,2),     &
512!hildarrayAModifier                          childarray(l,4,2)) = tempC%array4(i,j,k,l)
513                enddo
514                enddo
515                enddo
516                enddo
517            CASE (5)
518                do i = pttruetab(1),cetruetab(1)
519                do j = pttruetab(2),cetruetab(2)
520                do k = pttruetab(3),cetruetab(3)
521                do l = pttruetab(4),cetruetab(4)
522                do m = pttruetab(5),cetruetab(5)
523!hildarrayAModifier     if (restore%restore5D(i,j,k,l,m) == 0)  &
524!hildarrayAModifier         child%array5(childarray(i,1,2),     &
525!hildarrayAModifier                          childarray(j,2,2),     &
526!hildarrayAModifier                          childarray(k,3,2),     &
527!hildarrayAModifier                          childarray(l,4,2),     &
528!hildarrayAModifier                          childarray(m,5,2)) = tempC%array5(i,j,k,l,m)
529                enddo
530                enddo
531                enddo
532                enddo
533                enddo
534            CASE (6)
535                do i = pttruetab(1),cetruetab(1)
536                do j = pttruetab(2),cetruetab(2)
537                do k = pttruetab(3),cetruetab(3)
538                do l = pttruetab(4),cetruetab(4)
539                do m = pttruetab(5),cetruetab(5)
540                do n = pttruetab(6),cetruetab(6)
541!hildarrayAModifier     if (restore%restore6D(i,j,k,l,m,n) == 0)    &
542!hildarrayAModifier         child%array6(childarray(i,1,2),         &
543!hildarrayAModifier                          childarray(j,2,2),         &
544!hildarrayAModifier                          childarray(k,3,2),         &
545!hildarrayAModifier                          childarray(l,4,2),         &
546!hildarrayAModifier                          childarray(m,5,2),         &
547!hildarrayAModifier                          childarray(n,6,2)) = tempC%array6(i,j,k,l,m,n)
548                enddo
549                enddo
550                enddo
551                enddo
552                enddo
553                enddo
554            END SELECT
555!
556#else
557            select case (nbdim)
558            case (1)
559                do i = pttruetab(1),cetruetab(1)
560                    if (restore%restore1D(i) == 0)          &
561                        parray1(i) = tempC % array1(i)
562                enddo
563            case (2)
564                do j = pttruetab(2),cetruetab(2)
565                do i = pttruetab(1),cetruetab(1)
566                    if (restore%restore2D(i,j) == 0)        &
567                        parray2(i,j) = tempC % array2(i,j)
568                enddo
569                enddo
570            case (3)
571                do k = pttruetab(3),cetruetab(3)
572                do j = pttruetab(2),cetruetab(2)
573                do i = pttruetab(1),cetruetab(1)
574                    if (restore%restore3D(i,j,k) == 0)      &
575                        parray3(i,j,k) = tempC % array3(i,j,k)
576                enddo
577                enddo
578                enddo
579            case (4)
580                do l = pttruetab(4),cetruetab(4)
581                do k = pttruetab(3),cetruetab(3)
582                do j = pttruetab(2),cetruetab(2)
583                do i = pttruetab(1),cetruetab(1)
584                    if (restore%restore4D(i,j,k,l) == 0)    &
585                        parray4(i,j,k,l) = tempC % array4(i,j,k,l)
586                enddo
587                enddo
588                enddo
589                enddo
590            case (5)
591                do m = pttruetab(5),cetruetab(5)
592                do l = pttruetab(4),cetruetab(4)
593                do k = pttruetab(3),cetruetab(3)
594                do j = pttruetab(2),cetruetab(2)
595                do i = pttruetab(1),cetruetab(1)
596                    if (restore%restore5D(i,j,k,l,m) == 0)  &
597                        parray5(i,j,k,l,m) = tempC % array5(i,j,k,l,m)
598                enddo
599                enddo
600                enddo
601                enddo
602                enddo
603            case (6)
604                do n = pttruetab(6),cetruetab(6)
605                do m = pttruetab(5),cetruetab(5)
606                do l = pttruetab(4),cetruetab(4)
607                do k = pttruetab(3),cetruetab(3)
608                do j = pttruetab(2),cetruetab(2)
609                do i = pttruetab(1),cetruetab(1)
610                    if (restore%restore6D(i,j,k,l,m,n) == 0)    &
611                        parray6(i,j,k,l,m,n) = tempC % array6(i,j,k,l,m,n)
612                enddo
613                enddo
614                enddo
615                enddo
616                enddo
617                enddo
618            end select
619!
620#endif
621!
622        else    ! .not.to_restore
623!
624            if (memberin) then
625    !
626                if ( .not.in_bc ) then
627                    select case(nbdim)
628                    case(1)
629                        call procname(tempC % array1(            &
630                                childarray(1,1,1):childarray(1,2,1)), &
631                                childarray(1,1,2),childarray(1,2,2),.FALSE.,nb,ndir)
632                    case(2)
633                        call procname(                                &
634                                tempC % array2(            &
635                                childarray(1,1,1):childarray(1,2,1),  &
636                                childarray(2,1,1):childarray(2,2,1)), &
637                                childarray(1,1,2),childarray(1,2,2),  &
638                                childarray(2,1,2),childarray(2,2,2),.FALSE.,nb,ndir)
639                    case(3)
640                        call procname(                                &
641                                tempC % array3(            &
642                                childarray(1,1,1):childarray(1,2,1),  &
643                                childarray(2,1,1):childarray(2,2,1),  &
644                                childarray(3,1,1):childarray(3,2,1)), &
645                                childarray(1,1,2),childarray(1,2,2),  &
646                                childarray(2,1,2),childarray(2,2,2),  &
647                                childarray(3,1,2),childarray(3,2,2),.FALSE.,nb,ndir)
648                    case(4)
649                        call procname(                                &
650                                tempC % array4(            &
651                                childarray(1,1,1):childarray(1,2,1),  &
652                                childarray(2,1,1):childarray(2,2,1),  &
653                                childarray(3,1,1):childarray(3,2,1),  &
654                                childarray(4,1,1):childarray(4,2,1)), &
655                                childarray(1,1,2),childarray(1,2,2),  &
656                                childarray(2,1,2),childarray(2,2,2),  &
657                                childarray(3,1,2),childarray(3,2,2),  &
658                                childarray(4,1,2),childarray(4,2,2),.FALSE.,nb,ndir)
659                    case(5)
660                        call procname(                                &
661                                tempC % array5(            &
662                                childarray(1,1,1):childarray(1,2,1),  &
663                                childarray(2,1,1):childarray(2,2,1),  &
664                                childarray(3,1,1):childarray(3,2,1),  &
665                                childarray(4,1,1):childarray(4,2,1),  &
666                                childarray(5,1,1):childarray(5,2,1)), &
667                                childarray(1,1,2),childarray(1,2,2),  &
668                                childarray(2,1,2),childarray(2,2,2),  &
669                                childarray(3,1,2),childarray(3,2,2),  &
670                                childarray(4,1,2),childarray(4,2,2),  &
671                                childarray(5,1,2),childarray(5,2,2),.FALSE.,nb,ndir)
672                    case(6)
673                        call procname(                                &
674                                tempC % array6(            &
675                                childarray(1,1,1):childarray(1,2,1),  &
676                                childarray(2,1,1):childarray(2,2,1),  &
677                                childarray(3,1,1):childarray(3,2,1),  &
678                                childarray(4,1,1):childarray(4,2,1),  &
679                                childarray(5,1,1):childarray(5,2,1),  &
680                                childarray(6,1,1):childarray(6,2,1)), &
681                                childarray(1,1,2),childarray(1,2,2),  &
682                                childarray(2,1,2),childarray(2,2,2),  &
683                                childarray(3,1,2),childarray(3,2,2),  &
684                                childarray(4,1,2),childarray(4,2,2),  &
685                                childarray(5,1,2),childarray(5,2,2),  &
686                                childarray(6,1,2),childarray(6,2,2),.FALSE.,nb,ndir)
687                    end select
688                else    ! we are in_bc
689                    select case (nbdim)
690                    case (1)
691                        parray1(childarray(1,1,2):childarray(1,2,2)) =    &
692                         tempC%array1(childarray(1,1,1):childarray(1,2,1))
693                    case (2)
694                        parray2(childarray(1,1,2):childarray(1,2,2),      &
695                                      childarray(2,1,2):childarray(2,2,2)) =    &
696                         tempC%array2(childarray(1,1,1):childarray(1,2,1),      &
697                                      childarray(2,1,1):childarray(2,2,1))
698                    case (3)
699                        parray3(childarray(1,1,2):childarray(1,2,2),      &
700                                      childarray(2,1,2):childarray(2,2,2),      &
701                                      childarray(3,1,2):childarray(3,2,2)) =    &
702                         tempC%array3(childarray(1,1,1):childarray(1,2,1),      &
703                                      childarray(2,1,1):childarray(2,2,1),      &
704                                      childarray(3,1,1):childarray(3,2,1))
705                    case (4)
706                        parray4(childarray(1,1,2):childarray(1,2,2),      &
707                                      childarray(2,1,2):childarray(2,2,2),      &
708                                      childarray(3,1,2):childarray(3,2,2),      &
709                                      childarray(4,1,2):childarray(4,2,2)) =    &
710                         tempC%array4(childarray(1,1,1):childarray(1,2,1),      &
711                                      childarray(2,1,1):childarray(2,2,1),      &
712                                      childarray(3,1,1):childarray(3,2,1),      &
713                                      childarray(4,1,1):childarray(4,2,1))
714                    case (5)
715                        parray5(childarray(1,1,2):childarray(1,2,2),      &
716                                      childarray(2,1,2):childarray(2,2,2),      &
717                                      childarray(3,1,2):childarray(3,2,2),      &
718                                      childarray(4,1,2):childarray(4,2,2),      &
719                                      childarray(5,1,2):childarray(5,2,2)) =    &
720                         tempC%array5(childarray(1,1,1):childarray(1,2,1),      &
721                                      childarray(2,1,1):childarray(2,2,1),      &
722                                      childarray(3,1,1):childarray(3,2,1),      &
723                                      childarray(4,1,1):childarray(4,2,1),      &
724                                      childarray(5,1,1):childarray(5,2,1))
725                    case (6)
726                        parray6(childarray(1,1,2):childarray(1,2,2),      &
727                                      childarray(2,1,2):childarray(2,2,2),      &
728                                      childarray(3,1,2):childarray(3,2,2),      &
729                                      childarray(4,1,2):childarray(4,2,2),      &
730                                      childarray(5,1,2):childarray(5,2,2),      &
731                                      childarray(6,1,2):childarray(6,2,2)) =    &
732                         tempC%array6(childarray(1,1,1):childarray(1,2,1),      &
733                                      childarray(2,1,1):childarray(2,2,1),      &
734                                      childarray(3,1,1):childarray(3,2,1),      &
735                                      childarray(4,1,1):childarray(4,2,1),      &
736                                      childarray(5,1,1):childarray(5,2,1),      &
737                                      childarray(6,1,1):childarray(6,2,1))
738                    end select
739                endif  ! < (.not.in_bc)
740            endif  ! < memberin
741!
742        endif
743
744        call Agrif_array_deallocate(tempPextend,nbdim)
745        call Agrif_array_deallocate(tempC,nbdim)
746
747    endif
748!
749!   Deallocations
750#if defined AGRIF_MPI
751    if (member) then
752        call Agrif_array_deallocate(tempP,nbdim)
753    endif
754#endif
755!---------------------------------------------------------------------------------------------------
756end subroutine Agrif_InterpnD
757!===================================================================================================
758!
759!===================================================================================================
760!  subroutine Agrif_Parentbounds
761!
762!> Calculates the bounds of the parent grid for the interpolation of the child grid
763!---------------------------------------------------------------------------------------------------
764subroutine Agrif_Parentbounds ( type_interp, nbdim, indmin, indmax, &
765                                s_Parent_temp, s_Child_temp,        &
766                                s_Child, ds_Child,                  &
767                                s_Parent,ds_Parent,                 &
768                                pttruetab, cetruetab,               &
769                                pttab_Child, pttab_Parent, posvar, coords )
770!---------------------------------------------------------------------------------------------------
771    INTEGER, DIMENSION(6),     intent(in)  :: type_interp
772    INTEGER,                   intent(in)  :: nbdim
773    INTEGER, DIMENSION(nbdim), intent(out) :: indmin, indmax
774    REAL,    DIMENSION(nbdim), intent(out) :: s_Parent_temp, s_child_temp
775    REAL,    DIMENSION(nbdim), intent(in)  :: s_Child, ds_child
776    REAL,    DIMENSION(nbdim), intent(in)  :: s_Parent,ds_Parent
777    INTEGER, DIMENSION(nbdim), intent(in)  :: pttruetab, cetruetab
778    INTEGER, DIMENSION(nbdim), intent(in)  :: pttab_Child, pttab_Parent
779    INTEGER, DIMENSION(nbdim), intent(in)  :: posvar
780    INTEGER, DIMENSION(nbdim), intent(in)  :: coords
781!
782    INTEGER :: i
783    REAL,DIMENSION(nbdim) :: dim_newmin, dim_newmax
784!
785    dim_newmin = s_Child + (pttruetab - pttab_Child) * ds_Child
786    dim_newmax = s_Child + (cetruetab - pttab_Child) * ds_Child
787!
788    do i = 1,nbdim
789!
790        indmin(i) = pttab_Parent(i) + agrif_int((dim_newmin(i)-s_Parent(i))/ds_Parent(i))
791        indmax(i) = pttab_Parent(i) + agrif_ceiling((dim_newmax(i)-s_Parent(i))/ds_Parent(i))
792!
793!       Necessary for the Quadratic interpolation
794!
795        if ( (pttruetab(i) == cetruetab(i)) .and. (posvar(i) == 1) ) then
796        elseif ( coords(i) == 0 ) then  ! (interptab == 'N')
797        elseif ( (type_interp(i) == Agrif_ppm)     .or.     &
798                 (type_interp(i) == Agrif_eno)     .or.     &
799                 (type_interp(i) == Agrif_ppm_lim) .or.     &
800                 (type_interp(i) == Agrif_weno) ) then
801            indmin(i) = indmin(i) - 2
802            indmax(i) = indmax(i) + 2
803        elseif ( (type_interp(i) /= Agrif_constant) .and.   &
804                 (type_interp(i) /= Agrif_linear) ) then
805            indmin(i) = indmin(i) - 1
806            indmax(i) = indmax(i) + 1
807        endif
808!
809    enddo
810!
811    s_Parent_temp = s_Parent + (indmin - pttab_Parent) * ds_Parent
812    s_Child_temp  = s_Child + (pttruetab - pttab_Child) * ds_Child
813!---------------------------------------------------------------------------------------------------
814end subroutine Agrif_Parentbounds
815!===================================================================================================
816!
817!===================================================================================================
818!  subroutine Agrif_Interp_1D_Recursive
819!
820!> Subroutine for the interpolation of a 1D grid variable.
821!> It calls Agrif_InterpBase.
822!---------------------------------------------------------------------------------------------------
823subroutine Agrif_Interp_1D_recursive ( type_interp, tabin, tabout,  &
824                                       indmin, indmax,              &
825                                       pttab_child, petab_child,    &
826                                       s_child,  s_parent,          &
827                                       ds_child, ds_parent )
828!---------------------------------------------------------------------------------------------------
829    integer,            intent(in)  :: type_interp
830    integer,            intent(in)  :: indmin, indmax
831    integer,            intent(in)  :: pttab_child, petab_child
832    real,               intent(in)  :: s_child, s_parent
833    real,               intent(in)  :: ds_child, ds_parent
834    real, dimension(            &
835        indmin:indmax           &
836    ),                  intent(in)  :: tabin
837    real, dimension(            &
838        pttab_child:petab_child &
839    ),                  intent(out) :: tabout
840!---------------------------------------------------------------------------------------------------
841    call Agrif_InterpBase(type_interp,                      &
842                          tabin(indmin:indmax),             &
843                          tabout(pttab_child:petab_child),  &
844                          indmin, indmax,                   &
845                          pttab_child, petab_child,         &
846                          s_parent,    s_child,             &
847                         ds_parent,   ds_child)
848!---------------------------------------------------------------------------------------------------
849end subroutine Agrif_Interp_1D_recursive
850!===================================================================================================
851!
852!===================================================================================================
853!  subroutine Agrif_Interp_2D_Recursive
854!
855!> Subroutine for the interpolation of a 2D grid variable.
856!> It calls Agrif_Interp_1D_recursive and Agrif_InterpBase.
857!---------------------------------------------------------------------------------------------------
858subroutine Agrif_Interp_2D_recursive ( type_interp, tabin, tabout,  &
859                                       indmin, indmax,              &
860                                       pttab_child, petab_child,    &
861                                       s_child,  s_parent,          &
862                                       ds_child, ds_parent )
863!---------------------------------------------------------------------------------------------------
864    integer, dimension(2),              intent(in)  :: type_interp
865    integer, dimension(2),              intent(in)  :: indmin, indmax
866    integer, dimension(2),              intent(in)  :: pttab_child, petab_child
867    real,    dimension(2),              intent(in)  :: s_child, s_parent
868    real,    dimension(2),              intent(in)  :: ds_child, ds_parent
869    real,    dimension(                 &
870        indmin(1):indmax(1),            &
871        indmin(2):indmax(2)),           intent(in)  :: tabin
872    real,    dimension(                 &
873        pttab_child(1):petab_child(1),  &
874        pttab_child(2):petab_child(2)), intent(out) :: tabout
875!---------------------------------------------------------------------------------------------------
876    real, dimension(                    &
877        pttab_child(1):petab_child(1),  &
878        indmin(2):indmax(2))            :: tabtemp
879    real, dimension(                    &
880        pttab_child(2):petab_child(2),  &
881        pttab_child(1):petab_child(1))  :: tabout_trsp
882    real, dimension(                    &
883        indmin(2):indmax(2),            &
884        pttab_child(1):petab_child(1))  :: tabtemp_trsp
885    integer                             :: i, j, coeffraf
886!---------------------------------------------------------------------------------------------------
887!
888    coeffraf = nint ( ds_parent(1) / ds_child(1) )
889!
890    if ((type_interp(1) == Agrif_Linear) .and. (coeffraf /= 1)) then
891!---CDIR NEXPAND
892        if(.NOT. precomputedone(1))     &
893            call Linear1dPrecompute2d(                  &
894                    indmax(2)-indmin(2)+1,              &
895                    indmax(1)-indmin(1)+1,              &
896                    petab_child(1)-pttab_child(1)+1,    &
897                    s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
898!---CDIR NEXPAND
899        call Linear1dAfterCompute(tabin,tabtemp,size(tabin),size(tabtemp),1)
900!
901    elseif ((type_interp(1) == Agrif_PPM) .and. (coeffraf /= 1)) then
902!---CDIR NEXPAND
903        if(.NOT. precomputedone(1))     &
904            call PPM1dPrecompute2d(                     &
905                    indmax(2)-indmin(2)+1,              &
906                    indmax(1)-indmin(1)+1,              &
907                    petab_child(1)-pttab_child(1)+1,        &
908                    s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
909!---CDIR NEXPAND
910        call PPM1dAfterCompute(tabin,tabtemp,size(tabin),size(tabtemp),1)
911    else
912        do j = indmin(2),indmax(2)
913!
914!---CDIR NEXPAND
915            call Agrif_Interp_1D_recursive(type_interp(1),                  &
916                    tabin(indmin(1):indmax(1),j),               &
917                    tabtemp(pttab_child(1):petab_child(1),j),   &
918                    indmin(1),indmax(1),                    &
919                    pttab_child(1),petab_child(1),          &
920                    s_child(1), s_parent(1),                &
921                    ds_child(1),ds_parent(1))
922!
923        enddo
924    endif
925
926    coeffraf = nint(ds_parent(2)/ds_child(2))
927    tabtemp_trsp = TRANSPOSE(tabtemp)
928
929    if ((type_interp(2) == Agrif_Linear) .and. (coeffraf /= 1)) then
930!---CDIR NEXPAND
931        if(.NOT. precomputedone(2))     &
932            call Linear1dPrecompute2d(                  &
933                    petab_child(1)-pttab_child(1)+1,    &
934                    indmax(2)-indmin(2)+1,              &
935                    petab_child(2)-pttab_child(2)+1,    &
936                    s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
937!---CDIR NEXPAND
938        call Linear1dAfterCompute(tabtemp_trsp,tabout_trsp, &
939                size(tabtemp_trsp),size(tabout_trsp),2)
940
941    elseif ((type_interp(2) == Agrif_PPM) .and. (coeffraf /= 1)) then
942!---CDIR NEXPAND
943        if(.NOT. precomputedone(2))     &
944            call PPM1dPrecompute2d(                     &
945                    petab_child(1)-pttab_child(1)+1,    &
946                    indmax(2)-indmin(2)+1,              &
947                    petab_child(2)-pttab_child(2)+1,    &
948                    s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
949!---CDIR NEXPAND
950        call PPM1dAfterCompute(tabtemp_trsp, tabout_trsp,    &
951                               size(tabtemp_trsp), size(tabout_trsp), 2)
952    else
953        do i = pttab_child(1), petab_child(1)
954!
955!---CDIR NEXPAND
956            call Agrif_InterpBase(type_interp(2),                                   &
957                                  tabtemp_trsp(indmin(2):indmax(2), i),             &
958                                  tabout_trsp(pttab_child(2):petab_child(2), i),    &
959                                  indmin(2), indmax(2),                             &
960                                  pttab_child(2), petab_child(2),                   &
961                                  s_parent(2),    s_child(2),                       &
962                                 ds_parent(2),   ds_child(2) )
963        enddo
964    endif
965!
966    tabout = TRANSPOSE(tabout_trsp)
967!---------------------------------------------------------------------------------------------------
968end subroutine Agrif_Interp_2D_recursive
969!===================================================================================================
970!
971!===================================================================================================
972!  subroutine Agrif_Interp_3D_Recursive
973!
974!> Subroutine for the interpolation of a 3D grid variable.
975!> It calls #Agrif_Interp_2D_recursive and #Agrif_InterpBase.
976!---------------------------------------------------------------------------------------------------
977subroutine Agrif_Interp_3D_recursive ( type_interp, tabin, tabout,  &
978                                       indmin, indmax,              &
979                                       pttab_child, petab_child,    &
980                                       s_child,  s_parent,          &
981                                       ds_child, ds_parent )
982!---------------------------------------------------------------------------------------------------
983    integer, dimension(3),              intent(in)  :: type_interp
984    integer, dimension(3),              intent(in)  :: indmin, indmax
985    integer, dimension(3),              intent(in)  :: pttab_child, petab_child
986    real,    dimension(3),              intent(in)  :: s_child, s_parent
987    real,    dimension(3),              intent(in)  :: ds_child, ds_parent
988    real,    dimension(                 &
989        indmin(1):indmax(1),            &
990        indmin(2):indmax(2),            &
991        indmin(3):indmax(3)),           intent(in)  :: tabin
992    real,    dimension(                 &
993        pttab_child(1):petab_child(1),  &
994        pttab_child(2):petab_child(2),  &
995        pttab_child(3):petab_child(3)), intent(out) :: tabout
996!---------------------------------------------------------------------------------------------------
997    real, dimension(                    &
998        pttab_child(1):petab_child(1),  &
999        pttab_child(2):petab_child(2),  &
1000        indmin(3):indmax(3))            :: tabtemp
1001    integer                             :: i, j, k, coeffraf
1002    integer                             :: locind_child_left, kdeb
1003!
1004    coeffraf = nint ( ds_parent(1) / ds_child(1) )
1005    if ( (type_interp(1) == Agrif_Linear) .and. (coeffraf/=1) ) then
1006        call Linear1dPrecompute2d(indmax(2)-indmin(2)+1,            &
1007                                  indmax(1)-indmin(1)+1,            &
1008                                  petab_child(1)-pttab_child(1)+1,  &
1009                                  s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1010        precomputedone(1) = .TRUE.
1011    elseif ( (type_interp(1) == Agrif_PPM) .and. (coeffraf/=1) ) then
1012        call PPM1dPrecompute2d(indmax(2)-indmin(2)+1,           &
1013                               indmax(1)-indmin(1)+1,           &
1014                               petab_child(1)-pttab_child(1)+1, &
1015                               s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1016        precomputedone(1) = .TRUE.
1017    endif
1018
1019    coeffraf = nint ( ds_parent(2) / ds_child(2) )
1020    if ( (type_interp(2) == Agrif_Linear) .and. (coeffraf/=1) ) then
1021        call Linear1dPrecompute2d(petab_child(1)-pttab_child(1)+1,  &
1022                                  indmax(2)-indmin(2)+1,            &
1023                                  petab_child(2)-pttab_child(2)+1,  &
1024                                  s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1025        precomputedone(2) = .TRUE.
1026    elseif ( (type_interp(2) == Agrif_PPM) .and. (coeffraf/=1) ) then
1027        call PPM1dPrecompute2d(petab_child(1)-pttab_child(1)+1, &
1028                               indmax(2)-indmin(2)+1,           &
1029                               petab_child(2)-pttab_child(2)+1, &
1030                               s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1031        precomputedone(2) = .TRUE.
1032    endif
1033!
1034    do k = indmin(3), indmax(3)
1035        call Agrif_Interp_2D_recursive(type_interp(1:2),                            &
1036                                       tabin(indmin(1):indmax(1),                   &
1037                                             indmin(2):indmax(2), k),               &
1038                                       tabtemp(pttab_child(1):petab_child(1),       &
1039                                               pttab_child(2):petab_child(2), k),   &
1040                                       indmin(1:2), indmax(1:2),                    &
1041                                       pttab_child(1:2), petab_child(1:2),          &
1042                                       s_child(1:2),     s_parent(1:2),             &
1043                                       ds_child(1:2),   ds_parent(1:2) )
1044    enddo
1045!
1046    precomputedone(1) = .FALSE.
1047    precomputedone(2) = .FALSE.
1048    coeffraf = nint(ds_parent(3)/ds_child(3))
1049!
1050    if ( coeffraf == 1 ) then
1051        locind_child_left = 1 + agrif_int((s_child(3)-s_parent(3))/ds_parent(3))
1052        kdeb = indmin(3)+locind_child_left-2
1053        do k = pttab_child(3),petab_child(3)
1054            kdeb = kdeb + 1
1055            do j = pttab_child(2), petab_child(2)
1056            do i = pttab_child(1), petab_child(1)
1057                tabout(i,j,k) = tabtemp(i,j,kdeb)
1058            enddo
1059            enddo
1060        enddo
1061    else
1062        do j = pttab_child(2), petab_child(2)
1063        do i = pttab_child(1), petab_child(1)
1064            call Agrif_InterpBase(type_interp(3),                                   &
1065                                  tabtemp(i,j,indmin(3):indmax(3)),                 &
1066                                  tabout(i,j,pttab_child(3):petab_child(3)),        &
1067                                  indmin(3), indmax(3),                             &
1068                                  pttab_child(3), petab_child(3),                   &
1069                                  s_parent(3),    s_child(3),                       &
1070                                 ds_parent(3),   ds_child(3) )
1071        enddo
1072        enddo
1073    endif
1074!---------------------------------------------------------------------------------------------------
1075end subroutine Agrif_Interp_3D_recursive
1076!===================================================================================================
1077!
1078!===================================================================================================
1079!  subroutine Agrif_Interp_4D_Recursive
1080!
1081!> Subroutine for the interpolation of a 4D grid variable.
1082!> It calls #Agrif_Interp_3D_recursive and #Agrif_InterpBase.
1083!---------------------------------------------------------------------------------------------------
1084subroutine Agrif_Interp_4D_recursive ( type_interp, tabin, tabout,  &
1085                                       indmin, indmax,              &
1086                                       pttab_child, petab_child,    &
1087                                       s_child,  s_parent,          &
1088                                       ds_child, ds_parent )
1089!---------------------------------------------------------------------------------------------------
1090    integer, dimension(4),              intent(in)  :: type_interp
1091    integer, dimension(4),              intent(in)  :: indmin, indmax
1092    integer, dimension(4),              intent(in)  :: pttab_child, petab_child
1093    real,    dimension(4),              intent(in)  :: s_child, s_parent
1094    real,    dimension(4),              intent(in)  :: ds_child, ds_parent
1095    real,    dimension(                 &
1096        indmin(1):indmax(1),            &
1097        indmin(2):indmax(2),            &
1098        indmin(3):indmax(3),            &
1099        indmin(4):indmax(4)),           intent(in)  :: tabin
1100    real,    dimension(                 &
1101        pttab_child(1):petab_child(1),  &
1102        pttab_child(2):petab_child(2),  &
1103        pttab_child(3):petab_child(3),  &
1104        pttab_child(4):petab_child(4)), intent(out) :: tabout
1105!---------------------------------------------------------------------------------------------------
1106    real, dimension(                    &
1107        pttab_child(1):petab_child(1),  &
1108        pttab_child(2):petab_child(2),  &
1109        pttab_child(3):petab_child(3),  &
1110        indmin(4):indmax(4))            :: tabtemp
1111    integer                             :: i, j, k, l
1112!
1113    do l = indmin(4), indmax(4)
1114        call Agrif_Interp_3D_recursive(type_interp(1:3),                            &
1115                                       tabin(indmin(1):indmax(1),                   &
1116                                             indmin(2):indmax(2),                   &
1117                                             indmin(3):indmax(3), l),               &
1118                                       tabtemp(pttab_child(1):petab_child(1),       &
1119                                               pttab_child(2):petab_child(2),       &
1120                                               pttab_child(3):petab_child(3), l),   &
1121                                       indmin(1:3), indmax(1:3),                    &
1122                                       pttab_child(1:3), petab_child(1:3),          &
1123                                       s_child(1:3),    s_parent(1:3),              &
1124                                       ds_child(1:3),  ds_parent(1:3) )
1125    enddo
1126!
1127    do k = pttab_child(3), petab_child(3)
1128    do j = pttab_child(2), petab_child(2)
1129    do i = pttab_child(1), petab_child(1)
1130        call Agrif_InterpBase(type_interp(4),                                       &
1131                              tabtemp(i,j,k,indmin(4):indmax(4)),                   &
1132                              tabout(i,j,k,pttab_child(4):petab_child(4)),          &
1133                              indmin(4), indmax(4),                                 &
1134                              pttab_child(4), petab_child(4),                       &
1135                              s_parent(4),    s_child(4),                           &
1136                             ds_parent(4),   ds_child(4) )
1137    enddo
1138    enddo
1139    enddo
1140!---------------------------------------------------------------------------------------------------
1141end subroutine Agrif_Interp_4D_recursive
1142!===================================================================================================
1143!
1144!===================================================================================================
1145!  subroutine Agrif_Interp_5D_Recursive
1146!
1147!> Subroutine for the interpolation of a 5D grid variable.
1148!> It calls #Agrif_Interp_4D_recursive and #Agrif_InterpBase.
1149!---------------------------------------------------------------------------------------------------
1150subroutine Agrif_Interp_5D_recursive ( type_interp, tabin, tabout,  &
1151                                       indmin, indmax,              &
1152                                       pttab_child, petab_child,    &
1153                                       s_child,  s_parent,          &
1154                                       ds_child, ds_parent )
1155!---------------------------------------------------------------------------------------------------
1156    integer, dimension(5),              intent(in)  :: type_interp
1157    integer, dimension(5),              intent(in)  :: indmin, indmax
1158    integer, dimension(5),              intent(in)  :: pttab_child, petab_child
1159    real,    dimension(5),              intent(in)  :: s_child, s_parent
1160    real,    dimension(5),              intent(in)  :: ds_child, ds_parent
1161    real,    dimension(                 &
1162        indmin(1):indmax(1),            &
1163        indmin(2):indmax(2),            &
1164        indmin(3):indmax(3),            &
1165        indmin(4):indmax(4),            &
1166        indmin(5):indmax(5)),           intent(in)  :: tabin
1167    real,    dimension(                 &
1168        pttab_child(1):petab_child(1),  &
1169        pttab_child(2):petab_child(2),  &
1170        pttab_child(3):petab_child(3),  &
1171        pttab_child(4):petab_child(4),  &
1172        pttab_child(5):petab_child(5)), intent(out) :: tabout
1173!---------------------------------------------------------------------------------------------------
1174    real, dimension(                    &
1175        pttab_child(1):petab_child(1),  &
1176        pttab_child(2):petab_child(2),  &
1177        pttab_child(3):petab_child(3),  &
1178        pttab_child(4):petab_child(4),  &
1179        indmin(5):indmax(5))            :: tabtemp
1180    integer                             :: i, j, k, l, m
1181!
1182    do m = indmin(5), indmax(5)
1183        call Agrif_Interp_4D_recursive(type_interp(1:4),                            &
1184                                       tabin(indmin(1):indmax(1),                   &
1185                                             indmin(2):indmax(2),                   &
1186                                             indmin(3):indmax(3),                   &
1187                                             indmin(4):indmax(4),m),                &
1188                                       tabtemp(pttab_child(1):petab_child(1),       &
1189                                               pttab_child(2):petab_child(2),       &
1190                                               pttab_child(3):petab_child(3),       &
1191                                               pttab_child(4):petab_child(4), m),   &
1192                                       indmin(1:4),indmax(1:4),                     &
1193                                       pttab_child(1:4), petab_child(1:4),          &
1194                                       s_child(1:4),     s_parent(1:4),             &
1195                                       ds_child(1:4),   ds_parent(1:4) )
1196    enddo
1197!
1198    do l = pttab_child(4), petab_child(4)
1199    do k = pttab_child(3), petab_child(3)
1200    do j = pttab_child(2), petab_child(2)
1201    do i = pttab_child(1), petab_child(1)
1202        call Agrif_InterpBase(type_interp(5),                                       &
1203                              tabtemp(i,j,k,l,indmin(5):indmax(5)),                 &
1204                              tabout(i,j,k,l,pttab_child(5):petab_child(5)),        &
1205                              indmin(5), indmax(5),                                 &
1206                              pttab_child(5), petab_child(5),                       &
1207                              s_parent(5),   s_child(5),                            &
1208                              ds_parent(5), ds_child(5) )
1209    enddo
1210    enddo
1211    enddo
1212    enddo
1213!---------------------------------------------------------------------------------------------------
1214end subroutine Agrif_Interp_5D_recursive
1215!===================================================================================================
1216!
1217!===================================================================================================
1218!  subroutine Agrif_Interp_6D_Recursive
1219!
1220!> Subroutine for the interpolation of a 6D grid variable.
1221!> It calls #Agrif_Interp_5D_recursive and Agrif_InterpBase.
1222!---------------------------------------------------------------------------------------------------
1223subroutine Agrif_Interp_6D_recursive ( type_interp, tabin, tabout,  &
1224                                       indmin, indmax,              &
1225                                       pttab_child, petab_child,    &
1226                                       s_child,  s_parent,          &
1227                                       ds_child, ds_parent )
1228!---------------------------------------------------------------------------------------------------
1229    integer, dimension(6),                  intent(in)  :: type_interp
1230    integer, dimension(6),                  intent(in)  :: indmin, indmax
1231    integer, dimension(6),                  intent(in)  :: pttab_child, petab_child
1232    real,    dimension(6),                  intent(in)  :: s_child, s_parent
1233    real,    dimension(6),                  intent(in)  :: ds_child, ds_parent
1234    real,    dimension(                 &
1235        indmin(1):indmax(1),            &
1236        indmin(2):indmax(2),            &
1237        indmin(3):indmax(3),            &
1238        indmin(4):indmax(4),            &
1239        indmin(5):indmax(5),            &
1240        indmin(6):indmax(6)),               intent(in)  :: tabin
1241    real,    dimension(                 &
1242        pttab_child(1):petab_child(1),  &
1243        pttab_child(2):petab_child(2),  &
1244        pttab_child(3):petab_child(3),  &
1245        pttab_child(4):petab_child(4),  &
1246        pttab_child(5):petab_child(5),  &
1247        pttab_child(6):petab_child(6)),     intent(out) :: tabout
1248!---------------------------------------------------------------------------------------------------
1249    real, dimension(                    &
1250        pttab_child(1):petab_child(1),  &
1251        pttab_child(2):petab_child(2),  &
1252        pttab_child(3):petab_child(3),  &
1253        pttab_child(4):petab_child(4),  &
1254        pttab_child(5):petab_child(5),  &
1255        indmin(6):indmax(6))            :: tabtemp
1256    integer                             :: i, j, k, l, m, n
1257!
1258    do n = indmin(6), indmax(6)
1259        call Agrif_Interp_5D_recursive(type_interp(1:5),                            &
1260                                       tabin(indmin(1):indmax(1),                   &
1261                                             indmin(2):indmax(2),                   &
1262                                             indmin(3):indmax(3),                   &
1263                                             indmin(4):indmax(4),                   &
1264                                             indmin(5):indmax(5), n),               &
1265                                        tabtemp(pttab_child(1):petab_child(1),      &
1266                                                pttab_child(2):petab_child(2),      &
1267                                                pttab_child(3):petab_child(3),      &
1268                                                pttab_child(4):petab_child(4),      &
1269                                                pttab_child(5):petab_child(5), n),  &
1270                                        indmin(1:5),indmax(1:5),                    &
1271                                        pttab_child(1:5), petab_child(1:5),         &
1272                                        s_child(1:5), s_parent(1:5),                &
1273                                        ds_child(1:5),ds_parent(1:5) )
1274    enddo
1275!
1276    do m = pttab_child(5), petab_child(5)
1277    do l = pttab_child(4), petab_child(4)
1278    do k = pttab_child(3), petab_child(3)
1279    do j = pttab_child(2), petab_child(2)
1280    do i = pttab_child(1), petab_child(1)
1281        call Agrif_InterpBase(type_interp(6),                                       &
1282                              tabtemp(i,j,k,l,m,indmin(6):indmax(6)),               &
1283                              tabout(i,j,k,l,m,pttab_child(6):petab_child(6)),      &
1284                              indmin(6), indmax(6),                                 &
1285                              pttab_child(6), petab_child(6),                       &
1286                              s_parent(6),   s_child(6),                            &
1287                              ds_parent(6), ds_child(6) )
1288    enddo
1289    enddo
1290    enddo
1291    enddo
1292    enddo
1293!---------------------------------------------------------------------------------------------------
1294end subroutine Agrif_Interp_6D_recursive
1295!===================================================================================================
1296!
1297!===================================================================================================
1298!  subroutine Agrif_InterpBase
1299!
1300!> Calls the interpolation method chosen by the user (linear, lagrange, spline, etc.).
1301!---------------------------------------------------------------------------------------------------
1302subroutine Agrif_InterpBase ( type_interp, parenttab, childtab, indmin, indmax, &
1303                              pttab_child, petab_child,                         &
1304                              s_parent, s_child, ds_parent, ds_child )
1305!---------------------------------------------------------------------------------------------------
1306    INTEGER                                                 :: type_interp
1307    INTEGER                                                 :: indmin, indmax
1308    INTEGER                                                 :: pttab_child, petab_child
1309    REAL, DIMENSION(indmin:indmax),           INTENT(IN)    :: parenttab
1310    REAL, DIMENSION(pttab_child:petab_child), INTENT(OUT)   :: childtab
1311    REAL                                                    :: s_parent, s_child
1312    REAL                                                    :: ds_parent,ds_child
1313!
1314    if ( (indmin == indmax) .and. (pttab_child == petab_child) ) then
1315!
1316        childtab(pttab_child) = parenttab(indmin)
1317!
1318    elseif (type_interp == Agrif_LINEAR) then    !       Linear interpolation
1319!
1320        call Agrif_basicinterp_linear1D(parenttab,childtab,                   &
1321                indmax-indmin+1,petab_child-pttab_child+1,  &
1322                s_parent,s_child,ds_parent,ds_child)
1323!
1324    elseif ( type_interp == Agrif_PPM ) then     !       PPM interpolation
1325
1326        call PPM1d(parenttab,childtab,                      &
1327                indmax-indmin+1,petab_child-pttab_child+1,  &
1328                s_parent,s_child,ds_parent,ds_child)
1329!
1330    elseif ( type_interp == Agrif_PPM_LIM ) then     !       PPM interpolation
1331
1332        call PPM1d_lim(parenttab,childtab,                      &
1333                indmax-indmin+1,petab_child-pttab_child+1,  &
1334                s_parent,s_child,ds_parent,ds_child)
1335!
1336    elseif (type_interp == Agrif_LAGRANGE) then  !       Lagrange interpolation
1337!
1338        call lagrange1D(parenttab,childtab,                 &
1339                indmax-indmin+1,petab_child-pttab_child+1,  &
1340                s_parent,s_child,ds_parent,ds_child)
1341!
1342    elseif (type_interp == Agrif_ENO) then       !       Eno interpolation
1343!
1344        call ENO1d(parenttab,childtab,                      &
1345                indmax-indmin+1,petab_child-pttab_child+1,  &
1346                s_parent,s_child,ds_parent,ds_child)
1347!
1348    elseif (type_interp == Agrif_WENO) then      !       Weno interpolation
1349!
1350        call WENO1d(parenttab,childtab,                     &
1351                indmax-indmin+1,petab_child-pttab_child+1,  &
1352                s_parent,s_child,ds_parent,ds_child)
1353!
1354    elseif (type_interp == Agrif_LINEARCONSERV) then !   Linear conservative interpolation
1355!
1356        call Linear1dConserv(parenttab,childtab,            &
1357                indmax-indmin+1,petab_child-pttab_child+1,  &
1358                s_parent,s_child,ds_parent,ds_child)
1359!
1360    elseif (type_interp == Agrif_LINEARCONSERVLIM) then !Linear conservative interpolation
1361!
1362        call Linear1dConservLim(parenttab,childtab,         &
1363                indmax-indmin+1,petab_child-pttab_child+1,  &
1364                s_parent,s_child,ds_parent,ds_child)
1365!
1366    elseif (type_interp == Agrif_CONSTANT) then
1367!
1368        call Constant1d(parenttab,childtab,                 &
1369                indmax-indmin+1,petab_child-pttab_child+1,  &
1370                s_parent,s_child,ds_parent,ds_child)
1371!
1372    endif
1373!---------------------------------------------------------------------------------------------------
1374end subroutine Agrif_InterpBase
1375!===================================================================================================
1376!
1377!===================================================================================================
1378!  subroutine Agrif_Find_list_interp
1379!---------------------------------------------------------------------------------------------------
1380function Agrif_Find_list_interp ( list_interp, pttab, petab, pttab_Child, pttab_Parent,     &
1381                                    nbdim, indmin, indmax, indminglob,  indmaxglob,         &
1382                                    pttruetab, cetruetab, memberin                          &
1383#if defined AGRIF_MPI
1384                                   ,indminglob2, indmaxglob2, parentarray,                  &
1385                                    member, tab4t, memberinall, sendtoproc1, recvfromproc1  &
1386#endif
1387    ) result(find_list_interp)
1388!---------------------------------------------------------------------------------------------------
1389    type(Agrif_List_Interp_Loc),   pointer     :: list_interp
1390    integer,                       intent(in)  :: nbdim
1391    integer, dimension(nbdim),     intent(in)  :: pttab, petab, pttab_Child, pttab_Parent
1392    integer, dimension(nbdim),     intent(out) :: indmin, indmax
1393    integer, dimension(nbdim),     intent(out) :: indminglob, indmaxglob
1394    integer, dimension(nbdim),     intent(out) :: pttruetab, cetruetab
1395    logical,                       intent(out) :: memberin
1396#if defined AGRIF_MPI
1397    integer, dimension(nbdim),     intent(out) :: indminglob2, indmaxglob2
1398    integer, dimension(nbdim,2,2), intent(out) :: parentarray
1399    logical,                       intent(out) :: member
1400    integer, dimension(nbdim,0:Agrif_Nbprocs-1,8), intent(out) :: tab4t
1401    logical, dimension(0:Agrif_Nbprocs-1),         intent(out) :: memberinall
1402    logical, dimension(0:Agrif_Nbprocs-1),         intent(out) :: sendtoproc1, recvfromproc1
1403#endif
1404    logical :: find_list_interp
1405!
1406    integer :: i
1407    type(Agrif_List_Interp_Loc), pointer :: parcours
1408    type(Agrif_Interp_Loc),      pointer :: pil
1409
1410    find_list_interp = .false.
1411
1412    if ( .not. associated(list_interp) )    return
1413
1414    parcours => list_interp
1415    find_loop : do while ( associated(parcours) )
1416
1417        pil => parcours % interp_loc
1418
1419        do i = 1,nbdim
1420            if ( (pttab(i) /= pil % pttab(i)) .or. &
1421                 (petab(i) /= pil % petab(i)) .or. &
1422                 (pttab_child(i)  /= pil % pttab_child(i)) .or. &
1423                 (pttab_parent(i) /= pil % pttab_parent(i)) ) then
1424                parcours => parcours % suiv
1425                cycle find_loop
1426            endif
1427        enddo
1428
1429        indmin = pil % indmin(1:nbdim)
1430        indmax = pil % indmax(1:nbdim)
1431
1432        pttruetab = pil % pttruetab(1:nbdim)
1433        cetruetab = pil % cetruetab(1:nbdim)
1434
1435#if !defined AGRIF_MPI
1436        indminglob  = pil % indminglob(1:nbdim)
1437        indmaxglob  = pil % indmaxglob(1:nbdim)
1438#else
1439        indminglob  = pil % indminglob2(1:nbdim)
1440        indmaxglob  = pil % indmaxglob2(1:nbdim)
1441        indminglob2 = pil % indminglob2(1:nbdim)
1442        indmaxglob2 = pil % indmaxglob2(1:nbdim)
1443        parentarray = pil % parentarray(1:nbdim,:,:)
1444        member      = pil % member
1445        tab4t         = pil % tab4t(1:nbdim, 0:Agrif_Nbprocs-1, 1:8)
1446        memberinall   = pil % memberinall(0:Agrif_Nbprocs-1)
1447        sendtoproc1   = pil % sendtoproc1(0:Agrif_Nbprocs-1)
1448        recvfromproc1 = pil % recvfromproc1(0:Agrif_Nbprocs-1)
1449#endif
1450        memberin = pil % memberin
1451        find_list_interp = .true.
1452        exit find_loop
1453    enddo find_loop
1454!---------------------------------------------------------------------------------------------------
1455end function Agrif_Find_list_interp
1456!===================================================================================================
1457!
1458!===================================================================================================
1459!  subroutine Agrif_AddTo_list_interp
1460!---------------------------------------------------------------------------------------------------
1461subroutine Agrif_AddTo_list_interp ( list_interp, pttab, petab, pttab_Child, pttab_Parent,  &
1462                                     indmin, indmax, indminglob, indmaxglob,                &
1463                                     pttruetab, cetruetab,                                  &
1464                                     memberin, nbdim                                        &
1465#if defined AGRIF_MPI
1466                                    ,indminglob2, indmaxglob2,                              &
1467                                     parentarray,                                           &
1468                                     member,                                                &
1469                                     tab4t, memberinall, sendtoproc1, recvfromproc1         &
1470#endif
1471    )
1472!---------------------------------------------------------------------------------------------------
1473    type(Agrif_List_Interp_Loc), pointer    :: list_interp
1474    integer                                 :: nbdim
1475    integer, dimension(nbdim)               :: pttab, petab, pttab_Child, pttab_Parent
1476    integer, dimension(nbdim)               :: indmin,indmax
1477    integer, dimension(nbdim)               :: indminglob, indmaxglob
1478    integer, dimension(nbdim)               :: pttruetab, cetruetab
1479    logical                                 :: memberin
1480#if defined AGRIF_MPI
1481    integer, dimension(nbdim,2,2)           :: parentarray
1482    logical                                 :: member
1483    integer, dimension(nbdim)                       :: indminglob2,indmaxglob2
1484    integer, dimension(nbdim,0:Agrif_Nbprocs-1,8)   :: tab4t
1485    logical, dimension(0:Agrif_Nbprocs-1)           :: memberinall
1486    logical, dimension(0:Agrif_Nbprocs-1)           :: sendtoproc1
1487    logical, dimension(0:Agrif_Nbprocs-1)           :: recvfromproc1
1488#endif
1489!
1490    type(Agrif_List_Interp_Loc), pointer    :: parcours
1491    type(Agrif_Interp_Loc),      pointer    :: pil
1492!
1493    allocate(parcours)
1494    allocate(parcours % interp_loc)
1495
1496    pil => parcours % interp_loc
1497
1498    pil % pttab(1:nbdim) = pttab(1:nbdim)
1499    pil % petab(1:nbdim) = petab(1:nbdim)
1500    pil % pttab_child(1:nbdim) = pttab_child(1:nbdim)
1501    pil % pttab_parent(1:nbdim) = pttab_parent(1:nbdim)
1502
1503    pil % indmin(1:nbdim) = indmin(1:nbdim)
1504    pil % indmax(1:nbdim) = indmax(1:nbdim)
1505
1506    pil % memberin = memberin
1507#if !defined AGRIF_MPI
1508    pil % indminglob(1:nbdim) = indminglob(1:nbdim)
1509    pil % indmaxglob(1:nbdim) = indmaxglob(1:nbdim)
1510#else
1511    pil % indminglob2(1:nbdim) = indminglob2(1:nbdim)
1512    pil % indmaxglob2(1:nbdim) = indmaxglob2(1:nbdim)
1513    pil % parentarray(1:nbdim,:,:) = parentarray(1:nbdim,:,:)
1514    pil % member = member
1515    allocate(pil % tab4t(nbdim, 0:Agrif_Nbprocs-1, 8))
1516    allocate(pil % memberinall(0:Agrif_Nbprocs-1))
1517    allocate(pil % sendtoproc1(0:Agrif_Nbprocs-1))
1518    allocate(pil % recvfromproc1(0:Agrif_Nbprocs-1))
1519    pil % tab4t         = tab4t
1520    pil % memberinall   = memberinall
1521    pil % sendtoproc1   = sendtoproc1
1522    pil % recvfromproc1 = recvfromproc1
1523#endif
1524
1525    pil % pttruetab(1:nbdim) = pttruetab(1:nbdim)
1526    pil % cetruetab(1:nbdim) = cetruetab(1:nbdim)
1527
1528    parcours % suiv => list_interp
1529    list_interp => parcours
1530!---------------------------------------------------------------------------------------------------
1531end subroutine Agrif_Addto_list_interp
1532!===================================================================================================
1533!
1534end module Agrif_Interpolation
Note: See TracBrowser for help on using the repository browser.