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 @ 779

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

Agrif improvment for vectorization, see ticket #41

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 77.6 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,4) :: tab4t
634      LOGICAL, DIMENSION(0:Agrif_Nbprocs-1) :: memberinall
635      LOGICAL, DIMENSION(1) :: memberin1
636C
637#endif     
638C     
639C   
640C     Boundaries of the current grid where interpolation is done
641
642     
643     
644     
645      IF (Associated(child%var%list_interp)) THEN
646      Call Agrif_Find_list_interp(child%var%list_interp,pttab,petab,
647     &                          pttab_Child,pttab_Parent,nbdim,
648     &                          indmin,indmax,indminglob,
649     &       indmaxglob,indminglob2,indmaxglob2,parentarray,
650     &       pttruetab,cetruetab,member,memberin,find_list_interp
651#if defined AGRIF_MPI
652     &       ,tab4t,memberinall
653#endif
654     &    )
655      ELSE
656      find_list_interp = .FALSE.
657      ENDIF
658     
659      IF (.not.find_list_interp) THEN
660      Call Agrif_nbdim_Get_bound_dimension(child % var,
661     &                               lowerbound,upperbound,nbdim)
662
663      Call Agrif_Childbounds(nbdim,lowerbound,upperbound,
664     &                                   pttab,petab,
665     &                                   pttruetab,cetruetab,memberin)
666     
667C
668C
669
670      Call Agrif_Parentbounds(TYPEinterp,nbdim,indminglob,indmaxglob,
671     &                        s_Parent_temp,s_Child_temp,
672     &                        s_Child,ds_Child,
673     &                        s_Parent,ds_Parent,
674     &                        pttab,petab,
675     &                        pttab_Child,pttab_Parent,
676     &                        child%var%root_var%posvar,
677     &                        child % var % root_var % interptab)
678
679
680#ifdef AGRIF_MPI
681       IF (memberin) THEN
682        Call Agrif_Parentbounds(TYPEinterp,nbdim,indmin,indmax,
683     &                        s_Parent_temp,s_Child_temp,
684     &                        s_Child,ds_Child,
685     &                        s_Parent,ds_Parent,
686     &                        pttruetab,cetruetab,
687     &                        pttab_Child,pttab_Parent,
688     &                        child%var%root_var%posvar,
689     &                        child % var % root_var % interptab)
690      ENDIF
691
692      Call Agrif_nbdim_Get_bound_dimension(parent%var,
693     &                              lowerbound,upperbound,nbdim)
694
695      Call Agrif_ChildGrid_to_ParentGrid()
696C
697      Call Agrif_Childbounds(nbdim,
698     &                       lowerbound,upperbound,
699     &                       indminglob,indmaxglob,
700     &                       indminglob2,indmaxglob2,member)
701
702C
703      IF (member) THEN
704      Call Agrif_GlobtoLocInd2(parentarray,
705     &                     lowerbound,upperbound,
706     &                     indminglob2,indmaxglob2,
707     &                     nbdim,Agrif_Procrank,
708     &                     member)
709       endif
710       
711      Call Agrif_ParentGrid_to_ChildGrid()
712#else
713      parentarray(:,1,1) = indminglob
714      parentarray(:,2,1) = indmaxglob
715      parentarray(:,1,2) = indminglob
716      parentarray(:,2,2) = indmaxglob
717      indmin = indminglob
718      indmax = indmaxglob
719      member = .TRUE.
720#endif
721
722      ELSE
723     
724#if !defined AGRIF_MPI
725      parentarray(:,1,1) = indminglob
726      parentarray(:,2,1) = indmaxglob
727      parentarray(:,1,2) = indminglob
728      parentarray(:,2,2) = indmaxglob
729      indmin = indminglob
730      indmax = indmaxglob
731      member = .TRUE.
732      s_Parent_temp = s_Parent + (indminglob - pttab_Parent)*ds_Parent
733      s_Child_temp = s_Child + (pttab - pttab_Child) * ds_Child
734#else     
735      s_Parent_temp = s_Parent + (indmin - pttab_Parent)*ds_Parent
736      s_Child_temp = s_Child + (pttruetab - pttab_Child) * ds_Child
737#endif     
738           
739      ENDIF
740
741
742
743      IF (member) THEN
744      IF (.not.associated(tempP%var)) allocate(tempP%var)
745
746C
747      Call Agrif_nbdim_allocation(tempP%var,
748     &     parentarray(:,1,1),parentarray(:,2,1),nbdim)
749
750      Call Agrif_nbdim_Full_VarEQreal(tempP%var,0.,nbdim)
751
752
753
754      IF (present(procname)) THEN
755      Call Agrif_ChildGrid_to_ParentGrid()
756            SELECT CASE (nbdim)
757        CASE(1)
758          CALL procname(tempP%var%array1,
759     &                          parentarray(1,1,2),parentarray(1,2,2))
760        CASE(2)
761          CALL procname(tempP%var%array2,
762     &                          parentarray(1,1,2),parentarray(1,2,2),
763     &                          parentarray(2,1,2),parentarray(2,2,2))
764        CASE(3)
765          CALL procname(tempP%var%array3,
766     &                          parentarray(1,1,2),parentarray(1,2,2),
767     &                          parentarray(2,1,2),parentarray(2,2,2),
768     &                          parentarray(3,1,2),parentarray(3,2,2))
769        CASE(4)
770          CALL procname(tempP%var%array4,
771     &                          parentarray(1,1,2),parentarray(1,2,2),
772     &                          parentarray(2,1,2),parentarray(2,2,2),
773     &                          parentarray(3,1,2),parentarray(3,2,2),
774     &                          parentarray(4,1,2),parentarray(4,2,2))
775        CASE(5)
776          CALL procname(tempP%var%array5,
777     &                          parentarray(1,1,2),parentarray(1,2,2),
778     &                          parentarray(2,1,2),parentarray(2,2,2),
779     &                          parentarray(3,1,2),parentarray(3,2,2),
780     &                          parentarray(4,1,2),parentarray(4,2,2),
781     &                          parentarray(5,1,2),parentarray(5,2,2))
782        CASE(6)
783          CALL procname(tempP%var%array6,
784     &                          parentarray(1,1,2),parentarray(1,2,2),
785     &                          parentarray(2,1,2),parentarray(2,2,2),
786     &                          parentarray(3,1,2),parentarray(3,2,2),
787     &                          parentarray(4,1,2),parentarray(4,2,2),
788     &                          parentarray(5,1,2),parentarray(5,2,2),
789     &                          parentarray(6,1,2),parentarray(6,2,2))
790            END SELECT
791      Call Agrif_ParentGrid_to_ChildGrid()
792      ELSE
793
794      Call Agrif_nbdim_VarEQvar(tempP%var,
795     &        parentarray(:,1,1),parentarray(:,2,1),
796     &        parent%var,parentarray(:,1,2),parentarray(:,2,2),
797     &        nbdim)
798      ENDIF
799            endif
800
801#ifdef AGRIF_MPI
802      if (.not.find_list_interp) then
803      tab3(:,1) = indminglob2(:)
804      tab3(:,2) = indmaxglob2(:)
805      tab3(:,3) = indmin(:)
806      tab3(:,4) = indmax(:)
807C
808C
809      Call MPI_ALLGATHER(tab3,4*nbdim,MPI_INTEGER,tab4,4*nbdim,
810     &                   MPI_INTEGER,MPI_COMM_WORLD,code)
811
812      IF (.not.associated(tempPextend%var)) Allocate(tempPextend%var)
813
814      DO k=0,Agrif_Nbprocs-1
815       do j=1,4
816         do i=1,nbdim
817         tab4t(i,k,j) = tab4(i,j,k)
818                enddo
819      enddo
820      enddo
821     
822      memberin1(1) = memberin
823      CALL MPI_ALLGATHER(memberin1,1,MPI_LOGICAL,memberinall,
824     &                  1,MPI_LOGICAL,MPI_COMM_WORLD,code)
825     
826      endif     
827           
828      Call Get_External_Data(tempP,tempPextend,tab4t(:,:,1),
829     &            tab4t(:,:,2),
830     &            tab4t(:,:,3),tab4t(:,:,4),nbdim,member,memberin,
831     &            memberinall)
832#else
833      tempPextend%var => tempP%var
834#endif
835
836      if (.not.find_list_interp) then
837      Call Agrif_Addto_list_interp(child%var%list_interp,pttab,petab,
838     &                          pttab_Child,pttab_Parent,indmin,indmax,
839     &   indminglob,indmaxglob,indminglob2,indmaxglob2,parentarray,
840     &   pttruetab,cetruetab,member,memberin,nbdim
841#if defined AGRIF_MPI
842     &   ,tab4t,memberinall
843#endif
844     &    )
845      endif
846     
847C
848C
849      IF (memberin) THEN
850      IF (.not.associated(tempC%var)) allocate(tempC%var)
851C
852
853      Call Agrif_nbdim_allocation(tempC%var,pttruetab,cetruetab,nbdim)
854
855C
856C
857C     Special values on the parent grid
858      if (Agrif_UseSpecialValue) then
859C
860          noraftab(1:nbdim) =
861     &         child % var % root_var % interptab(1:nbdim) .EQ. 'N'
862C
863          IF (.not.associated(parentvalues%var))
864     &            Allocate(parentvalues%var)
865C
866          Call Agrif_nbdim_allocation
867     &               (parentvalues%var,indmin,indmax,nbdim)
868          Call Agrif_nbdim_Full_VarEQvar
869     &               (parentvalues%var,tempPextend%var,nbdim)
870C
871          Call Agrif_CheckMasknD(tempPextend,
872     &                           parentvalues,
873     &                           indmin(1:nbdim),indmax(1:nbdim),
874     &                           indmin(1:nbdim),indmax(1:nbdim),
875     &                           noraftab(1:nbdim),nbdim)
876C
877          Call Agrif_nbdim_deallocation(parentvalues%var,nbdim)
878C          Deallocate(parentvalues%var)
879C
880C
881      endif     
882
883C
884C
885C     Interpolation of the current grid
886
887      IF (memberin) THEN
888      if ( nbdim .EQ. 1 ) then
889         Call Agrif_Interp_1D_recursive(TypeInterp,
890     &           tempPextend%var%array1,tempC%var%array1,
891     &           indmin,indmax,
892     &           pttruetab,cetruetab,
893     &           s_Child_temp,s_Parent_temp,
894     &           ds_Child,ds_Parent,nbdim)
895      elseif ( nbdim .EQ. 2 ) then
896
897         Call Agrif_Interp_2D_recursive(TypeInterp,
898     &           tempPextend%var%array2,tempC%var%array2,
899     &           indmin,indmax,
900     &           pttruetab,cetruetab,
901     &           s_Child_temp,s_Parent_temp,
902     &           ds_Child,ds_Parent,nbdim)
903      elseif ( nbdim .EQ. 3 ) then
904
905         Call Agrif_Interp_3D_recursive(TypeInterp,
906     &           tempPextend%var%array3,tempC%var%array3,
907     &           indmin,indmax,
908     &           pttruetab,cetruetab,
909     &           s_Child_temp,s_Parent_temp,
910     &           ds_Child,ds_Parent,nbdim)
911      elseif ( nbdim .EQ. 4 ) then
912         Call Agrif_Interp_4D_recursive(TypeInterp,
913     &           tempPextend%var%array4,tempC%var%array4,
914     &           indmin,indmax,
915     &           pttruetab,cetruetab,
916     &           s_Child_temp,s_Parent_temp,
917     &           ds_Child,ds_Parent,nbdim)
918      elseif ( nbdim .EQ. 5 ) then
919         Call Agrif_Interp_5D_recursive(TypeInterp,
920     &           tempPextend%var%array5,tempC%var%array5,
921     &           indmin,indmax,
922     &           pttruetab,cetruetab,
923     &           s_Child_temp,s_Parent_temp,
924     &           ds_Child,ds_Parent,nbdim)
925      elseif ( nbdim .EQ. 6 ) then
926         Call Agrif_Interp_6D_recursive(TypeInterp,
927     &           tempPextend%var%array6,tempC%var%array6,
928     &           indmin,indmax,
929     &           pttruetab,cetruetab,
930     &           s_Child_temp,s_Parent_temp,
931     &           ds_Child,ds_Parent,nbdim)
932       endif
933
934
935C
936C
937C     Special values on the child grid 
938      if (Agrif_UseSpecialValueFineGrid) then
939C
940#ifdef AGRIF_MPI
941C       
942          Call GiveAgrif_SpecialValueToTab_mpi(child%var,tempC%var,
943     &                 childarray,
944     &                 pttruetab,cetruetab,
945     &                 Agrif_SpecialValueFineGrid,nbdim)
946C
947#else
948C
949          Call GiveAgrif_SpecialValueToTab(child%var,tempC%var,
950     &                  pttruetab,cetruetab,
951     &                  Agrif_SpecialValueFineGrid,nbdim)
952C
953#endif
954C
955C       
956      endif
957C 
958
959      Call Agrif_nbdim_Get_bound_dimension(child % var,
960     &                               lowerbound,upperbound,nbdim)
961
962#ifdef AGRIF_MPI
963      Call Agrif_GlobtoLocInd2(childarray,
964     &                     lowerbound,upperbound,
965     &                     pttruetab,cetruetab,
966     &                     nbdim,Agrif_Procrank,
967     &                     memberout)
968
969#else
970       childarray(:,1,1) = pttruetab
971       childarray(:,2,1) = cetruetab
972       childarray(:,1,2) = pttruetab
973       childarray(:,2,2) = cetruetab
974ccccccccccccccc       memberout = .TRUE.
975#endif
976
977      endif
978
979C
980      if (torestore) then
981C
982#ifdef AGRIF_MPI
983C
984        SELECT CASE (nbdim)
985        CASE (1)
986             do i = pttruetab(1),cetruetab(1)         
987ChildarrayAModifier                if (restore%var%restore1D(i) == 0) 
988ChildarrayAModifier     &                child%var%array1(childarray(i,1,2)
989ChildarrayAModifier     &                                 ) = 
990ChildarrayAModifier     &                tempC%var%array1(i)
991             enddo
992        CASE (2)
993             do i = pttruetab(1),cetruetab(1)
994             do j = pttruetab(2),cetruetab(2)
995ChildarrayAModifier                   if (restore%var%restore2D(i,j) == 0) 
996ChildarrayAModifier     &                child%var%array2(childarray(i,1,2),
997ChildarrayAModifier     &                                 childarray(j,2,2)) = 
998ChildarrayAModifier     &                tempC%var%array2(i,j)
999             enddo
1000             enddo
1001        CASE (3)
1002             do i = pttruetab(1),cetruetab(1)
1003             do j = pttruetab(2),cetruetab(2) 
1004             do k = pttruetab(3),cetruetab(3)
1005ChildarrayAModifier                      if (restore%var%restore3D(i,j,k) == 0) 
1006ChildarrayAModifier     &                child%var%array3(childarray(i,1,2),
1007ChildarrayAModifier     &                                 childarray(j,2,2),
1008ChildarrayAModifier     &                                 childarray(k,3,2)) = 
1009ChildarrayAModifier     &                tempC%var%array3(i,j,k)
1010             enddo
1011             enddo
1012             enddo
1013        CASE (4)
1014             do i = pttruetab(1),cetruetab(1)
1015             do j = pttruetab(2),cetruetab(2)
1016             do k = pttruetab(3),cetruetab(3)
1017             do l = pttruetab(4),cetruetab(4)
1018ChildarrayAModifier                         if (restore%var%restore4D(i,j,k,l) == 0) 
1019ChildarrayAModifier     &                      child%var%array4(childarray(i,1,2),
1020ChildarrayAModifier     &                                       childarray(j,2,2),
1021ChildarrayAModifier     &                                       childarray(k,3,2),
1022ChildarrayAModifier     &                                       childarray(l,4,2)) = 
1023ChildarrayAModifier     &                      tempC%var%array4(i,j,k,l)
1024             enddo
1025             enddo
1026             enddo
1027             enddo
1028        CASE (5)
1029             do i = pttruetab(1),cetruetab(1)
1030             do j = pttruetab(2),cetruetab(2)
1031             do k = pttruetab(3),cetruetab(3)
1032             do l = pttruetab(4),cetruetab(4)
1033             do m = pttruetab(5),cetruetab(5)
1034ChildarrayAModifier              if (restore%var%restore5D(i,j,k,l,m) == 0) 
1035ChildarrayAModifier     &                child%var%array5(childarray(i,1,2),
1036ChildarrayAModifier     &                                 childarray(j,2,2),
1037ChildarrayAModifier     &                                 childarray(k,3,2),
1038ChildarrayAModifier     &                                 childarray(l,4,2),
1039ChildarrayAModifier     &                                 childarray(m,5,2)) = 
1040ChildarrayAModifier     &                tempC%var%array5(i,j,k,l,m)
1041             enddo
1042             enddo
1043             enddo
1044             enddo
1045             enddo
1046        CASE (6)
1047             do i = pttruetab(1),cetruetab(1)
1048             do j = pttruetab(2),cetruetab(2)
1049             do k = pttruetab(3),cetruetab(3)
1050             do l = pttruetab(4),cetruetab(4)
1051             do m = pttruetab(5),cetruetab(5)
1052             do n = pttruetab(6),cetruetab(6)
1053ChildarrayAModifier              if (restore%var%restore6D(i,j,k,l,m,n) == 0) 
1054ChildarrayAModifier     &                child%var%array6(childarray(i,1,2),
1055ChildarrayAModifier     &                                 childarray(j,2,2),
1056ChildarrayAModifier     &                                 childarray(k,3,2),
1057ChildarrayAModifier     &                                 childarray(l,4,2),
1058ChildarrayAModifier     &                                 childarray(m,5,2),
1059ChildarrayAModifier     &                                 childarray(n,6,2)) = 
1060ChildarrayAModifier     &                tempC%var%array6(i,j,k,l,m,n)
1061             enddo
1062             enddo
1063             enddo
1064             enddo
1065             enddo
1066             enddo
1067        END SELECT
1068C
1069#else
1070        SELECT CASE (nbdim)
1071        CASE (1)
1072           do i = pttruetab(1),cetruetab(1)         
1073            if (restore%var%restore1D(i) == 0)
1074     &            child % var % array1(i) = 
1075     &            tempC % var % array1(i)   
1076          enddo
1077        CASE (2)
1078           do j = pttruetab(2),cetruetab(2)
1079             do i = pttruetab(1),cetruetab(1)           
1080              if (restore%var%restore2D(i,j) == 0)     
1081     &              child % var % array2(i,j) = 
1082     &              tempC % var % array2(i,j)   
1083              enddo
1084             enddo
1085        CASE (3)
1086           do k = pttruetab(3),cetruetab(3)
1087           do j = pttruetab(2),cetruetab(2)
1088             do i = pttruetab(1),cetruetab(1) 
1089              if (restore%var%restore3D(i,j,k) == 0)
1090     &                  child % var % array3(i,j,k) =
1091     &                  tempC % var % array3(i,j,k)   
1092                  enddo
1093              enddo
1094             enddo
1095        CASE (4)
1096           do l = pttruetab(4),cetruetab(4)
1097           do k = pttruetab(3),cetruetab(3)
1098          do j = pttruetab(2),cetruetab(2)
1099             do i = pttruetab(1),cetruetab(1)
1100                if (restore%var%restore4D(i,j,k,l) == 0)
1101     &                 child % var % array4(i,j,k,l) = 
1102     &                 tempC % var % array4(i,j,k,l)   
1103             enddo
1104             enddo
1105              enddo
1106             enddo
1107        CASE (5)
1108           do m = pttruetab(5),cetruetab(5)
1109          do l = pttruetab(4),cetruetab(4)
1110         do k = pttruetab(3),cetruetab(3)
1111           do j = pttruetab(2),cetruetab(2)
1112             do i = pttruetab(1),cetruetab(1)
1113                if (restore%var%restore5D(i,j,k,l,m) == 0)
1114     &                  child % var % array5(i,j,k,l,m) = 
1115     &                  tempC % var % array5(i,j,k,l,m)   
1116             enddo
1117             enddo
1118                  enddo
1119              enddo
1120             enddo
1121        CASE (6)
1122           do n = pttruetab(6),cetruetab(6)
1123          do m = pttruetab(5),cetruetab(5)
1124          do l = pttruetab(4),cetruetab(4)
1125         do k = pttruetab(3),cetruetab(3)
1126          do j = pttruetab(2),cetruetab(2)
1127             do i = pttruetab(1),cetruetab(1)
1128                if (restore%var%restore6D(i,j,k,l,m,n) == 0)
1129     &                      child % var % array6(i,j,k,l,m,n) = 
1130     &                      tempC % var % array6(i,j,k,l,m,n)   
1131             enddo
1132            enddo
1133                      enddo
1134                  enddo
1135              enddo
1136             enddo
1137        END SELECT
1138C
1139#endif
1140C       
1141        else
1142C
1143C
1144          IF (memberin) THEN
1145          SELECT CASE (nbdim)
1146          CASE (1)
1147            child%var%array1(childarray(1,1,2):childarray(1,2,2)) =
1148     &       tempC%var%array1(childarray(1,1,1):childarray(1,2,1))
1149          CASE (2)
1150            child%var%array2(childarray(1,1,2):childarray(1,2,2),
1151     &                       childarray(2,1,2):childarray(2,2,2)) =
1152     &      tempC%var%array2(childarray(1,1,1):childarray(1,2,1),
1153     &                       childarray(2,1,1):childarray(2,2,1))
1154          CASE (3)
1155            child%var%array3(childarray(1,1,2):childarray(1,2,2),
1156     &                       childarray(2,1,2):childarray(2,2,2),
1157     &                       childarray(3,1,2):childarray(3,2,2)) =
1158     &      tempC%var%array3(childarray(1,1,1):childarray(1,2,1),
1159     &                       childarray(2,1,1):childarray(2,2,1),
1160     &                       childarray(3,1,1):childarray(3,2,1))
1161          CASE (4)
1162            child%var%array4(childarray(1,1,2):childarray(1,2,2),
1163     &                       childarray(2,1,2):childarray(2,2,2),
1164     &                       childarray(3,1,2):childarray(3,2,2),
1165     &                       childarray(4,1,2):childarray(4,2,2)) =
1166     &      tempC%var%array4(childarray(1,1,1):childarray(1,2,1),
1167     &                       childarray(2,1,1):childarray(2,2,1),
1168     &                       childarray(3,1,1):childarray(3,2,1),
1169     &                       childarray(4,1,1):childarray(4,2,1))
1170          CASE (5)
1171            child%var%array5(childarray(1,1,2):childarray(1,2,2),
1172     &                       childarray(2,1,2):childarray(2,2,2),
1173     &                       childarray(3,1,2):childarray(3,2,2),
1174     &                       childarray(4,1,2):childarray(4,2,2),
1175     &                       childarray(5,1,2):childarray(5,2,2)) =
1176     &      tempC%var%array5(childarray(1,1,1):childarray(1,2,1),
1177     &                       childarray(2,1,1):childarray(2,2,1),
1178     &                       childarray(3,1,1):childarray(3,2,1),
1179     &                       childarray(4,1,1):childarray(4,2,1),
1180     &                       childarray(5,1,1):childarray(5,2,1))
1181          CASE (6)
1182            child%var%array6(childarray(1,1,2):childarray(1,2,2),
1183     &                       childarray(2,1,2):childarray(2,2,2),
1184     &                       childarray(3,1,2):childarray(3,2,2),
1185     &                       childarray(4,1,2):childarray(4,2,2),
1186     &                       childarray(5,1,2):childarray(5,2,2),
1187     &                       childarray(6,1,2):childarray(6,2,2)) =
1188     &      tempC%var%array6(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     &                       childarray(6,1,1):childarray(6,2,1))
1194          END SELECT
1195          ENDIF
1196C
1197C       
1198      endif
1199
1200        Call Agrif_nbdim_deallocation(tempPextend%var,nbdim)
1201C        deallocate(tempPextend%var)
1202
1203      Call Agrif_nbdim_deallocation(tempC%var,nbdim)
1204     
1205C      Deallocate(tempC % var)
1206      ELSE
1207     
1208C      deallocate(tempPextend%var)
1209
1210      ENDIF
1211C
1212C             
1213C     Deallocations
1214#ifdef AGRIF_MPI       
1215      IF (member) THEN
1216      Call Agrif_nbdim_deallocation(tempP%var,nbdim)
1217C      Deallocate(tempP % var)
1218      endif
1219#endif
1220C
1221C
1222     
1223C
1224C
1225      End Subroutine Agrif_InterpnD 
1226C
1227C
1228C
1229C                 
1230C
1231C     **************************************************************************
1232CCC   Subroutine Agrif_Parentbounds
1233C     **************************************************************************
1234C
1235      Subroutine Agrif_Parentbounds(TYPEinterp,nbdim,indmin,indmax,
1236     &                              s_Parent_temp,
1237     &                              s_Child_temp,s_Child,ds_Child,
1238     &                              s_Parent,ds_Parent,
1239     &                              pttruetab,cetruetab,pttab_Child,
1240     &                              pttab_Parent,posvar,interptab)
1241C
1242CCC   Description:
1243CCC   Subroutine calculating the bounds of the parent grid for the interpolation
1244CCC   of the child grid     
1245C
1246C
1247C     Declarations:
1248C
1249C
1250C     Arguments
1251      INTEGER :: nbdim
1252      INTEGER, DIMENSION(6) :: TypeInterp
1253      INTEGER,DIMENSION(nbdim) :: indmin,indmax
1254      REAL,DIMENSION(nbdim) :: s_Parent_temp,s_child_temp
1255      REAL,DIMENSION(nbdim) :: s_Child,ds_child
1256      REAL,DIMENSION(nbdim) :: s_Parent,ds_Parent
1257      INTEGER,DIMENSION(nbdim) :: pttruetab,cetruetab
1258      INTEGER,DIMENSION(nbdim) :: pttab_Child,pttab_Parent
1259      INTEGER,DIMENSION(nbdim) :: posvar
1260      CHARACTER(6), DIMENSION(nbdim) :: interptab
1261C
1262C     Local variables
1263      INTEGER :: i
1264      REAL,DIMENSION(nbdim) :: dim_newmin,dim_newmax     
1265C
1266      dim_newmin = s_Child + (pttruetab - pttab_Child) * ds_Child
1267      dim_newmax = s_Child + (cetruetab - pttab_Child) * ds_Child
1268     
1269      DO i = 1,nbdim         
1270C     
1271        indmin(i) = pttab_Parent(i) + 
1272     &         agrif_int((dim_newmin(i)-s_Parent(i))/ds_Parent(i))
1273C
1274        indmax(i) = pttab_Parent(i) + 
1275     &                agrif_ceiling((dim_newmax(i)-
1276     &                s_Parent(i))/ds_Parent(i))
1277     
1278C
1279C
1280C       Necessary for the Quadratic interpolation
1281C 
1282
1283        IF ((pttruetab(i) == cetruetab(i)) .AND. 
1284     &                           (posvar(i) == 1)) THEN
1285        ELSEIF (interptab(i) .EQ. 'N') THEN
1286        ELSEIF ( TYPEinterp(i) .eq. Agrif_ppm .or.
1287     &      TYPEinterp(i) .eq. Agrif_eno  .or.
1288     &      TYPEinterp(i) .eq. Agrif_weno) THEN           
1289           indmin(i) = indmin(i) - 2 
1290           indmax(i) = indmax(i) + 2                 
1291        ELSE IF (( TYPEinterp(i) .ne. Agrif_constant )
1292     &        .AND.( TYPEinterp(i) .ne. Agrif_linear )) THEN
1293           indmin(i) = indmin(i) - 1 
1294           indmax(i) = indmax(i) + 1
1295        ENDIF
1296       
1297
1298C       
1299       ENDDO 
1300C
1301        s_Parent_temp = s_Parent + (indmin - pttab_Parent) * ds_Parent
1302C     
1303        s_Child_temp = s_Child + (pttruetab - pttab_Child) * ds_Child
1304C
1305C
1306      Return
1307C
1308C
1309      End Subroutine Agrif_Parentbounds
1310C
1311C
1312C
1313C     **************************************************************************
1314CCC   Subroutine Agrif_Interp_1D_Recursive 
1315C     **************************************************************************
1316C
1317      Subroutine Agrif_Interp_1D_recursive(TypeInterp,tabin,tabout,
1318     &           indmin,indmax, 
1319     &           pttab_child,petab_child,
1320     &           s_child,s_parent,ds_child,ds_parent,nbdim)     
1321C
1322CCC   Description:
1323CCC   Subroutine for the interpolation of a 1D grid variable. 
1324CCC   It calls Agrif_InterpBase. 
1325C
1326C     Declarations:
1327C
1328     
1329C
1330C     Arguments
1331      INTEGER :: nbdim
1332      INTEGER,DIMENSION(1) :: TypeInterp
1333      INTEGER, DIMENSION(nbdim) :: indmin,indmax
1334      INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child
1335      REAL, DIMENSION(nbdim) :: s_child,s_parent
1336      REAL, DIMENSION(nbdim) :: ds_child,ds_parent
1337      REAL, INTENT(IN),DIMENSION(indmin(nbdim):indmax(nbdim)) :: tabin       
1338      REAL, INTENT(OUT),
1339     &      DIMENSION(pttab_child(nbdim):petab_child(nbdim)) :: tabout
1340      INTEGER :: coeffraf
1341C
1342C
1343C     Commentaire perso : nbdim vaut toujours 1 ici. 
1344C
1345      coeffraf = nint(ds_parent(nbdim)/ds_child(nbdim))
1346     
1347      Call Agrif_InterpBase(TypeInterp(1),
1348     &                  tabin(indmin(nbdim):indmax(nbdim)),
1349     &                  tabout(pttab_child(nbdim):petab_child(nbdim)),
1350     &                  indmin(nbdim),indmax(nbdim),           
1351     &                  pttab_child(nbdim),petab_child(nbdim),
1352     &                  s_parent(nbdim),s_child(nbdim),
1353     &                  ds_parent(nbdim),ds_child(nbdim),coeffraf)
1354               
1355C               
1356      Return
1357C
1358C
1359      End Subroutine Agrif_Interp_1D_recursive
1360C
1361C
1362C     
1363C     **************************************************************************
1364CCC   Subroutine Agrif_Interp_2D_Recursive 
1365C     **************************************************************************
1366C
1367      Subroutine Agrif_Interp_2D_recursive(TypeInterp,
1368     &           tabin,tabout,
1369     &           indmin,indmax,   
1370     &           pttab_child,petab_child,
1371     &            s_child, s_parent,
1372     &           ds_child,ds_parent,
1373     &           nbdim)
1374C
1375CCC   Description:
1376CCC   Subroutine for the interpolation of a 2D grid variable. 
1377CCC   It calls Agrif_Interp_1D_recursive and Agrif_InterpBase.   
1378C
1379C     Declarations:
1380C
1381     
1382C     
1383      INTEGER                   :: nbdim
1384      INTEGER,DIMENSION(2)      :: TypeInterp
1385      INTEGER, DIMENSION(nbdim) :: indmin,indmax
1386      INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child
1387      REAL   , DIMENSION(nbdim) ::  s_child, s_parent
1388      REAL   , DIMENSION(nbdim) :: ds_child,ds_parent
1389      REAL   ,INTENT(IN), DIMENSION(
1390     &                indmin(nbdim-1):indmax(nbdim-1),
1391     &                indmin(nbdim):indmax(nbdim)
1392     &                ) :: tabin       
1393      REAL   ,INTENT(OUT), DIMENSION(
1394     &                pttab_child(nbdim-1):petab_child(nbdim-1),
1395     &                pttab_child(nbdim):petab_child(nbdim)
1396     &                ) :: tabout
1397C
1398C     Local variables     
1399      REAL, DIMENSION(pttab_child(nbdim-1):petab_child(nbdim-1),
1400     &                 indmin(nbdim):indmax(nbdim)) :: tabtemp
1401      INTEGER i,j
1402      INTEGER :: coeffraf
1403      REAL   , DIMENSION(
1404     &                pttab_child(nbdim):petab_child(nbdim),
1405     &                pttab_child(nbdim-1):petab_child(nbdim-1)
1406     &                ) :: tabout_trsp
1407      REAL, DIMENSION(indmin(nbdim):indmax(nbdim),
1408     &        pttab_child(nbdim-1):petab_child(nbdim-1)) :: tabtemp_trsp
1409
1410C
1411C
1412C
1413C
1414C     Commentaire perso : nbdim vaut toujours 2 ici.
1415C
1416      coeffraf = nint ( ds_parent(1) / ds_child(1) )
1417      IF((TypeInterp(1) == Agrif_Linear) .AND. (coeffraf /= 1 ) )THEN
1418
1419!---CDIR NEXPAND
1420          IF(.NOT. precomputedone(1) ) call linear1Dprecompute2D(
1421     &          indmax(2)-indmin(2)+1,   
1422     &          indmax(1)-indmin(1)+1,   
1423     &          petab_child(1)-pttab_child(1)+1,
1424     &          s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1425!---CDIR NEXPAND
1426           call linear1daftercompute(tabin,tabtemp,
1427     &          size(tabin), size(tabtemp), 
1428     &          s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1429 
1430      ELSEIF((TypeInterp(1) == Agrif_PPM) .AND. (coeffraf /= 1 ) )THEN
1431!---CDIR NEXPAND
1432          IF(.NOT. precomputedone(1) ) call ppm1Dprecompute2D(
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 ppm1daftercompute(tabin,tabtemp,
1439     &          size(tabin), size(tabtemp), 
1440     &          s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1441
1442      ELSE
1443
1444      do j = indmin(nbdim),indmax(nbdim)
1445C       
1446!---CDIR NEXPAND
1447        Call Agrif_Interp_1D_recursive(TypeInterp(1),
1448     &         tabin(indmin(nbdim-1):indmax(nbdim-1),j),
1449     &         tabtemp(pttab_child(nbdim-1):petab_child(nbdim-1),j),
1450     &         indmin(1:nbdim-1),indmax(1:nbdim-1),
1451     &         pttab_child(1:nbdim-1),petab_child(1:nbdim-1),
1452     &         s_child(1:nbdim-1),s_parent(1:nbdim-1),
1453     &         ds_child(1:nbdim-1),ds_parent(1:nbdim-1),nbdim-1)
1454C       
1455      enddo
1456      ENDIF   
1457
1458      coeffraf = nint(ds_parent(nbdim)/ds_child(nbdim))
1459     
1460      tabtemp_trsp = TRANSPOSE(tabtemp)
1461
1462      IF((TypeInterp(2) == Agrif_Linear) .AND. (coeffraf /= 1 ) )THEN
1463
1464!---CDIR NEXPAND
1465          IF(.NOT. precomputedone(2) ) call linear1Dprecompute2D(
1466     &          petab_child(1)-pttab_child(1)+1,
1467     &          indmax(2)-indmin(2)+1,
1468     &          petab_child(2)-pttab_child(2)+1,
1469     &          s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1470!---CDIR NEXPAND
1471           call linear1daftercompute(tabtemp_trsp,tabout_trsp,
1472     &          size(tabtemp_trsp), size(tabout_trsp),
1473     &          s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1474
1475      ELSEIF((TypeInterp(2) == Agrif_PPM) .AND. (coeffraf /= 1 ) )THEN
1476
1477!---CDIR NEXPAND
1478           IF(.NOT. precomputedone(1) )call ppm1Dprecompute2D(
1479     &          petab_child(1)-pttab_child(1)+1,
1480     &          indmax(2)-indmin(2)+1,
1481     &          petab_child(2)-pttab_child(2)+1,
1482     &          s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1483!---CDIR NEXPAND
1484           call ppm1daftercompute(tabtemp_trsp,tabout_trsp,
1485     &          size(tabtemp_trsp), size(tabout_trsp),
1486     &          s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1487
1488      ELSE
1489      do i=pttab_child(nbdim-1),petab_child(nbdim-1)
1490C
1491!---CDIR NEXPAND
1492        Call Agrif_InterpBase(TypeInterp(2),
1493     &           tabtemp_trsp(indmin(nbdim):indmax(nbdim),i),
1494     &           tabout_trsp(pttab_child(nbdim):petab_child(nbdim),i),
1495     &           indmin(nbdim),indmax(nbdim),
1496     &           pttab_child(nbdim),petab_child(nbdim),
1497     &           s_parent(nbdim),s_child(nbdim),
1498     &           ds_parent(nbdim),ds_child(nbdim),coeffraf)
1499
1500C       
1501      enddo
1502      ENDIF
1503     
1504      tabout = TRANSPOSE(tabout_trsp)
1505C
1506      Return
1507C
1508C
1509      End Subroutine Agrif_Interp_2D_recursive
1510C
1511C
1512C     
1513C     **************************************************************************
1514CCC   Subroutine Agrif_Interp_3D_Recursive 
1515C     **************************************************************************
1516C
1517      Subroutine Agrif_Interp_3D_recursive(TypeInterp,tabin,tabout,
1518     &           indmin,indmax,   
1519     &           pttab_child,petab_child,
1520     &           s_child,s_parent,ds_child,ds_parent,nbdim)
1521C
1522CCC   Description:
1523CCC   Subroutine for the interpolation of a 3D grid variable. 
1524CCC   It calls Agrif_Interp_2D_recursive and Agrif_InterpBase.   
1525C
1526C     Declarations:
1527C
1528     
1529C     
1530      INTEGER :: nbdim
1531      INTEGER,DIMENSION(3) :: TypeInterp
1532      INTEGER, DIMENSION(nbdim) :: indmin,indmax
1533      INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child
1534      REAL, DIMENSION(nbdim) :: s_child,s_parent,ds_child,ds_parent
1535      REAL,INTENT(IN), DIMENSION(indmin(nbdim-2):indmax(nbdim-2),
1536     &                indmin(nbdim-1):indmax(nbdim-1),
1537     &                indmin(nbdim)  :indmax(nbdim)) :: tabin       
1538      REAL,INTENT(OUT),
1539     &        DIMENSION(pttab_child(nbdim-2):petab_child(nbdim-2),
1540     &                pttab_child(nbdim-1):petab_child(nbdim-1),
1541     &                pttab_child(nbdim):petab_child(nbdim)) :: tabout
1542C
1543C     Local variables     
1544      REAL, DIMENSION(pttab_child(nbdim-2):petab_child(nbdim-2),
1545     &                 pttab_child(nbdim-1):petab_child(nbdim-1),
1546     &                 indmin(nbdim):indmax(nbdim)) :: tabtemp
1547      INTEGER i,j,k
1548      INTEGER :: coeffraf, locind_child_left, kdeb
1549C
1550C
1551      coeffraf = nint ( ds_parent(1) / ds_child(1) )
1552      IF((TypeInterp(1) == Agrif_Linear) .AND. (coeffraf/=1))THEN
1553        Call linear1Dprecompute2D(
1554     &          indmax(2)-indmin(2)+1,
1555     &          indmax(1)-indmin(1)+1,
1556     &          petab_child(1)-pttab_child(1)+1,
1557     &          s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1558      precomputedone(1) = .TRUE. 
1559      ELSEIF((TypeInterp(1) == Agrif_PPM) .AND. (coeffraf/=1))THEN
1560        Call ppm1Dprecompute2D(
1561     &          indmax(2)-indmin(2)+1,
1562     &          indmax(1)-indmin(1)+1,
1563     &          petab_child(1)-pttab_child(1)+1,
1564     &          s_parent(1),s_child(1),ds_parent(1),ds_child(1),1)
1565      precomputedone(1) = .TRUE.
1566      ENDIF
1567
1568      coeffraf = nint ( ds_parent(2) / ds_child(2) )
1569      IF((TypeInterp(2) == Agrif_Linear) .AND. (coeffraf/=1)) THEN
1570         Call linear1Dprecompute2D(
1571     &          petab_child(1)-pttab_child(1)+1,
1572     &          indmax(2)-indmin(2)+1,
1573     &          petab_child(2)-pttab_child(2)+1,
1574     &          s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1575      precomputedone(2) = .TRUE. 
1576      ELSEIF((TypeInterp(2) == Agrif_PPM) .AND. (coeffraf/=1)) THEN
1577         Call ppm1Dprecompute2D(
1578     &          petab_child(1)-pttab_child(1)+1,
1579     &          indmax(2)-indmin(2)+1,
1580     &          petab_child(2)-pttab_child(2)+1,
1581     &          s_parent(2),s_child(2),ds_parent(2),ds_child(2),2)
1582      precomputedone(2) = .TRUE.
1583      ENDIF
1584
1585      do k = indmin(nbdim),indmax(nbdim)
1586C       
1587        Call Agrif_Interp_2D_recursive(TypeInterp(1:2),
1588     &         tabin(indmin(nbdim-2):indmax(nbdim-2),
1589     &         indmin(nbdim-1):indmax(nbdim-1),k),
1590     &         tabtemp(pttab_child(nbdim-2):petab_child(nbdim-2),
1591     &         pttab_child(nbdim-1):petab_child(nbdim-1),k),
1592     &         indmin(1:nbdim-1),indmax(1:nbdim-1),
1593     &         pttab_child(1:nbdim-1),petab_child(1:nbdim-1),
1594     &         s_child(1:nbdim-1),s_parent(1:nbdim-1),
1595     &         ds_child(1:nbdim-1),ds_parent(1:nbdim-1),nbdim-1)
1596C       
1597      enddo
1598     
1599      precomputedone(1) = .FALSE.
1600      precomputedone(2) = .FALSE.
1601      coeffraf = nint ( ds_parent(3) / ds_child(3) )
1602
1603      Call Agrif_Compute_nbdim_interp(s_parent(nbdim),s_child(nbdim),
1604     &  ds_parent(nbdim),ds_child(nbdim),coeffraf,locind_child_left)
1605           
1606      IF (coeffraf == 1) THEN
1607     
1608      kdeb = indmin(3)+locind_child_left-2
1609      do k=pttab_child(3),petab_child(3)
1610      kdeb = kdeb + 1
1611      do j = pttab_child(2),petab_child(2)
1612        do i = pttab_child(1),petab_child(1)
1613        tabout(i,j,k) = tabtemp(i,j,kdeb)
1614      enddo
1615      enddo
1616      enddo
1617             
1618      ELSE     
1619C
1620      do j=pttab_child(nbdim-1),petab_child(nbdim-1) 
1621C       
1622        do i=pttab_child(nbdim-2),petab_child(nbdim-2)
1623C
1624          Call Agrif_InterpBase(TypeInterp(3),
1625     &           tabtemp(i,j,indmin(nbdim):indmax(nbdim)),
1626     &           tabout(i,j,pttab_child(nbdim):petab_child(nbdim)),
1627     &           indmin(nbdim),indmax(nbdim),
1628     &           pttab_child(nbdim),petab_child(nbdim),
1629     &           s_parent(nbdim),s_child(nbdim),
1630     &           ds_parent(nbdim),ds_child(nbdim),coeffraf)
1631C
1632        enddo 
1633C       
1634      enddo
1635      ENDIF
1636C
1637      Return
1638C       
1639C
1640      End Subroutine Agrif_Interp_3D_recursive
1641C
1642C
1643C
1644C     **************************************************************************
1645CCC   Subroutine Agrif_Interp_4D_Recursive 
1646C     **************************************************************************
1647C
1648      Subroutine Agrif_Interp_4D_recursive(TypeInterp,tabin,tabout,
1649     &           indmin,indmax,   
1650     &           pttab_child,petab_child,
1651     &           s_child,s_parent,ds_child,ds_parent,nbdim)
1652C
1653CCC   Description:
1654CCC   Subroutine for the interpolation of a 4D grid variable. 
1655CCC   It calls Agrif_Interp_3D_recursive and Agrif_InterpBase.   
1656C
1657C     Declarations:
1658C
1659     
1660C     
1661      INTEGER :: nbdim
1662      INTEGER,DIMENSION(4) :: TypeInterp
1663      INTEGER, DIMENSION(nbdim) :: indmin,indmax
1664      INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child
1665      REAL, DIMENSION(nbdim) :: s_child,s_parent,ds_child,ds_parent
1666      REAL,INTENT(IN), DIMENSION(indmin(nbdim-3):indmax(nbdim-3),
1667     &                indmin(nbdim-2):indmax(nbdim-2),
1668     &                indmin(nbdim-1):indmax(nbdim-1),
1669     &                indmin(nbdim):indmax(nbdim)) :: tabin       
1670      REAL,INTENT(OUT),
1671     &       DIMENSION(pttab_child(nbdim-3):petab_child(nbdim-3),
1672     &                pttab_child(nbdim-2):petab_child(nbdim-2),
1673     &                pttab_child(nbdim-1):petab_child(nbdim-1),
1674     &                pttab_child(nbdim):petab_child(nbdim)) :: tabout
1675C
1676C     Local variables     
1677      REAL, DIMENSION(pttab_child(nbdim-3):petab_child(nbdim-3),
1678     &                 pttab_child(nbdim-2):petab_child(nbdim-2),
1679     &                 pttab_child(nbdim-1):petab_child(nbdim-1), 
1680     &                 indmin(nbdim):indmax(nbdim)) :: tabtemp
1681      INTEGER i,j,k,l
1682      INTEGER :: coeffraf
1683C
1684C
1685      do l = indmin(nbdim),indmax(nbdim)
1686C       
1687        Call Agrif_Interp_3D_recursive(TypeInterp(1:3),
1688     &         tabin(indmin(nbdim-3):indmax(nbdim-3),
1689     &               indmin(nbdim-2):indmax(nbdim-2),
1690     &               indmin(nbdim-1):indmax(nbdim-1),l),
1691     &         tabtemp(pttab_child(nbdim-3):petab_child(nbdim-3),
1692     &         pttab_child(nbdim-2):petab_child(nbdim-2),
1693     &         pttab_child(nbdim-1):petab_child(nbdim-1),l),
1694     &         indmin(1:nbdim-1),indmax(1:nbdim-1),
1695     &         pttab_child(1:nbdim-1),petab_child(1:nbdim-1),
1696     &         s_child(1:nbdim-1),s_parent(1:nbdim-1),
1697     &         ds_child(1:nbdim-1),ds_parent(1:nbdim-1),nbdim-1)
1698C       
1699      enddo
1700C
1701      coeffraf = nint(ds_parent(nbdim)/ds_child(nbdim))
1702     
1703      do k = pttab_child(nbdim-1),petab_child(nbdim-1)
1704C
1705        do j = pttab_child(nbdim-2),petab_child(nbdim-2) 
1706C       
1707          do i = pttab_child(nbdim-3),petab_child(nbdim-3)
1708C
1709            Call Agrif_InterpBase(TypeInterp(4),
1710     &           tabtemp(i,j,k,indmin(nbdim):indmax(nbdim)),
1711     &           tabout(i,j,k,pttab_child(nbdim):petab_child(nbdim)),
1712     &           indmin(nbdim),indmax(nbdim),
1713     &           pttab_child(nbdim),petab_child(nbdim),
1714     &           s_parent(nbdim),s_child(nbdim),
1715     &           ds_parent(nbdim),ds_child(nbdim),coeffraf)
1716C
1717          enddo
1718C
1719        enddo 
1720C       
1721      enddo
1722C
1723      Return
1724C
1725C       
1726      End Subroutine Agrif_Interp_4D_recursive
1727C
1728C
1729C
1730C     **************************************************************************
1731CCC   Subroutine Agrif_Interp_5D_Recursive 
1732C     **************************************************************************
1733C
1734      Subroutine Agrif_Interp_5D_recursive(TypeInterp,tabin,tabout,
1735     &           indmin,indmax,   
1736     &           pttab_child,petab_child,
1737     &           s_child,s_parent,ds_child,ds_parent,nbdim)
1738C
1739CCC   Description:
1740CCC   Subroutine for the interpolation of a 5D grid variable. 
1741CCC   It calls Agrif_Interp_4D_recursive and Agrif_InterpBase.   
1742C
1743C     Declarations:
1744C
1745     
1746C     
1747      INTEGER :: nbdim
1748      INTEGER,DIMENSION(5) :: TypeInterp
1749      INTEGER, DIMENSION(nbdim) :: indmin,indmax
1750      INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child
1751      REAL, DIMENSION(nbdim) :: s_child,s_parent,ds_child,ds_parent
1752      REAL,INTENT(IN), DIMENSION(indmin(nbdim-4):indmax(nbdim-4),
1753     &                indmin(nbdim-3):indmax(nbdim-3),
1754     &                indmin(nbdim-2):indmax(nbdim-2),
1755     &                indmin(nbdim-1):indmax(nbdim-1),
1756     &                indmin(nbdim):indmax(nbdim)) :: tabin 
1757      REAL,INTENT(OUT),
1758     &    DIMENSION(pttab_child(nbdim-4):petab_child(nbdim-4),
1759     &                pttab_child(nbdim-3):petab_child(nbdim-3),
1760     &                pttab_child(nbdim-2):petab_child(nbdim-2),
1761     &                pttab_child(nbdim-1):petab_child(nbdim-1),
1762     &                pttab_child(nbdim):petab_child(nbdim)) :: tabout
1763C
1764C     Local variables     
1765      REAL, DIMENSION(pttab_child(nbdim-4):petab_child(nbdim-4),
1766     &                 pttab_child(nbdim-3):petab_child(nbdim-3),
1767     &                 pttab_child(nbdim-2):petab_child(nbdim-2),
1768     &                 pttab_child(nbdim-1):petab_child(nbdim-1),   
1769     &                 indmin(nbdim):indmax(nbdim)) :: tabtemp
1770      INTEGER i,j,k,l,m
1771      INTEGER :: coeffraf
1772C
1773C
1774      do m = indmin(nbdim),indmax(nbdim)
1775C       
1776        Call Agrif_Interp_4D_recursive(TypeInterp(1:4),
1777     &         tabin(indmin(nbdim-4):indmax(nbdim-4),
1778     &               indmin(nbdim-3):indmax(nbdim-3),
1779     &               indmin(nbdim-2):indmax(nbdim-2),
1780     &               indmin(nbdim-1):indmax(nbdim-1),m),
1781     &         tabtemp(pttab_child(nbdim-4):petab_child(nbdim-4),
1782     &                 pttab_child(nbdim-3):petab_child(nbdim-3),
1783     &                 pttab_child(nbdim-2):petab_child(nbdim-2),
1784     &                 pttab_child(nbdim-1):petab_child(nbdim-1),m),
1785     &         indmin(1:nbdim-1),indmax(1:nbdim-1),
1786     &         pttab_child(1:nbdim-1),petab_child(1:nbdim-1),
1787     &         s_child(1:nbdim-1),s_parent(1:nbdim-1),
1788     &         ds_child(1:nbdim-1),ds_parent(1:nbdim-1),nbdim-1)
1789C       
1790      enddo
1791     
1792      coeffraf = nint(ds_parent(nbdim)/ds_child(nbdim))
1793C
1794      do l = pttab_child(nbdim-1),petab_child(nbdim-1) 
1795C
1796        do k = pttab_child(nbdim-2),petab_child(nbdim-2)
1797C
1798          do j = pttab_child(nbdim-3),petab_child(nbdim-3) 
1799C       
1800            do i = pttab_child(nbdim-4),petab_child(nbdim-4)
1801C
1802              Call Agrif_InterpBase(TypeInterp(5),
1803     &             tabtemp(i,j,k,l,indmin(nbdim):indmax(nbdim)),
1804     &                    tabout(i,j,k,l,
1805     &             pttab_child(nbdim):petab_child(nbdim)),
1806     &             indmin(nbdim),indmax(nbdim),
1807     &             pttab_child(nbdim),petab_child(nbdim),
1808     &             s_parent(nbdim),s_child(nbdim),
1809     &             ds_parent(nbdim),ds_child(nbdim),coeffraf)
1810C
1811            enddo
1812C
1813          enddo
1814C
1815        enddo 
1816C       
1817      enddo
1818C
1819C
1820      Return
1821C
1822C       
1823      End Subroutine Agrif_Interp_5D_recursive
1824C
1825C
1826C
1827C     **************************************************************************
1828CCC   Subroutine Agrif_Interp_6D_Recursive 
1829C     **************************************************************************
1830C
1831      Subroutine Agrif_Interp_6D_recursive(TypeInterp,tabin,tabout,
1832     &           indmin,indmax,   
1833     &           pttab_child,petab_child,
1834     &           s_child,s_parent,ds_child,ds_parent,nbdim)
1835C
1836CCC   Description:
1837CCC   Subroutine for the interpolation of a 6D grid variable. 
1838CCC   It calls Agrif_Interp_4D_recursive and Agrif_InterpBase.   
1839C
1840C     Declarations:
1841C
1842     
1843C     
1844      INTEGER :: nbdim
1845      INTEGER,DIMENSION(6) :: TypeInterp
1846      INTEGER, DIMENSION(nbdim) :: indmin,indmax
1847      INTEGER, DIMENSION(nbdim) :: pttab_child,petab_child
1848      REAL, DIMENSION(nbdim) :: s_child,s_parent,ds_child,ds_parent
1849      REAL,INTENT(IN), DIMENSION(indmin(nbdim-5):indmax(nbdim-5),
1850     &                indmin(nbdim-4):indmax(nbdim-4),
1851     &                indmin(nbdim-3):indmax(nbdim-3), 
1852     &                indmin(nbdim-2):indmax(nbdim-2),
1853     &                indmin(nbdim-1):indmax(nbdim-1),
1854     &                indmin(nbdim):indmax(nbdim)) :: tabin       
1855      REAL,INTENT(OUT),
1856     &    DIMENSION(pttab_child(nbdim-5):petab_child(nbdim-5),
1857     &                pttab_child(nbdim-4):petab_child(nbdim-4),
1858     &                pttab_child(nbdim-3):petab_child(nbdim-3),
1859     &                pttab_child(nbdim-2):petab_child(nbdim-2),
1860     &                pttab_child(nbdim-1):petab_child(nbdim-1),
1861     &                pttab_child(nbdim):petab_child(nbdim)) :: tabout
1862C
1863C     Local variables     
1864      REAL, DIMENSION(pttab_child(nbdim-5):petab_child(nbdim-5),
1865     &                 pttab_child(nbdim-4):petab_child(nbdim-4),
1866     &                 pttab_child(nbdim-3):petab_child(nbdim-3),
1867     &                 pttab_child(nbdim-2):petab_child(nbdim-2),   
1868     &                 pttab_child(nbdim-1):petab_child(nbdim-1),   
1869     &                 indmin(nbdim):indmax(nbdim)) :: tabtemp
1870      INTEGER i,j,k,l,m,n
1871      INTEGER :: coeffraf
1872C
1873C       
1874C
1875      do n = indmin(nbdim),indmax(nbdim)
1876C       
1877        Call Agrif_Interp_5D_recursive(TypeInterp(1:5),
1878     &         tabin(indmin(nbdim-5):indmax(nbdim-5),
1879     &               indmin(nbdim-4):indmax(nbdim-4),
1880     &               indmin(nbdim-3):indmax(nbdim-3),
1881     &               indmin(nbdim-2):indmax(nbdim-2),
1882     &               indmin(nbdim-1):indmax(nbdim-1),n),
1883     &         tabtemp(pttab_child(nbdim-5):petab_child(nbdim-5),
1884     &                 pttab_child(nbdim-4):petab_child(nbdim-4),
1885     &                 pttab_child(nbdim-3):petab_child(nbdim-3),
1886     &                 pttab_child(nbdim-2):petab_child(nbdim-2),
1887     &                 pttab_child(nbdim-1):petab_child(nbdim-1),n),
1888     &         indmin(1:nbdim-1),indmax(1:nbdim-1),
1889     &         pttab_child(1:nbdim-1),petab_child(1:nbdim-1),
1890     &         s_child(1:nbdim-1),s_parent(1:nbdim-1),
1891     &         ds_child(1:nbdim-1),ds_parent(1:nbdim-1),nbdim-1)
1892C       
1893      enddo
1894     
1895      coeffraf = nint(ds_parent(nbdim)/ds_child(nbdim))
1896C
1897      do m = pttab_child(nbdim-1),petab_child(nbdim-1) 
1898      do l = pttab_child(nbdim-2),petab_child(nbdim-2) 
1899C
1900        do k = pttab_child(nbdim-3),petab_child(nbdim-3)
1901C
1902          do j = pttab_child(nbdim-4),petab_child(nbdim-4) 
1903C       
1904            do i = pttab_child(nbdim-5),petab_child(nbdim-5)
1905C
1906              Call Agrif_InterpBase(TypeInterp(6),
1907     &             tabtemp(i,j,k,l,m,indmin(nbdim):indmax(nbdim)),
1908     &                    tabout(i,j,k,l,m,
1909     &                    pttab_child(nbdim):petab_child(nbdim)),
1910     &             indmin(nbdim),indmax(nbdim),
1911     &             pttab_child(nbdim),petab_child(nbdim),
1912     &             s_parent(nbdim),s_child(nbdim),
1913     &             ds_parent(nbdim),ds_child(nbdim),coeffraf)
1914C
1915            enddo
1916C
1917          enddo
1918C
1919        enddo 
1920C       
1921      enddo
1922      enddo
1923C               
1924C
1925      Return
1926C
1927C       
1928      End Subroutine Agrif_Interp_6D_recursive
1929C
1930C
1931C
1932C     **************************************************************************
1933CCC   Subroutine Agrif_InterpBase 
1934C     **************************************************************************
1935C 
1936      Subroutine Agrif_InterpBase(TypeInterp,
1937     &                           parenttab,childtab,
1938     &                           indmin,indmax,pttab_child,petab_child,
1939     &                           s_parent,s_child,ds_parent,ds_child,
1940     &                           coeffraf)   
1941C
1942CCC   Description:
1943CCC   Subroutine calling the interpolation method chosen by the user (linear, 
1944CCC   lagrange or spline). 
1945C
1946C     Declarations:
1947C
1948     
1949C
1950      INTEGER                :: TypeInterp
1951      INTEGER :: indmin,indmax
1952      INTEGER :: pttab_child,petab_child
1953      REAL,INTENT(IN),DIMENSION(indmin:indmax)           :: parenttab       
1954      REAL,INTENT(OUT),DIMENSION(pttab_child:petab_child) :: childtab       
1955      REAL    :: s_parent,s_child,ds_parent,ds_child 
1956      INTEGER :: coeffraf
1957C 
1958C
1959       IF ((indmin == indmax).AND.(pttab_child == petab_child)) THEN
1960         childtab(pttab_child) = parenttab(indmin)
1961       ELSEIF (TYPEinterp .EQ. AGRIF_LINEAR) then   
1962C
1963C         Linear interpolation 
1964 
1965          Call linear1D
1966     &         (parenttab,childtab,
1967     &          indmax-indmin+1,petab_child-pttab_child+1,
1968     &          s_parent,s_child,ds_parent,ds_child)
1969C         
1970      elseif ( TYPEinterp .EQ. AGRIF_PPM ) then
1971
1972          Call ppm1D
1973     &         (parenttab,childtab,
1974     &         indmax-indmin+1,petab_child-pttab_child+1,
1975     &         s_parent,s_child,ds_parent,ds_child)
1976C
1977
1978        elseif (TYPEinterp .EQ. AGRIF_LAGRANGE) then
1979C         
1980C         Lagrange interpolation   
1981          Call lagrange1D
1982     &        (parenttab,childtab,
1983     &         indmax-indmin+1,petab_child-pttab_child+1,
1984     &         s_parent,s_child,ds_parent,ds_child)
1985C           
1986        elseif (TYPEinterp .EQ. AGRIF_ENO) then
1987C         
1988C         Eno interpolation
1989          Call eno1D
1990     &         (parenttab,childtab,
1991     &         indmax-indmin+1,petab_child-pttab_child+1,
1992     &         s_parent,s_child,ds_parent,ds_child)
1993C       
1994        elseif (TYPEinterp .EQ. AGRIF_WENO) then
1995C         
1996C         Eno interpolation
1997          Call weno1D
1998     &         (parenttab,childtab,
1999     &         indmax-indmin+1,petab_child-pttab_child+1,
2000     &         s_parent,s_child,ds_parent,ds_child)
2001C           
2002        Else if (TYPEinterp .EQ. AGRIF_LINEARCONSERV) then
2003C         
2004C         Linear conservative interpolation
2005         
2006          Call linear1Dconserv
2007     &         (parenttab,childtab,
2008     &         indmax-indmin+1,petab_child-pttab_child+1,
2009     &         s_parent,s_child,ds_parent,ds_child)   
2010C             
2011        Else if (TYPEinterp .EQ. AGRIF_LINEARCONSERVLIM) then
2012C         
2013C         Linear conservative interpolation
2014         
2015          Call linear1Dconservlim
2016     &         (parenttab,childtab,
2017     &         indmax-indmin+1,petab_child-pttab_child+1,
2018     &         s_parent,s_child,ds_parent,ds_child)         
2019C             
2020        elseif (TYPEinterp .EQ. AGRIF_CONSTANT) then
2021C         
2022          Call constant1D
2023     &         (parenttab,childtab,
2024     &         indmax-indmin+1,petab_child-pttab_child+1,
2025     &         s_parent,s_child,ds_parent,ds_child)
2026C             
2027      endif
2028C
2029C     
2030      End Subroutine Agrif_InterpBase 
2031C
2032
2033      Subroutine Agrif_Compute_nbdim_interp(s_parent,s_child,
2034     &  ds_parent,ds_child,coeffraf,locind_child_left)
2035      real :: s_parent,s_child,ds_parent,ds_child
2036      integer :: coeffraf,locind_child_left
2037     
2038      coeffraf = nint(ds_parent/ds_child)
2039      locind_child_left = 1 + agrif_int((s_child-s_parent)/ds_parent)
2040      End Subroutine Agrif_Compute_nbdim_interp
2041C         
2042
2043      Subroutine Agrif_Find_list_interp(list_interp,pttab,petab,
2044     &                          pttab_Child,pttab_Parent,nbdim,
2045     &                          indmin,indmax,indminglob,
2046     &      indmaxglob,indminglob2,indmaxglob2,parentarray,
2047     &       pttruetab,cetruetab,member,memberin,
2048     &      find_list_interp
2049#if defined AGRIF_MPI
2050     &     ,tab4t,memberinall
2051#endif
2052     &    )     
2053      TYPE(Agrif_List_Interp_Loc), Pointer :: list_interp
2054      INTEGER :: nbdim
2055      INTEGER,DIMENSION(nbdim)   :: pttab,petab,pttab_Child,pttab_Parent
2056      LOGICAL :: find_list_interp
2057      Type(Agrif_List_Interp_loc), Pointer :: parcours
2058      INTEGER,DIMENSION(nbdim)   :: indmin,indmax 
2059      INTEGER,DIMENSION(nbdim)   :: indminglob,indmaxglob
2060      INTEGER,DIMENSION(nbdim)   :: pttruetab,cetruetab     
2061      INTEGER,DIMENSION(nbdim)   :: indminglob2,indmaxglob2 
2062      INTEGER,DIMENSION(nbdim,2,2) :: parentarray
2063      LOGICAL :: member, memberin
2064      INTEGER :: i
2065#ifdef AGRIF_MPI
2066C
2067      INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1,4) :: tab4t
2068      LOGICAL, DIMENSION(0:Agrif_Nbprocs-1) :: memberinall
2069#endif
2070                   
2071      find_list_interp = .FALSE.
2072     
2073      parcours => list_interp
2074      Find_loop :   Do While (associated(parcours))
2075        Do i=1,nbdim
2076          IF ((pttab(i) /= parcours%interp_loc%pttab(i)).OR.
2077     &        (petab(i) /= parcours%interp_loc%petab(i)).OR.
2078     &        (pttab_child(i) /= parcours%interp_loc%pttab_child(i)).OR.
2079     &        (pttab_parent(i) /= parcours%interp_loc%pttab_parent(i)))
2080     &               THEN
2081            parcours=>parcours%suiv
2082            Cycle Find_loop
2083          ENDIF
2084        EndDo
2085C        print *,'ok trouve'
2086        indmin = parcours%interp_loc%indmin(1:nbdim)
2087        indmax = parcours%interp_loc%indmax(1:nbdim)
2088       
2089        pttruetab = parcours%interp_loc%pttruetab(1:nbdim)
2090        cetruetab = parcours%interp_loc%cetruetab(1:nbdim)
2091               
2092#if !defined AGRIF_MPI
2093        indminglob = parcours%interp_loc%indminglob(1:nbdim)
2094        indmaxglob = parcours%interp_loc%indmaxglob(1:nbdim)
2095#else
2096        indminglob2 = parcours%interp_loc%indminglob2(1:nbdim)
2097        indmaxglob2 = parcours%interp_loc%indmaxglob2(1:nbdim)
2098        parentarray = parcours%interp_loc%parentarray(1:nbdim,:,:)
2099        member = parcours%interp_loc%member
2100        tab4t = parcours%interp_loc%tab4t(1:nbdim,0:Agrif_Nbprocs-1,1:4)
2101        memberinall = parcours%interp_loc%memberinall(0:Agrif_Nbprocs-1)
2102#endif       
2103        memberin = parcours%interp_loc%memberin
2104        find_list_interp = .TRUE.   
2105        Exit Find_loop
2106      End Do Find_loop 
2107                             
2108      End Subroutine Agrif_Find_list_interp   
2109     
2110      Subroutine Agrif_AddTo_list_interp(list_interp,pttab,petab,
2111     &                          pttab_Child,pttab_Parent,indmin,indmax,
2112     &                          indminglob,indmaxglob,
2113     &                          indminglob2,indmaxglob2,
2114     &                          parentarray,pttruetab,cetruetab,
2115     &                          member,memberin,nbdim
2116#if defined AGRIF_MPI
2117     &      ,tab4t,memberinall
2118#endif
2119     &    )
2120         
2121      TYPE(Agrif_List_Interp_Loc), Pointer :: list_interp
2122      INTEGER :: nbdim
2123      INTEGER,DIMENSION(nbdim)   :: pttab,petab,pttab_Child,pttab_Parent
2124      INTEGER,DIMENSION(nbdim)   :: indmin,indmax
2125      INTEGER,DIMENSION(nbdim)   :: indminglob,indmaxglob
2126      INTEGER,DIMENSION(nbdim)   :: indminglob2,indmaxglob2
2127      INTEGER,DIMENSION(nbdim)   :: pttruetab,cetruetab
2128      INTEGER,DIMENSION(nbdim,2,2) :: parentarray
2129      LOGICAL :: member, memberin
2130#ifdef AGRIF_MPI
2131C
2132      INTEGER,DIMENSION(nbdim,0:Agrif_Nbprocs-1,4) :: tab4t
2133      LOGICAL,DIMENSION(0:Agrif_Nbprocs-1) :: memberinall
2134#endif                   
2135      Type(Agrif_List_Interp_loc), Pointer :: parcours
2136           
2137       Allocate(parcours)
2138      Allocate(parcours%interp_loc)
2139     
2140      parcours%interp_loc%pttab(1:nbdim) = pttab(1:nbdim)
2141      parcours%interp_loc%petab(1:nbdim) = petab(1:nbdim)
2142      parcours%interp_loc%pttab_child(1:nbdim) = pttab_child(1:nbdim)
2143      parcours%interp_loc%pttab_parent(1:nbdim) = pttab_parent(1:nbdim)
2144 
2145 
2146      parcours%interp_loc%indmin(1:nbdim) = indmin(1:nbdim)
2147      parcours%interp_loc%indmax(1:nbdim) = indmax(1:nbdim)
2148
2149      parcours%interp_loc%memberin = memberin
2150#if !defined AGRIF_MPI
2151      parcours%interp_loc%indminglob(1:nbdim) = indminglob(1:nbdim)
2152      parcours%interp_loc%indmaxglob(1:nbdim) = indmaxglob(1:nbdim)
2153#else
2154      parcours%interp_loc%indminglob2(1:nbdim) = indminglob2(1:nbdim)
2155      parcours%interp_loc%indmaxglob2(1:nbdim) = indmaxglob2(1:nbdim)
2156      parcours%interp_loc%parentarray(1:nbdim,:,:) 
2157     &       = parentarray(1:nbdim,:,:)
2158      parcours%interp_loc%member = member
2159      Allocate(parcours%interp_loc%tab4t(nbdim,0:Agrif_Nbprocs-1,4))
2160      Allocate(parcours%interp_loc%memberinall(0:Agrif_Nbprocs-1))
2161      parcours%interp_loc%tab4t=tab4t   
2162      parcours%interp_loc%memberinall=memberinall   
2163#endif     
2164
2165      parcours%interp_loc%pttruetab(1:nbdim) = pttruetab(1:nbdim)
2166      parcours%interp_loc%cetruetab(1:nbdim) = cetruetab(1:nbdim)
2167     
2168      parcours%suiv => list_interp
2169     
2170      list_interp => parcours
2171      End Subroutine Agrif_Addto_list_interp 
2172               
2173      End Module Agrif_Interpolation
Note: See TracBrowser for help on using the repository browser.