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

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

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

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

  • Property svn:keywords set to Id
File size: 78.1 KB
Line 
1!
2! $Id$
3!
4!     AGRIF (Adaptive Grid Refinement In Fortran)
5!
6!     Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
7!                        Christophe Vouland (Christophe.Vouland@imag.fr)
8!
9!     This program is free software; you can redistribute it and/or modify
10!     it under the terms of the GNU General Public License as published by
11!     the Free Software Foundation; either version 2 of the License, or
12!     (at your option) any later version.
13!
14!     This program is distributed in the hope that it will be useful,
15!     but WITHOUT ANY WARRANTY; without even the implied warranty of
16!     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
17!     GNU General Public License for more details.
18!
19!     You should have received a copy of the GNU General Public License
20!     along with this program; if not, write to the Free Software
21!     Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA.
22!
23!
24!> Module 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
804            if (Agrif_UseSpecialValue) then
805               indmin(i) = indmin(i)-MaxSearch
806               indmax(i) = indmax(i)+MaxSearch
807            endif
808
809        elseif ( (type_interp(i) /= Agrif_constant) .and.   &
810                 (type_interp(i) /= Agrif_linear) ) then
811            indmin(i) = indmin(i) - 1
812            indmax(i) = indmax(i) + 1
813
814            if (Agrif_UseSpecialValue) then
815               indmin(i) = indmin(i)-MaxSearch
816               indmax(i) = indmax(i)+MaxSearch
817            endif
818
819        elseif ( (type_interp(i) == Agrif_constant) .or.   &
820                 (type_interp(i) == Agrif_linear) ) then
821            if (Agrif_UseSpecialValue) then
822               indmin(i) = indmin(i)-MaxSearch
823               indmax(i) = indmax(i)+MaxSearch
824            endif
825
826        endif
827!
828    enddo
829!
830    s_Parent_temp = s_Parent + (indmin - pttab_Parent) * ds_Parent
831    s_Child_temp  = s_Child + (pttruetab - pttab_Child) * ds_Child
832!---------------------------------------------------------------------------------------------------
833end subroutine Agrif_Parentbounds
834!===================================================================================================
835!
836!===================================================================================================
837!  subroutine Agrif_Interp_1D_Recursive
838!
839!> Subroutine for the interpolation of a 1D grid variable.
840!> It calls Agrif_InterpBase.
841!---------------------------------------------------------------------------------------------------
842subroutine Agrif_Interp_1D_recursive ( type_interp, tabin, tabout,  &
843                                       indmin, indmax,              &
844                                       pttab_child, petab_child,    &
845                                       s_child,  s_parent,          &
846                                       ds_child, ds_parent )
847!---------------------------------------------------------------------------------------------------
848    integer,            intent(in)  :: type_interp
849    integer,            intent(in)  :: indmin, indmax
850    integer,            intent(in)  :: pttab_child, petab_child
851    real,               intent(in)  :: s_child, s_parent
852    real,               intent(in)  :: ds_child, ds_parent
853    real, dimension(            &
854        indmin:indmax           &
855    ),                  intent(in)  :: tabin
856    real, dimension(            &
857        pttab_child:petab_child &
858    ),                  intent(out) :: tabout
859!---------------------------------------------------------------------------------------------------
860    call Agrif_InterpBase(type_interp,                      &
861                          tabin(indmin:indmax),             &
862                          tabout(pttab_child:petab_child),  &
863                          indmin, indmax,                   &
864                          pttab_child, petab_child,         &
865                          s_parent,    s_child,             &
866                         ds_parent,   ds_child)
867!---------------------------------------------------------------------------------------------------
868end subroutine Agrif_Interp_1D_recursive
869!===================================================================================================
870!
871!===================================================================================================
872!  subroutine Agrif_Interp_2D_Recursive
873!
874!> Subroutine for the interpolation of a 2D grid variable.
875!> It calls Agrif_Interp_1D_recursive and Agrif_InterpBase.
876!---------------------------------------------------------------------------------------------------
877subroutine Agrif_Interp_2D_recursive ( type_interp, tabin, tabout,  &
878                                       indmin, indmax,              &
879                                       pttab_child, petab_child,    &
880                                       s_child,  s_parent,          &
881                                       ds_child, ds_parent )
882!---------------------------------------------------------------------------------------------------
883    integer, dimension(2),              intent(in)  :: type_interp
884    integer, dimension(2),              intent(in)  :: indmin, indmax
885    integer, dimension(2),              intent(in)  :: pttab_child, petab_child
886    real,    dimension(2),              intent(in)  :: s_child, s_parent
887    real,    dimension(2),              intent(in)  :: ds_child, ds_parent
888    real,    dimension(                 &
889        indmin(1):indmax(1),            &
890        indmin(2):indmax(2)),           intent(in)  :: tabin
891    real,    dimension(                 &
892        pttab_child(1):petab_child(1),  &
893        pttab_child(2):petab_child(2)), intent(out) :: tabout
894!---------------------------------------------------------------------------------------------------
895    real, dimension(                    &
896        pttab_child(1):petab_child(1),  &
897        indmin(2):indmax(2))            :: tabtemp
898    real, dimension(                    &
899        pttab_child(2):petab_child(2),  &
900        pttab_child(1):petab_child(1))  :: tabout_trsp
901    real, dimension(                    &
902        indmin(2):indmax(2),            &
903        pttab_child(1):petab_child(1))  :: tabtemp_trsp
904    integer                             :: i, j, coeffraf
905!---------------------------------------------------------------------------------------------------
906!
907    coeffraf = nint ( ds_parent(1) / ds_child(1) )
908!
909    if ((type_interp(1) == Agrif_Linear) .and. (coeffraf /= 1)) then
910!---CDIR NEXPAND
911        if(.NOT. precomputedone(1))     &
912            call Linear1dPrecompute2d(                  &
913                    indmax(2)-indmin(2)+1,              &
914                    indmax(1)-indmin(1)+1,              &
915                    petab_child(1)-pttab_child(1)+1,    &
916                    s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
917!---CDIR NEXPAND
918        call Linear1dAfterCompute(tabin,tabtemp,size(tabin),size(tabtemp),1)
919!
920    elseif ((type_interp(1) == Agrif_PPM) .and. (coeffraf /= 1)) then
921!---CDIR NEXPAND
922        if(.NOT. precomputedone(1))     &
923            call PPM1dPrecompute2d(                     &
924                    indmax(2)-indmin(2)+1,              &
925                    indmax(1)-indmin(1)+1,              &
926                    petab_child(1)-pttab_child(1)+1,        &
927                    s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
928!---CDIR NEXPAND
929        call PPM1dAfterCompute(tabin,tabtemp,size(tabin),size(tabtemp),1)
930    else
931        do j = indmin(2),indmax(2)
932!
933!---CDIR NEXPAND
934            call Agrif_Interp_1D_recursive(type_interp(1),                  &
935                    tabin(indmin(1):indmax(1),j),               &
936                    tabtemp(pttab_child(1):petab_child(1),j),   &
937                    indmin(1),indmax(1),                    &
938                    pttab_child(1),petab_child(1),          &
939                    s_child(1), s_parent(1),                &
940                    ds_child(1),ds_parent(1))
941!
942        enddo
943    endif
944
945    coeffraf = nint(ds_parent(2)/ds_child(2))
946    tabtemp_trsp = TRANSPOSE(tabtemp)
947
948    if ((type_interp(2) == Agrif_Linear) .and. (coeffraf /= 1)) then
949!---CDIR NEXPAND
950        if(.NOT. precomputedone(2))     &
951            call Linear1dPrecompute2d(                  &
952                    petab_child(1)-pttab_child(1)+1,    &
953                    indmax(2)-indmin(2)+1,              &
954                    petab_child(2)-pttab_child(2)+1,    &
955                    s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
956!---CDIR NEXPAND
957        call Linear1dAfterCompute(tabtemp_trsp,tabout_trsp, &
958                size(tabtemp_trsp),size(tabout_trsp),2)
959
960    elseif ((type_interp(2) == Agrif_PPM) .and. (coeffraf /= 1)) then
961!---CDIR NEXPAND
962        if(.NOT. precomputedone(2))     &
963            call PPM1dPrecompute2d(                     &
964                    petab_child(1)-pttab_child(1)+1,    &
965                    indmax(2)-indmin(2)+1,              &
966                    petab_child(2)-pttab_child(2)+1,    &
967                    s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
968!---CDIR NEXPAND
969        call PPM1dAfterCompute(tabtemp_trsp, tabout_trsp,    &
970                               size(tabtemp_trsp), size(tabout_trsp), 2)
971    else
972        do i = pttab_child(1), petab_child(1)
973!
974!---CDIR NEXPAND
975            call Agrif_InterpBase(type_interp(2),                                   &
976                                  tabtemp_trsp(indmin(2):indmax(2), i),             &
977                                  tabout_trsp(pttab_child(2):petab_child(2), i),    &
978                                  indmin(2), indmax(2),                             &
979                                  pttab_child(2), petab_child(2),                   &
980                                  s_parent(2),    s_child(2),                       &
981                                 ds_parent(2),   ds_child(2) )
982        enddo
983    endif
984!
985    tabout = TRANSPOSE(tabout_trsp)
986!---------------------------------------------------------------------------------------------------
987end subroutine Agrif_Interp_2D_recursive
988!===================================================================================================
989!
990!===================================================================================================
991!  subroutine Agrif_Interp_3D_Recursive
992!
993!> Subroutine for the interpolation of a 3D grid variable.
994!> It calls #Agrif_Interp_2D_recursive and #Agrif_InterpBase.
995!---------------------------------------------------------------------------------------------------
996subroutine Agrif_Interp_3D_recursive ( type_interp, tabin, tabout,  &
997                                       indmin, indmax,              &
998                                       pttab_child, petab_child,    &
999                                       s_child,  s_parent,          &
1000                                       ds_child, ds_parent )
1001!---------------------------------------------------------------------------------------------------
1002    integer, dimension(3),              intent(in)  :: type_interp
1003    integer, dimension(3),              intent(in)  :: indmin, indmax
1004    integer, dimension(3),              intent(in)  :: pttab_child, petab_child
1005    real,    dimension(3),              intent(in)  :: s_child, s_parent
1006    real,    dimension(3),              intent(in)  :: ds_child, ds_parent
1007    real,    dimension(                 &
1008        indmin(1):indmax(1),            &
1009        indmin(2):indmax(2),            &
1010        indmin(3):indmax(3)),           intent(in)  :: tabin
1011    real,    dimension(                 &
1012        pttab_child(1):petab_child(1),  &
1013        pttab_child(2):petab_child(2),  &
1014        pttab_child(3):petab_child(3)), intent(out) :: tabout
1015!---------------------------------------------------------------------------------------------------
1016    real, dimension(                    &
1017        pttab_child(1):petab_child(1),  &
1018        pttab_child(2):petab_child(2),  &
1019        indmin(3):indmax(3))            :: tabtemp
1020    integer                             :: i, j, k, coeffraf
1021    integer                             :: locind_child_left, kdeb
1022!
1023    coeffraf = nint ( ds_parent(1) / ds_child(1) )
1024    if ( (type_interp(1) == Agrif_Linear) .and. (coeffraf/=1) ) then
1025        call Linear1dPrecompute2d(indmax(2)-indmin(2)+1,            &
1026                                  indmax(1)-indmin(1)+1,            &
1027                                  petab_child(1)-pttab_child(1)+1,  &
1028                                  s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1029        precomputedone(1) = .TRUE.
1030    elseif ( (type_interp(1) == Agrif_PPM) .and. (coeffraf/=1) ) then
1031        call PPM1dPrecompute2d(indmax(2)-indmin(2)+1,           &
1032                               indmax(1)-indmin(1)+1,           &
1033                               petab_child(1)-pttab_child(1)+1, &
1034                               s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1035        precomputedone(1) = .TRUE.
1036    endif
1037
1038    coeffraf = nint ( ds_parent(2) / ds_child(2) )
1039    if ( (type_interp(2) == Agrif_Linear) .and. (coeffraf/=1) ) then
1040        call Linear1dPrecompute2d(petab_child(1)-pttab_child(1)+1,  &
1041                                  indmax(2)-indmin(2)+1,            &
1042                                  petab_child(2)-pttab_child(2)+1,  &
1043                                  s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1044        precomputedone(2) = .TRUE.
1045    elseif ( (type_interp(2) == Agrif_PPM) .and. (coeffraf/=1) ) then
1046        call PPM1dPrecompute2d(petab_child(1)-pttab_child(1)+1, &
1047                               indmax(2)-indmin(2)+1,           &
1048                               petab_child(2)-pttab_child(2)+1, &
1049                               s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1050        precomputedone(2) = .TRUE.
1051    endif
1052!
1053    do k = indmin(3), indmax(3)
1054        call Agrif_Interp_2D_recursive(type_interp(1:2),                            &
1055                                       tabin(indmin(1):indmax(1),                   &
1056                                             indmin(2):indmax(2), k),               &
1057                                       tabtemp(pttab_child(1):petab_child(1),       &
1058                                               pttab_child(2):petab_child(2), k),   &
1059                                       indmin(1:2), indmax(1:2),                    &
1060                                       pttab_child(1:2), petab_child(1:2),          &
1061                                       s_child(1:2),     s_parent(1:2),             &
1062                                       ds_child(1:2),   ds_parent(1:2) )
1063    enddo
1064!
1065    precomputedone(1) = .FALSE.
1066    precomputedone(2) = .FALSE.
1067    coeffraf = nint(ds_parent(3)/ds_child(3))
1068!
1069    if ( coeffraf == 1 ) then
1070        locind_child_left = 1 + agrif_int((s_child(3)-s_parent(3))/ds_parent(3))
1071        kdeb = indmin(3)+locind_child_left-2
1072        do k = pttab_child(3),petab_child(3)
1073            kdeb = kdeb + 1
1074            do j = pttab_child(2), petab_child(2)
1075            do i = pttab_child(1), petab_child(1)
1076                tabout(i,j,k) = tabtemp(i,j,kdeb)
1077            enddo
1078            enddo
1079        enddo
1080    else
1081        do j = pttab_child(2), petab_child(2)
1082        do i = pttab_child(1), petab_child(1)
1083            call Agrif_InterpBase(type_interp(3),                                   &
1084                                  tabtemp(i,j,indmin(3):indmax(3)),                 &
1085                                  tabout(i,j,pttab_child(3):petab_child(3)),        &
1086                                  indmin(3), indmax(3),                             &
1087                                  pttab_child(3), petab_child(3),                   &
1088                                  s_parent(3),    s_child(3),                       &
1089                                 ds_parent(3),   ds_child(3) )
1090        enddo
1091        enddo
1092    endif
1093!---------------------------------------------------------------------------------------------------
1094end subroutine Agrif_Interp_3D_recursive
1095!===================================================================================================
1096!
1097!===================================================================================================
1098!  subroutine Agrif_Interp_4D_Recursive
1099!
1100!> Subroutine for the interpolation of a 4D grid variable.
1101!> It calls #Agrif_Interp_3D_recursive and #Agrif_InterpBase.
1102!---------------------------------------------------------------------------------------------------
1103subroutine Agrif_Interp_4D_recursive ( type_interp, tabin, tabout,  &
1104                                       indmin, indmax,              &
1105                                       pttab_child, petab_child,    &
1106                                       s_child,  s_parent,          &
1107                                       ds_child, ds_parent )
1108!---------------------------------------------------------------------------------------------------
1109    integer, dimension(4),              intent(in)  :: type_interp
1110    integer, dimension(4),              intent(in)  :: indmin, indmax
1111    integer, dimension(4),              intent(in)  :: pttab_child, petab_child
1112    real,    dimension(4),              intent(in)  :: s_child, s_parent
1113    real,    dimension(4),              intent(in)  :: ds_child, ds_parent
1114    real,    dimension(                 &
1115        indmin(1):indmax(1),            &
1116        indmin(2):indmax(2),            &
1117        indmin(3):indmax(3),            &
1118        indmin(4):indmax(4)),           intent(in)  :: tabin
1119    real,    dimension(                 &
1120        pttab_child(1):petab_child(1),  &
1121        pttab_child(2):petab_child(2),  &
1122        pttab_child(3):petab_child(3),  &
1123        pttab_child(4):petab_child(4)), intent(out) :: tabout
1124!---------------------------------------------------------------------------------------------------
1125    real, dimension(                    &
1126        pttab_child(1):petab_child(1),  &
1127        pttab_child(2):petab_child(2),  &
1128        pttab_child(3):petab_child(3),  &
1129        indmin(4):indmax(4))            :: tabtemp
1130    integer                             :: i, j, k, l
1131!
1132    do l = indmin(4), indmax(4)
1133        call Agrif_Interp_3D_recursive(type_interp(1:3),                            &
1134                                       tabin(indmin(1):indmax(1),                   &
1135                                             indmin(2):indmax(2),                   &
1136                                             indmin(3):indmax(3), l),               &
1137                                       tabtemp(pttab_child(1):petab_child(1),       &
1138                                               pttab_child(2):petab_child(2),       &
1139                                               pttab_child(3):petab_child(3), l),   &
1140                                       indmin(1:3), indmax(1:3),                    &
1141                                       pttab_child(1:3), petab_child(1:3),          &
1142                                       s_child(1:3),    s_parent(1:3),              &
1143                                       ds_child(1:3),  ds_parent(1:3) )
1144    enddo
1145!
1146    do k = pttab_child(3), petab_child(3)
1147    do j = pttab_child(2), petab_child(2)
1148    do i = pttab_child(1), petab_child(1)
1149        call Agrif_InterpBase(type_interp(4),                                       &
1150                              tabtemp(i,j,k,indmin(4):indmax(4)),                   &
1151                              tabout(i,j,k,pttab_child(4):petab_child(4)),          &
1152                              indmin(4), indmax(4),                                 &
1153                              pttab_child(4), petab_child(4),                       &
1154                              s_parent(4),    s_child(4),                           &
1155                             ds_parent(4),   ds_child(4) )
1156    enddo
1157    enddo
1158    enddo
1159!---------------------------------------------------------------------------------------------------
1160end subroutine Agrif_Interp_4D_recursive
1161!===================================================================================================
1162!
1163!===================================================================================================
1164!  subroutine Agrif_Interp_5D_Recursive
1165!
1166!> Subroutine for the interpolation of a 5D grid variable.
1167!> It calls #Agrif_Interp_4D_recursive and #Agrif_InterpBase.
1168!---------------------------------------------------------------------------------------------------
1169subroutine Agrif_Interp_5D_recursive ( type_interp, tabin, tabout,  &
1170                                       indmin, indmax,              &
1171                                       pttab_child, petab_child,    &
1172                                       s_child,  s_parent,          &
1173                                       ds_child, ds_parent )
1174!---------------------------------------------------------------------------------------------------
1175    integer, dimension(5),              intent(in)  :: type_interp
1176    integer, dimension(5),              intent(in)  :: indmin, indmax
1177    integer, dimension(5),              intent(in)  :: pttab_child, petab_child
1178    real,    dimension(5),              intent(in)  :: s_child, s_parent
1179    real,    dimension(5),              intent(in)  :: ds_child, ds_parent
1180    real,    dimension(                 &
1181        indmin(1):indmax(1),            &
1182        indmin(2):indmax(2),            &
1183        indmin(3):indmax(3),            &
1184        indmin(4):indmax(4),            &
1185        indmin(5):indmax(5)),           intent(in)  :: tabin
1186    real,    dimension(                 &
1187        pttab_child(1):petab_child(1),  &
1188        pttab_child(2):petab_child(2),  &
1189        pttab_child(3):petab_child(3),  &
1190        pttab_child(4):petab_child(4),  &
1191        pttab_child(5):petab_child(5)), intent(out) :: tabout
1192!---------------------------------------------------------------------------------------------------
1193    real, dimension(                    &
1194        pttab_child(1):petab_child(1),  &
1195        pttab_child(2):petab_child(2),  &
1196        pttab_child(3):petab_child(3),  &
1197        pttab_child(4):petab_child(4),  &
1198        indmin(5):indmax(5))            :: tabtemp
1199    integer                             :: i, j, k, l, m
1200!
1201    do m = indmin(5), indmax(5)
1202        call Agrif_Interp_4D_recursive(type_interp(1:4),                            &
1203                                       tabin(indmin(1):indmax(1),                   &
1204                                             indmin(2):indmax(2),                   &
1205                                             indmin(3):indmax(3),                   &
1206                                             indmin(4):indmax(4),m),                &
1207                                       tabtemp(pttab_child(1):petab_child(1),       &
1208                                               pttab_child(2):petab_child(2),       &
1209                                               pttab_child(3):petab_child(3),       &
1210                                               pttab_child(4):petab_child(4), m),   &
1211                                       indmin(1:4),indmax(1:4),                     &
1212                                       pttab_child(1:4), petab_child(1:4),          &
1213                                       s_child(1:4),     s_parent(1:4),             &
1214                                       ds_child(1:4),   ds_parent(1:4) )
1215    enddo
1216!
1217    do l = pttab_child(4), petab_child(4)
1218    do k = pttab_child(3), petab_child(3)
1219    do j = pttab_child(2), petab_child(2)
1220    do i = pttab_child(1), petab_child(1)
1221        call Agrif_InterpBase(type_interp(5),                                       &
1222                              tabtemp(i,j,k,l,indmin(5):indmax(5)),                 &
1223                              tabout(i,j,k,l,pttab_child(5):petab_child(5)),        &
1224                              indmin(5), indmax(5),                                 &
1225                              pttab_child(5), petab_child(5),                       &
1226                              s_parent(5),   s_child(5),                            &
1227                              ds_parent(5), ds_child(5) )
1228    enddo
1229    enddo
1230    enddo
1231    enddo
1232!---------------------------------------------------------------------------------------------------
1233end subroutine Agrif_Interp_5D_recursive
1234!===================================================================================================
1235!
1236!===================================================================================================
1237!  subroutine Agrif_Interp_6D_Recursive
1238!
1239!> Subroutine for the interpolation of a 6D grid variable.
1240!> It calls #Agrif_Interp_5D_recursive and Agrif_InterpBase.
1241!---------------------------------------------------------------------------------------------------
1242subroutine Agrif_Interp_6D_recursive ( type_interp, tabin, tabout,  &
1243                                       indmin, indmax,              &
1244                                       pttab_child, petab_child,    &
1245                                       s_child,  s_parent,          &
1246                                       ds_child, ds_parent )
1247!---------------------------------------------------------------------------------------------------
1248    integer, dimension(6),                  intent(in)  :: type_interp
1249    integer, dimension(6),                  intent(in)  :: indmin, indmax
1250    integer, dimension(6),                  intent(in)  :: pttab_child, petab_child
1251    real,    dimension(6),                  intent(in)  :: s_child, s_parent
1252    real,    dimension(6),                  intent(in)  :: ds_child, ds_parent
1253    real,    dimension(                 &
1254        indmin(1):indmax(1),            &
1255        indmin(2):indmax(2),            &
1256        indmin(3):indmax(3),            &
1257        indmin(4):indmax(4),            &
1258        indmin(5):indmax(5),            &
1259        indmin(6):indmax(6)),               intent(in)  :: tabin
1260    real,    dimension(                 &
1261        pttab_child(1):petab_child(1),  &
1262        pttab_child(2):petab_child(2),  &
1263        pttab_child(3):petab_child(3),  &
1264        pttab_child(4):petab_child(4),  &
1265        pttab_child(5):petab_child(5),  &
1266        pttab_child(6):petab_child(6)),     intent(out) :: tabout
1267!---------------------------------------------------------------------------------------------------
1268    real, dimension(                    &
1269        pttab_child(1):petab_child(1),  &
1270        pttab_child(2):petab_child(2),  &
1271        pttab_child(3):petab_child(3),  &
1272        pttab_child(4):petab_child(4),  &
1273        pttab_child(5):petab_child(5),  &
1274        indmin(6):indmax(6))            :: tabtemp
1275    integer                             :: i, j, k, l, m, n
1276!
1277    do n = indmin(6), indmax(6)
1278        call Agrif_Interp_5D_recursive(type_interp(1:5),                            &
1279                                       tabin(indmin(1):indmax(1),                   &
1280                                             indmin(2):indmax(2),                   &
1281                                             indmin(3):indmax(3),                   &
1282                                             indmin(4):indmax(4),                   &
1283                                             indmin(5):indmax(5), n),               &
1284                                        tabtemp(pttab_child(1):petab_child(1),      &
1285                                                pttab_child(2):petab_child(2),      &
1286                                                pttab_child(3):petab_child(3),      &
1287                                                pttab_child(4):petab_child(4),      &
1288                                                pttab_child(5):petab_child(5), n),  &
1289                                        indmin(1:5),indmax(1:5),                    &
1290                                        pttab_child(1:5), petab_child(1:5),         &
1291                                        s_child(1:5), s_parent(1:5),                &
1292                                        ds_child(1:5),ds_parent(1:5) )
1293    enddo
1294!
1295    do m = pttab_child(5), petab_child(5)
1296    do l = pttab_child(4), petab_child(4)
1297    do k = pttab_child(3), petab_child(3)
1298    do j = pttab_child(2), petab_child(2)
1299    do i = pttab_child(1), petab_child(1)
1300        call Agrif_InterpBase(type_interp(6),                                       &
1301                              tabtemp(i,j,k,l,m,indmin(6):indmax(6)),               &
1302                              tabout(i,j,k,l,m,pttab_child(6):petab_child(6)),      &
1303                              indmin(6), indmax(6),                                 &
1304                              pttab_child(6), petab_child(6),                       &
1305                              s_parent(6),   s_child(6),                            &
1306                              ds_parent(6), ds_child(6) )
1307    enddo
1308    enddo
1309    enddo
1310    enddo
1311    enddo
1312!---------------------------------------------------------------------------------------------------
1313end subroutine Agrif_Interp_6D_recursive
1314!===================================================================================================
1315!
1316!===================================================================================================
1317!  subroutine Agrif_InterpBase
1318!
1319!> Calls the interpolation method chosen by the user (linear, lagrange, spline, etc.).
1320!---------------------------------------------------------------------------------------------------
1321subroutine Agrif_InterpBase ( type_interp, parenttab, childtab, indmin, indmax, &
1322                              pttab_child, petab_child,                         &
1323                              s_parent, s_child, ds_parent, ds_child )
1324!---------------------------------------------------------------------------------------------------
1325    INTEGER                                                 :: type_interp
1326    INTEGER                                                 :: indmin, indmax
1327    INTEGER                                                 :: pttab_child, petab_child
1328    REAL, DIMENSION(indmin:indmax),           INTENT(IN)    :: parenttab
1329    REAL, DIMENSION(pttab_child:petab_child), INTENT(OUT)   :: childtab
1330    REAL                                                    :: s_parent, s_child
1331    REAL                                                    :: ds_parent,ds_child
1332!
1333    if ( (indmin == indmax) .and. (pttab_child == petab_child) ) then
1334!
1335        childtab(pttab_child) = parenttab(indmin)
1336!
1337    elseif (type_interp == Agrif_LINEAR) then    !       Linear interpolation
1338!
1339        call Agrif_basicinterp_linear1D(parenttab,childtab,                   &
1340                indmax-indmin+1,petab_child-pttab_child+1,  &
1341                s_parent,s_child,ds_parent,ds_child)
1342!
1343    elseif ( type_interp == Agrif_PPM ) then     !       PPM interpolation
1344
1345        call PPM1d(parenttab,childtab,                      &
1346                indmax-indmin+1,petab_child-pttab_child+1,  &
1347                s_parent,s_child,ds_parent,ds_child)
1348!
1349    elseif ( type_interp == Agrif_PPM_LIM ) then     !       PPM interpolation
1350
1351        call PPM1d_lim(parenttab,childtab,                      &
1352                indmax-indmin+1,petab_child-pttab_child+1,  &
1353                s_parent,s_child,ds_parent,ds_child)
1354!
1355    elseif (type_interp == Agrif_LAGRANGE) then  !       Lagrange interpolation
1356!
1357        call lagrange1D(parenttab,childtab,                 &
1358                indmax-indmin+1,petab_child-pttab_child+1,  &
1359                s_parent,s_child,ds_parent,ds_child)
1360!
1361    elseif (type_interp == Agrif_ENO) then       !       Eno interpolation
1362!
1363        call ENO1d(parenttab,childtab,                      &
1364                indmax-indmin+1,petab_child-pttab_child+1,  &
1365                s_parent,s_child,ds_parent,ds_child)
1366!
1367    elseif (type_interp == Agrif_WENO) then      !       Weno interpolation
1368!
1369        call WENO1d(parenttab,childtab,                     &
1370                indmax-indmin+1,petab_child-pttab_child+1,  &
1371                s_parent,s_child,ds_parent,ds_child)
1372!
1373    elseif (type_interp == Agrif_LINEARCONSERV) then !   Linear conservative interpolation
1374!
1375        call Linear1dConserv(parenttab,childtab,            &
1376                indmax-indmin+1,petab_child-pttab_child+1,  &
1377                s_parent,s_child,ds_parent,ds_child)
1378!
1379    elseif (type_interp == Agrif_LINEARCONSERVLIM) then !Linear conservative interpolation
1380!
1381        call Linear1dConservLim(parenttab,childtab,         &
1382                indmax-indmin+1,petab_child-pttab_child+1,  &
1383                s_parent,s_child,ds_parent,ds_child)
1384!
1385    elseif (type_interp == Agrif_CONSTANT) then
1386!
1387        call Constant1d(parenttab,childtab,                 &
1388                indmax-indmin+1,petab_child-pttab_child+1,  &
1389                s_parent,s_child,ds_parent,ds_child)
1390!
1391    endif
1392!---------------------------------------------------------------------------------------------------
1393end subroutine Agrif_InterpBase
1394!===================================================================================================
1395!
1396!===================================================================================================
1397!  subroutine Agrif_Find_list_interp
1398!---------------------------------------------------------------------------------------------------
1399function Agrif_Find_list_interp ( list_interp, pttab, petab, pttab_Child, pttab_Parent,     &
1400                                    nbdim, indmin, indmax, indminglob,  indmaxglob,         &
1401                                    pttruetab, cetruetab, memberin                          &
1402#if defined AGRIF_MPI
1403                                   ,indminglob2, indmaxglob2, parentarray,                  &
1404                                    member, tab4t, memberinall, sendtoproc1, recvfromproc1  &
1405#endif
1406    ) result(find_list_interp)
1407!---------------------------------------------------------------------------------------------------
1408    type(Agrif_List_Interp_Loc),   pointer     :: list_interp
1409    integer,                       intent(in)  :: nbdim
1410    integer, dimension(nbdim),     intent(in)  :: pttab, petab, pttab_Child, pttab_Parent
1411    integer, dimension(nbdim),     intent(out) :: indmin, indmax
1412    integer, dimension(nbdim),     intent(out) :: indminglob, indmaxglob
1413    integer, dimension(nbdim),     intent(out) :: pttruetab, cetruetab
1414    logical,                       intent(out) :: memberin
1415#if defined AGRIF_MPI
1416    integer, dimension(nbdim),     intent(out) :: indminglob2, indmaxglob2
1417    integer, dimension(nbdim,2,2), intent(out) :: parentarray
1418    logical,                       intent(out) :: member
1419    integer, dimension(nbdim,0:Agrif_Nbprocs-1,8), intent(out) :: tab4t
1420    logical, dimension(0:Agrif_Nbprocs-1),         intent(out) :: memberinall
1421    logical, dimension(0:Agrif_Nbprocs-1),         intent(out) :: sendtoproc1, recvfromproc1
1422#endif
1423    logical :: find_list_interp
1424!
1425    integer :: i
1426    type(Agrif_List_Interp_Loc), pointer :: parcours
1427    type(Agrif_Interp_Loc),      pointer :: pil
1428
1429    find_list_interp = .false.
1430
1431    if ( .not. associated(list_interp) )    return
1432
1433    parcours => list_interp
1434    find_loop : do while ( associated(parcours) )
1435
1436        pil => parcours % interp_loc
1437
1438        do i = 1,nbdim
1439            if ( (pttab(i) /= pil % pttab(i)) .or. &
1440                 (petab(i) /= pil % petab(i)) .or. &
1441                 (pttab_child(i)  /= pil % pttab_child(i)) .or. &
1442                 (pttab_parent(i) /= pil % pttab_parent(i)) ) then
1443                parcours => parcours % suiv
1444                cycle find_loop
1445            endif
1446        enddo
1447
1448        indmin = pil % indmin(1:nbdim)
1449        indmax = pil % indmax(1:nbdim)
1450
1451        pttruetab = pil % pttruetab(1:nbdim)
1452        cetruetab = pil % cetruetab(1:nbdim)
1453
1454#if !defined AGRIF_MPI
1455        indminglob  = pil % indminglob(1:nbdim)
1456        indmaxglob  = pil % indmaxglob(1:nbdim)
1457#else
1458        indminglob  = pil % indminglob2(1:nbdim)
1459        indmaxglob  = pil % indmaxglob2(1:nbdim)
1460        indminglob2 = pil % indminglob2(1:nbdim)
1461        indmaxglob2 = pil % indmaxglob2(1:nbdim)
1462        parentarray = pil % parentarray(1:nbdim,:,:)
1463        member      = pil % member
1464        tab4t         = pil % tab4t(1:nbdim, 0:Agrif_Nbprocs-1, 1:8)
1465        memberinall   = pil % memberinall(0:Agrif_Nbprocs-1)
1466        sendtoproc1   = pil % sendtoproc1(0:Agrif_Nbprocs-1)
1467        recvfromproc1 = pil % recvfromproc1(0:Agrif_Nbprocs-1)
1468#endif
1469        memberin = pil % memberin
1470        find_list_interp = .true.
1471        exit find_loop
1472    enddo find_loop
1473!---------------------------------------------------------------------------------------------------
1474end function Agrif_Find_list_interp
1475!===================================================================================================
1476!
1477!===================================================================================================
1478!  subroutine Agrif_AddTo_list_interp
1479!---------------------------------------------------------------------------------------------------
1480subroutine Agrif_AddTo_list_interp ( list_interp, pttab, petab, pttab_Child, pttab_Parent,  &
1481                                     indmin, indmax, indminglob, indmaxglob,                &
1482                                     pttruetab, cetruetab,                                  &
1483                                     memberin, nbdim                                        &
1484#if defined AGRIF_MPI
1485                                    ,indminglob2, indmaxglob2,                              &
1486                                     parentarray,                                           &
1487                                     member,                                                &
1488                                     tab4t, memberinall, sendtoproc1, recvfromproc1         &
1489#endif
1490    )
1491!---------------------------------------------------------------------------------------------------
1492    type(Agrif_List_Interp_Loc), pointer    :: list_interp
1493    integer                                 :: nbdim
1494    integer, dimension(nbdim)               :: pttab, petab, pttab_Child, pttab_Parent
1495    integer, dimension(nbdim)               :: indmin,indmax
1496    integer, dimension(nbdim)               :: indminglob, indmaxglob
1497    integer, dimension(nbdim)               :: pttruetab, cetruetab
1498    logical                                 :: memberin
1499#if defined AGRIF_MPI
1500    integer, dimension(nbdim,2,2)           :: parentarray
1501    logical                                 :: member
1502    integer, dimension(nbdim)                       :: indminglob2,indmaxglob2
1503    integer, dimension(nbdim,0:Agrif_Nbprocs-1,8)   :: tab4t
1504    logical, dimension(0:Agrif_Nbprocs-1)           :: memberinall
1505    logical, dimension(0:Agrif_Nbprocs-1)           :: sendtoproc1
1506    logical, dimension(0:Agrif_Nbprocs-1)           :: recvfromproc1
1507#endif
1508!
1509    type(Agrif_List_Interp_Loc), pointer    :: parcours
1510    type(Agrif_Interp_Loc),      pointer    :: pil
1511!
1512    allocate(parcours)
1513    allocate(parcours % interp_loc)
1514
1515    pil => parcours % interp_loc
1516
1517    pil % pttab(1:nbdim) = pttab(1:nbdim)
1518    pil % petab(1:nbdim) = petab(1:nbdim)
1519    pil % pttab_child(1:nbdim) = pttab_child(1:nbdim)
1520    pil % pttab_parent(1:nbdim) = pttab_parent(1:nbdim)
1521
1522    pil % indmin(1:nbdim) = indmin(1:nbdim)
1523    pil % indmax(1:nbdim) = indmax(1:nbdim)
1524
1525    pil % memberin = memberin
1526#if !defined AGRIF_MPI
1527    pil % indminglob(1:nbdim) = indminglob(1:nbdim)
1528    pil % indmaxglob(1:nbdim) = indmaxglob(1:nbdim)
1529#else
1530    pil % indminglob2(1:nbdim) = indminglob2(1:nbdim)
1531    pil % indmaxglob2(1:nbdim) = indmaxglob2(1:nbdim)
1532    pil % parentarray(1:nbdim,:,:) = parentarray(1:nbdim,:,:)
1533    pil % member = member
1534    allocate(pil % tab4t(nbdim, 0:Agrif_Nbprocs-1, 8))
1535    allocate(pil % memberinall(0:Agrif_Nbprocs-1))
1536    allocate(pil % sendtoproc1(0:Agrif_Nbprocs-1))
1537    allocate(pil % recvfromproc1(0:Agrif_Nbprocs-1))
1538    pil % tab4t         = tab4t
1539    pil % memberinall   = memberinall
1540    pil % sendtoproc1   = sendtoproc1
1541    pil % recvfromproc1 = recvfromproc1
1542#endif
1543
1544    pil % pttruetab(1:nbdim) = pttruetab(1:nbdim)
1545    pil % cetruetab(1:nbdim) = cetruetab(1:nbdim)
1546
1547    parcours % suiv => list_interp
1548    list_interp => parcours
1549!---------------------------------------------------------------------------------------------------
1550end subroutine Agrif_Addto_list_interp
1551!===================================================================================================
1552!
1553end module Agrif_Interpolation
Note: See TracBrowser for help on using the repository browser.