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.F in trunk/AGRIF/AGRIF_FILES – NEMO

source: trunk/AGRIF/AGRIF_FILES/modinterp.F @ 898

Last change on this file since 898 was 898, checked in by rblod, 16 years ago

Correct some bugs in agrif optimization and add MPP optimization, see #42

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